1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine initialize(itime,ldt,up,vp,wp, & |
---|
5 | usigold,vsigold,wsigold,xt,yt,zt,icbt) |
---|
6 | ! i i o o o |
---|
7 | ! o o o i i i o |
---|
8 | !***************************************************************************** |
---|
9 | ! * |
---|
10 | ! Calculation of trajectories utilizing a zero-acceleration scheme. The time* |
---|
11 | ! step is determined by the Courant-Friedrichs-Lewy (CFL) criterion. This * |
---|
12 | ! means that the time step must be so small that the displacement within * |
---|
13 | ! this time step is smaller than 1 grid distance. Additionally, a temporal * |
---|
14 | ! CFL criterion is introduced: the time step must be smaller than the time * |
---|
15 | ! interval of the wind fields used for interpolation. * |
---|
16 | ! For random walk simulations, these are the only time step criteria. * |
---|
17 | ! For the other options, the time step is also limited by the Lagrangian * |
---|
18 | ! time scale. * |
---|
19 | ! * |
---|
20 | ! Author: A. Stohl * |
---|
21 | ! * |
---|
22 | ! 16 December 1997 * |
---|
23 | ! * |
---|
24 | ! Literature: * |
---|
25 | ! * |
---|
26 | !***************************************************************************** |
---|
27 | ! * |
---|
28 | ! Variables: * |
---|
29 | ! h [m] Mixing height * |
---|
30 | ! lwindinterv [s] time interval between two wind fields * |
---|
31 | ! itime [s] current temporal position * |
---|
32 | ! ldt [s] Suggested time step for next integration * |
---|
33 | ! ladvance [s] Total integration time period * |
---|
34 | ! rannumb(maxrand) normally distributed random variables * |
---|
35 | ! up,vp,wp random velocities due to turbulence * |
---|
36 | ! usig,vsig,wsig uncertainties of wind velocities due to interpolation * |
---|
37 | ! usigold,vsigold,wsigold like usig, etc., but for the last time step * |
---|
38 | ! xt,yt,zt Next time step's spatial position of trajectory * |
---|
39 | ! * |
---|
40 | ! * |
---|
41 | ! Constants: * |
---|
42 | ! cfl factor, by which the time step has to be smaller than * |
---|
43 | ! the spatial CFL-criterion * |
---|
44 | ! cflt factor, by which the time step has to be smaller than * |
---|
45 | ! the temporal CFL-criterion * |
---|
46 | ! * |
---|
47 | !***************************************************************************** |
---|
48 | |
---|
49 | use par_mod |
---|
50 | use com_mod |
---|
51 | use interpol_mod |
---|
52 | use hanna_mod |
---|
53 | use random_mod, only: ran3 |
---|
54 | |
---|
55 | implicit none |
---|
56 | |
---|
57 | integer :: itime |
---|
58 | integer :: ldt,nrand |
---|
59 | integer(kind=2) :: icbt |
---|
60 | real :: zt,dz,dz1,dz2,up,vp,wp,usigold,vsigold,wsigold |
---|
61 | real(kind=dp) :: xt,yt |
---|
62 | save idummy |
---|
63 | |
---|
64 | integer :: idummy = -7 |
---|
65 | |
---|
66 | icbt=1 ! initialize particle to "no reflection" |
---|
67 | |
---|
68 | nrand=int(ran3(idummy)*real(maxrand-1))+1 |
---|
69 | |
---|
70 | |
---|
71 | !****************************** |
---|
72 | ! 2. Interpolate necessary data |
---|
73 | !****************************** |
---|
74 | |
---|
75 | ! Compute maximum mixing height around particle position |
---|
76 | !******************************************************* |
---|
77 | |
---|
78 | ix=int(xt) |
---|
79 | jy=int(yt) |
---|
80 | ixp=ix+1 |
---|
81 | jyp=jy+1 |
---|
82 | |
---|
83 | h=max(hmix(ix ,jy,1,memind(1)), & |
---|
84 | hmix(ixp,jy ,1,memind(1)), & |
---|
85 | hmix(ix ,jyp,1,memind(1)), & |
---|
86 | hmix(ixp,jyp,1,memind(1)), & |
---|
87 | hmix(ix ,jy ,1,memind(2)), & |
---|
88 | hmix(ixp,jy ,1,memind(2)), & |
---|
89 | hmix(ix ,jyp,1,memind(2)), & |
---|
90 | hmix(ixp,jyp,1,memind(2))) |
---|
91 | |
---|
92 | zeta=zt/h |
---|
93 | |
---|
94 | |
---|
95 | !************************************************************* |
---|
96 | ! If particle is in the PBL, interpolate once and then make a |
---|
97 | ! time loop until end of interval is reached |
---|
98 | !************************************************************* |
---|
99 | |
---|
100 | if (zeta.le.1.) then |
---|
101 | |
---|
102 | call interpol_all(itime,real(xt),real(yt),zt) |
---|
103 | |
---|
104 | |
---|
105 | ! Vertical interpolation of u,v,w,rho and drhodz |
---|
106 | !*********************************************** |
---|
107 | |
---|
108 | ! Vertical distance to the level below and above current position |
---|
109 | ! both in terms of (u,v) and (w) fields |
---|
110 | !**************************************************************** |
---|
111 | |
---|
112 | dz1=zt-height(indz) |
---|
113 | dz2=height(indzp)-zt |
---|
114 | dz=1./(dz1+dz2) |
---|
115 | |
---|
116 | u=(dz1*uprof(indzp)+dz2*uprof(indz))*dz |
---|
117 | v=(dz1*vprof(indzp)+dz2*vprof(indz))*dz |
---|
118 | w=(dz1*wprof(indzp)+dz2*wprof(indz))*dz |
---|
119 | |
---|
120 | |
---|
121 | ! Compute the turbulent disturbances |
---|
122 | |
---|
123 | ! Determine the sigmas and the timescales |
---|
124 | !**************************************** |
---|
125 | |
---|
126 | if (turbswitch) then |
---|
127 | call hanna(zt) |
---|
128 | else |
---|
129 | call hanna1(zt) |
---|
130 | endif |
---|
131 | |
---|
132 | |
---|
133 | ! Determine the new diffusivity velocities |
---|
134 | !***************************************** |
---|
135 | |
---|
136 | if (nrand+2.gt.maxrand) nrand=1 |
---|
137 | up=rannumb(nrand)*sigu |
---|
138 | vp=rannumb(nrand+1)*sigv |
---|
139 | wp=rannumb(nrand+2) |
---|
140 | if (.not.turbswitch) then ! modified by mc |
---|
141 | wp=wp*sigw |
---|
142 | else if (cblflag.eq.1) then ! modified by mc |
---|
143 | if(-h/ol.gt.5) then |
---|
144 | !if (ol.lt.0.) then |
---|
145 | !if (ol.gt.0.) then !by mc : only for test correct is lt.0 |
---|
146 | call initialize_cbl_vel(idummy,zt,ust,wst,h,sigw,wp,ol) |
---|
147 | else |
---|
148 | wp=wp*sigw |
---|
149 | end if |
---|
150 | end if |
---|
151 | |
---|
152 | |
---|
153 | ! Determine time step for next integration |
---|
154 | !***************************************** |
---|
155 | |
---|
156 | if (turbswitch) then |
---|
157 | ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), & |
---|
158 | 0.5/abs(dsigwdz),600.)*ctl) |
---|
159 | else |
---|
160 | ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5),600.)*ctl) |
---|
161 | endif |
---|
162 | ldt=max(ldt,mintime) |
---|
163 | |
---|
164 | |
---|
165 | usig=(usigprof(indzp)+usigprof(indz))/2. |
---|
166 | vsig=(vsigprof(indzp)+vsigprof(indz))/2. |
---|
167 | wsig=(wsigprof(indzp)+wsigprof(indz))/2. |
---|
168 | |
---|
169 | else |
---|
170 | |
---|
171 | |
---|
172 | |
---|
173 | !********************************************************** |
---|
174 | ! For all particles that are outside the PBL, make a single |
---|
175 | ! time step. Only horizontal turbulent disturbances are |
---|
176 | ! calculated. Vertical disturbances are reset. |
---|
177 | !********************************************************** |
---|
178 | |
---|
179 | |
---|
180 | ! Interpolate the wind |
---|
181 | !********************* |
---|
182 | |
---|
183 | call interpol_wind(itime,real(xt),real(yt),zt) |
---|
184 | |
---|
185 | |
---|
186 | ! Compute everything for above the PBL |
---|
187 | |
---|
188 | ! Assume constant turbulent perturbations |
---|
189 | !**************************************** |
---|
190 | |
---|
191 | ldt=abs(lsynctime) |
---|
192 | |
---|
193 | if (nrand+1.gt.maxrand) nrand=1 |
---|
194 | up=rannumb(nrand)*0.3 |
---|
195 | vp=rannumb(nrand+1)*0.3 |
---|
196 | nrand=nrand+2 |
---|
197 | wp=0. |
---|
198 | sigw=0. |
---|
199 | |
---|
200 | endif |
---|
201 | |
---|
202 | !**************************************************************** |
---|
203 | ! Add mesoscale random disturbances |
---|
204 | ! This is done only once for the whole lsynctime interval to save |
---|
205 | ! computation time |
---|
206 | !**************************************************************** |
---|
207 | |
---|
208 | |
---|
209 | ! It is assumed that the average interpolation error is 1/2 sigma |
---|
210 | ! of the surrounding points, autocorrelation time constant is |
---|
211 | ! 1/2 of time interval between wind fields |
---|
212 | !**************************************************************** |
---|
213 | |
---|
214 | if (nrand+2.gt.maxrand) nrand=1 |
---|
215 | usigold=rannumb(nrand)*usig |
---|
216 | vsigold=rannumb(nrand+1)*vsig |
---|
217 | wsigold=rannumb(nrand+2)*wsig |
---|
218 | |
---|
219 | end subroutine initialize |
---|