1 | !********************************************************************** |
---|
2 | ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 * |
---|
3 | ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, * |
---|
4 | ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann * |
---|
5 | ! * |
---|
6 | ! This file is part of FLEXPART. * |
---|
7 | ! * |
---|
8 | ! FLEXPART is free software: you can redistribute it and/or modify * |
---|
9 | ! it under the terms of the GNU General Public License as published by* |
---|
10 | ! the Free Software Foundation, either version 3 of the License, or * |
---|
11 | ! (at your option) any later version. * |
---|
12 | ! * |
---|
13 | ! FLEXPART is distributed in the hope that it will be useful, * |
---|
14 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
15 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
16 | ! GNU General Public License for more details. * |
---|
17 | ! * |
---|
18 | ! You should have received a copy of the GNU General Public License * |
---|
19 | ! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
20 | !********************************************************************** |
---|
21 | |
---|
22 | subroutine advance(itime,nrelpoint,ldt,up,vp,wp, & |
---|
23 | usigold,vsigold,wsigold,nstop,xt,yt,zt,prob,icbt) |
---|
24 | ! i i i/oi/oi/o |
---|
25 | ! i/o i/o i/o o i/oi/oi/o i/o i/o |
---|
26 | !***************************************************************************** |
---|
27 | ! * |
---|
28 | ! Calculation of turbulent particle trajectories utilizing a * |
---|
29 | ! zero-acceleration scheme, which is corrected by a numerically more * |
---|
30 | ! accurate Petterssen scheme whenever possible. * |
---|
31 | ! * |
---|
32 | ! Particle positions are read in, incremented, and returned to the calling * |
---|
33 | ! program. * |
---|
34 | ! * |
---|
35 | ! In different regions of the atmosphere (PBL vs. free troposphere), * |
---|
36 | ! different parameters are needed for advection, parameterizing turbulent * |
---|
37 | ! velocities, etc. For efficiency, different interpolation routines have * |
---|
38 | ! been written for these different cases, with the disadvantage that there * |
---|
39 | ! exist several routines doing almost the same. They all share the * |
---|
40 | ! included file 'interpol_mod'. The following * |
---|
41 | ! interpolation routines are used: * |
---|
42 | ! * |
---|
43 | ! interpol_all(_nests) interpolates everything (called inside the PBL) * |
---|
44 | ! interpol_misslev(_nests) if a particle moves vertically in the PBL, * |
---|
45 | ! additional parameters are interpolated if it * |
---|
46 | ! crosses a model level * |
---|
47 | ! interpol_wind(_nests) interpolates the wind and determines the * |
---|
48 | ! standard deviation of the wind (called outside * |
---|
49 | ! PBL) also interpolates potential vorticity * |
---|
50 | ! interpol_wind_short(_nests) only interpolates the wind (needed for the * |
---|
51 | ! Petterssen scheme) * |
---|
52 | ! interpol_vdep(_nests) interpolates deposition velocities * |
---|
53 | ! * |
---|
54 | ! * |
---|
55 | ! Author: A. Stohl * |
---|
56 | ! * |
---|
57 | ! 16 December 1997 * |
---|
58 | ! * |
---|
59 | ! Changes: * |
---|
60 | ! * |
---|
61 | ! 8 April 2000: Deep convection parameterization * |
---|
62 | ! * |
---|
63 | ! May 2002: Petterssen scheme introduced * |
---|
64 | ! * |
---|
65 | !***************************************************************************** |
---|
66 | ! * |
---|
67 | ! Variables: * |
---|
68 | ! icbt 1 if particle not transferred to forbidden state, * |
---|
69 | ! else -1 * |
---|
70 | ! dawsave accumulated displacement in along-wind direction * |
---|
71 | ! dcwsave accumulated displacement in cross-wind direction * |
---|
72 | ! dxsave accumulated displacement in longitude * |
---|
73 | ! dysave accumulated displacement in latitude * |
---|
74 | ! h [m] Mixing height * |
---|
75 | ! lwindinterv [s] time interval between two wind fields * |
---|
76 | ! itime [s] time at which this subroutine is entered * |
---|
77 | ! itimec [s] actual time, which is incremented in this subroutine * |
---|
78 | ! href [m] height for which dry deposition velocity is calculated * |
---|
79 | ! ladvance [s] Total integration time period * |
---|
80 | ! ldirect 1 forward, -1 backward * |
---|
81 | ! ldt [s] Time step for the next integration * |
---|
82 | ! lsynctime [s] Synchronisation interval of FLEXPART * |
---|
83 | ! ngrid index which grid is to be used * |
---|
84 | ! nrand index for a variable to be picked from rannumb * |
---|
85 | ! nstop if > 1 particle has left domain and must be stopped * |
---|
86 | ! prob probability of absorption due to dry deposition * |
---|
87 | ! rannumb(maxrand) normally distributed random variables * |
---|
88 | ! rhoa air density * |
---|
89 | ! rhograd vertical gradient of the air density * |
---|
90 | ! up,vp,wp random velocities due to turbulence (along wind, cross * |
---|
91 | ! wind, vertical wind * |
---|
92 | ! usig,vsig,wsig mesoscale wind fluctuations * |
---|
93 | ! usigold,vsigold,wsigold like usig, etc., but for the last time step * |
---|
94 | ! vdepo Deposition velocities for all species * |
---|
95 | ! xt,yt,zt Particle position * |
---|
96 | ! * |
---|
97 | !***************************************************************************** |
---|
98 | |
---|
99 | use point_mod |
---|
100 | use par_mod |
---|
101 | use com_mod |
---|
102 | use interpol_mod |
---|
103 | use hanna_mod |
---|
104 | use cmapf_mod |
---|
105 | use random_mod, only: ran3 |
---|
106 | |
---|
107 | implicit none |
---|
108 | |
---|
109 | real(kind=dp) :: xt,yt |
---|
110 | real :: zt,xts,yts,weight |
---|
111 | integer :: itime,itimec,nstop,ldt,i,j,k,nrand,loop,memindnext,mind |
---|
112 | integer :: ngr,nix,njy,ks,nsp,nrelpoint |
---|
113 | real :: dz,dz1,dz2,xlon,ylat,xpol,ypol,gridsize |
---|
114 | real :: ru,rv,rw,dt,ux,vy,cosfact,xtn,ytn,tropop |
---|
115 | real :: prob(maxspec),up,vp,wp,dxsave,dysave,dawsave |
---|
116 | real :: dcwsave |
---|
117 | real :: usigold,vsigold,wsigold,r,rs |
---|
118 | real :: uold,vold,wold,vdepo(maxspec) |
---|
119 | real :: h1(2) |
---|
120 | !real uprof(nzmax),vprof(nzmax),wprof(nzmax) |
---|
121 | !real usigprof(nzmax),vsigprof(nzmax),wsigprof(nzmax) |
---|
122 | !real rhoprof(nzmax),rhogradprof(nzmax) |
---|
123 | real :: rhoa,rhograd,delz,dtf,rhoaux,dtftlw,uxscale,wpscale |
---|
124 | integer(kind=2) :: icbt |
---|
125 | real,parameter :: eps=nxmax/3.e5,eps2=1.e-9,eps3=tiny(1.0) |
---|
126 | real :: ptot_lhh,Q_lhh,phi_lhh,ath,bth !modified by mc |
---|
127 | real :: old_wp_buf,dcas,dcas1,del_test !added by mc |
---|
128 | integer :: i_well,jj,flagrein !test well mixed: modified by mc |
---|
129 | |
---|
130 | |
---|
131 | !!! CHANGE: TEST OF THE WELL-MIXED CRITERION |
---|
132 | ! integer,parameter :: iclass=10 |
---|
133 | ! real(kind=dp) :: zacc,tacc,t(iclass),th(0:iclass),hsave |
---|
134 | ! logical dump |
---|
135 | ! save zacc,tacc,t,th,hsave,dump |
---|
136 | !!! CHANGE |
---|
137 | |
---|
138 | integer :: idummy = -7 |
---|
139 | real :: settling = 0. |
---|
140 | |
---|
141 | |
---|
142 | !!! CHANGE: TEST OF THE WELL-MIXED CRITERION |
---|
143 | !if (idummy.eq.-7) then |
---|
144 | !open(550,file='WELLMIXEDTEST') |
---|
145 | !do 17 i=0,iclass |
---|
146 | !7 th(i)=real(i)/real(iclass) |
---|
147 | !endif |
---|
148 | !!! CHANGE |
---|
149 | |
---|
150 | |
---|
151 | nstop=0 |
---|
152 | do i=1,nmixz |
---|
153 | indzindicator(i)=.true. |
---|
154 | end do |
---|
155 | |
---|
156 | |
---|
157 | if (DRYDEP) then ! reset probability for deposition |
---|
158 | do ks=1,nspec |
---|
159 | depoindicator(ks)=.true. |
---|
160 | prob(ks)=0. |
---|
161 | end do |
---|
162 | endif |
---|
163 | |
---|
164 | dxsave=0. ! reset position displacements |
---|
165 | dysave=0. ! due to mean wind |
---|
166 | dawsave=0. ! and turbulent wind |
---|
167 | dcwsave=0. |
---|
168 | |
---|
169 | itimec=itime |
---|
170 | |
---|
171 | nrand=int(ran3(idummy)*real(maxrand-1))+1 |
---|
172 | |
---|
173 | |
---|
174 | ! Determine whether lat/long grid or polarstereographic projection |
---|
175 | ! is to be used |
---|
176 | ! Furthermore, determine which nesting level to be used |
---|
177 | !***************************************************************** |
---|
178 | |
---|
179 | if (nglobal.and.(yt.gt.switchnorthg)) then |
---|
180 | ngrid=-1 |
---|
181 | else if (sglobal.and.(yt.lt.switchsouthg)) then |
---|
182 | ngrid=-2 |
---|
183 | else |
---|
184 | ngrid=0 |
---|
185 | do j=numbnests,1,-1 |
---|
186 | if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. & |
---|
187 | (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then |
---|
188 | ngrid=j |
---|
189 | goto 23 |
---|
190 | endif |
---|
191 | end do |
---|
192 | 23 continue |
---|
193 | endif |
---|
194 | |
---|
195 | |
---|
196 | !*************************** |
---|
197 | ! Interpolate necessary data |
---|
198 | !*************************** |
---|
199 | |
---|
200 | if (abs(itime-memtime(1)).lt.abs(itime-memtime(2))) then |
---|
201 | memindnext=1 |
---|
202 | else |
---|
203 | memindnext=2 |
---|
204 | endif |
---|
205 | |
---|
206 | ! Determine nested grid coordinates |
---|
207 | !********************************** |
---|
208 | |
---|
209 | if (ngrid.gt.0) then |
---|
210 | xtn=(xt-xln(ngrid))*xresoln(ngrid) |
---|
211 | ytn=(yt-yln(ngrid))*yresoln(ngrid) |
---|
212 | ix=int(xtn) |
---|
213 | jy=int(ytn) |
---|
214 | nix=nint(xtn) |
---|
215 | njy=nint(ytn) |
---|
216 | else |
---|
217 | ix=int(xt) |
---|
218 | jy=int(yt) |
---|
219 | nix=nint(xt) |
---|
220 | njy=nint(yt) |
---|
221 | endif |
---|
222 | ixp=ix+1 |
---|
223 | jyp=jy+1 |
---|
224 | |
---|
225 | |
---|
226 | ! Determine the lower left corner and its distance to the current position |
---|
227 | !************************************************************************* |
---|
228 | |
---|
229 | ddx=xt-real(ix) |
---|
230 | ddy=yt-real(jy) |
---|
231 | rddx=1.-ddx |
---|
232 | rddy=1.-ddy |
---|
233 | p1=rddx*rddy |
---|
234 | p2=ddx*rddy |
---|
235 | p3=rddx*ddy |
---|
236 | p4=ddx*ddy |
---|
237 | |
---|
238 | ! Calculate variables for time interpolation |
---|
239 | !******************************************* |
---|
240 | |
---|
241 | dt1=real(itime-memtime(1)) |
---|
242 | dt2=real(memtime(2)-itime) |
---|
243 | dtt=1./(dt1+dt2) |
---|
244 | |
---|
245 | ! eso: Temporary fix for particle exactly at north pole |
---|
246 | if (jyp >= nymax) then |
---|
247 | ! write(*,*) 'WARNING: advance.f90 jyp >= nymax. xt,yt:',xt,yt |
---|
248 | jyp=jyp-1 |
---|
249 | end if |
---|
250 | |
---|
251 | ! Compute maximum mixing height around particle position |
---|
252 | !******************************************************* |
---|
253 | |
---|
254 | h=0. |
---|
255 | if (ngrid.le.0) then |
---|
256 | do k=1,2 |
---|
257 | mind=memind(k) ! eso: compatibility with 3-field version |
---|
258 | if (interpolhmix) then |
---|
259 | h1(k)=p1*hmix(ix ,jy ,1,mind) & |
---|
260 | + p2*hmix(ixp,jy ,1,mind) & |
---|
261 | + p3*hmix(ix ,jyp,1,mind) & |
---|
262 | + p4*hmix(ixp,jyp,1,mind) |
---|
263 | else |
---|
264 | do j=jy,jyp |
---|
265 | do i=ix,ixp |
---|
266 | if (hmix(i,j,1,mind).gt.h) h=hmix(i,j,1,mind) |
---|
267 | end do |
---|
268 | end do |
---|
269 | endif |
---|
270 | end do |
---|
271 | tropop=tropopause(nix,njy,1,1) |
---|
272 | else |
---|
273 | do k=1,2 |
---|
274 | mind=memind(k) |
---|
275 | do j=jy,jyp |
---|
276 | do i=ix,ixp |
---|
277 | if (hmixn(i,j,1,mind,ngrid).gt.h) h=hmixn(i,j,1,mind,ngrid) |
---|
278 | end do |
---|
279 | end do |
---|
280 | end do |
---|
281 | tropop=tropopausen(nix,njy,1,1,ngrid) |
---|
282 | endif |
---|
283 | |
---|
284 | if (interpolhmix) h=(h1(1)*dt2+h1(2)*dt1)*dtt |
---|
285 | zeta=zt/h |
---|
286 | |
---|
287 | |
---|
288 | |
---|
289 | !************************************************************* |
---|
290 | ! If particle is in the PBL, interpolate once and then make a |
---|
291 | ! time loop until end of interval is reached |
---|
292 | !************************************************************* |
---|
293 | |
---|
294 | if (zeta.le.1.) then |
---|
295 | |
---|
296 | ! BEGIN TIME LOOP |
---|
297 | !================ |
---|
298 | |
---|
299 | loop=0 |
---|
300 | 100 loop=loop+1 |
---|
301 | if (method.eq.1) then |
---|
302 | ldt=min(ldt,abs(lsynctime-itimec+itime)) |
---|
303 | itimec=itimec+ldt*ldirect |
---|
304 | else |
---|
305 | ldt=abs(lsynctime) |
---|
306 | itimec=itime+lsynctime |
---|
307 | endif |
---|
308 | dt=real(ldt) |
---|
309 | |
---|
310 | zeta=zt/h |
---|
311 | |
---|
312 | |
---|
313 | if (loop.eq.1) then |
---|
314 | if (ngrid.le.0) then |
---|
315 | xts=real(xt) |
---|
316 | yts=real(yt) |
---|
317 | call interpol_all(itime,xts,yts,zt) |
---|
318 | else |
---|
319 | call interpol_all_nests(itime,xtn,ytn,zt) |
---|
320 | endif |
---|
321 | |
---|
322 | else |
---|
323 | |
---|
324 | |
---|
325 | ! Determine the level below the current position for u,v,rho |
---|
326 | !*********************************************************** |
---|
327 | |
---|
328 | do i=2,nz |
---|
329 | if (height(i).gt.zt) then |
---|
330 | indz=i-1 |
---|
331 | indzp=i |
---|
332 | goto 6 |
---|
333 | endif |
---|
334 | end do |
---|
335 | 6 continue |
---|
336 | |
---|
337 | ! If one of the levels necessary is not yet available, |
---|
338 | ! calculate it |
---|
339 | !***************************************************** |
---|
340 | |
---|
341 | do i=indz,indzp |
---|
342 | if (indzindicator(i)) then |
---|
343 | if (ngrid.le.0) then |
---|
344 | call interpol_misslev(i) |
---|
345 | else |
---|
346 | call interpol_misslev_nests(i) |
---|
347 | endif |
---|
348 | endif |
---|
349 | end do |
---|
350 | endif |
---|
351 | |
---|
352 | |
---|
353 | ! Vertical interpolation of u,v,w,rho and drhodz |
---|
354 | !*********************************************** |
---|
355 | |
---|
356 | ! Vertical distance to the level below and above current position |
---|
357 | ! both in terms of (u,v) and (w) fields |
---|
358 | !**************************************************************** |
---|
359 | |
---|
360 | dz=1./(height(indzp)-height(indz)) |
---|
361 | dz1=(zt-height(indz))*dz |
---|
362 | dz2=(height(indzp)-zt)*dz |
---|
363 | |
---|
364 | u=dz1*uprof(indzp)+dz2*uprof(indz) |
---|
365 | v=dz1*vprof(indzp)+dz2*vprof(indz) |
---|
366 | w=dz1*wprof(indzp)+dz2*wprof(indz) |
---|
367 | rhoa=dz1*rhoprof(indzp)+dz2*rhoprof(indz) |
---|
368 | rhograd=dz1*rhogradprof(indzp)+dz2*rhogradprof(indz) |
---|
369 | |
---|
370 | |
---|
371 | ! Compute the turbulent disturbances |
---|
372 | ! Determine the sigmas and the timescales |
---|
373 | !**************************************** |
---|
374 | |
---|
375 | if (turbswitch) then |
---|
376 | call hanna(zt) |
---|
377 | else |
---|
378 | call hanna1(zt) |
---|
379 | endif |
---|
380 | |
---|
381 | |
---|
382 | !***************************************** |
---|
383 | ! Determine the new diffusivity velocities |
---|
384 | !***************************************** |
---|
385 | |
---|
386 | ! Horizontal components |
---|
387 | !********************** |
---|
388 | |
---|
389 | if (nrand+1.gt.maxrand) nrand=1 |
---|
390 | if (dt/tlu.lt..5) then |
---|
391 | up=(1.-dt/tlu)*up+rannumb(nrand)*sigu*sqrt(2.*dt/tlu) |
---|
392 | else |
---|
393 | ru=exp(-dt/tlu) |
---|
394 | up=ru*up+rannumb(nrand)*sigu*sqrt(1.-ru**2) |
---|
395 | endif |
---|
396 | if (dt/tlv.lt..5) then |
---|
397 | vp=(1.-dt/tlv)*vp+rannumb(nrand+1)*sigv*sqrt(2.*dt/tlv) |
---|
398 | else |
---|
399 | rv=exp(-dt/tlv) |
---|
400 | vp=rv*vp+rannumb(nrand+1)*sigv*sqrt(1.-rv**2) |
---|
401 | endif |
---|
402 | nrand=nrand+2 |
---|
403 | |
---|
404 | |
---|
405 | if (nrand+ifine.gt.maxrand) nrand=1 |
---|
406 | rhoaux=rhograd/rhoa |
---|
407 | dtf=dt*fine |
---|
408 | |
---|
409 | dtftlw=dtf/tlw |
---|
410 | |
---|
411 | ! Loop over ifine short time steps for vertical component |
---|
412 | !******************************************************** |
---|
413 | |
---|
414 | do i=1,ifine |
---|
415 | |
---|
416 | ! Determine the drift velocity and density correction velocity |
---|
417 | !************************************************************* |
---|
418 | |
---|
419 | if (turbswitch) then |
---|
420 | if (dtftlw.lt..5) then |
---|
421 | !************************************************************* |
---|
422 | !************** CBL options added by mc see routine cblf90 *** |
---|
423 | if (cblflag.eq.1) then !modified by mc |
---|
424 | if (-h/ol.gt.5) then !modified by mc |
---|
425 | !if (ol.lt.0.) then !modified by mc |
---|
426 | !if (ol.gt.0.) then !modified by mc : for test |
---|
427 | !print *,zt,wp,ath,bth,tlw,dtf,'prima' |
---|
428 | flagrein=0 |
---|
429 | nrand=nrand+1 |
---|
430 | old_wp_buf=wp |
---|
431 | call cbl(wp,zt,ust,wst,h,rhoa,rhograd,sigw,dsigwdz,tlw,ptot_lhh,Q_lhh,phi_lhh,ath,bth,ol,flagrein) !inside the routine for inverse time |
---|
432 | wp=(wp+ath*dtf+bth*rannumb(nrand)*sqrt(dtf))*real(icbt) |
---|
433 | ! wp=(wp+ath*dtf+bth*gasdev2(mydum)*sqrt(dtf))*real(icbt) |
---|
434 | delz=wp*dtf |
---|
435 | if (flagrein.eq.1) then |
---|
436 | call re_initialize_particle(zt,ust,wst,h,sigw,old_wp_buf,nrand,ol) |
---|
437 | wp=old_wp_buf |
---|
438 | delz=wp*dtf |
---|
439 | nan_count=nan_count+1 |
---|
440 | end if |
---|
441 | !print *,zt,wp,ath,bth,tlw,dtf,rannumb(nrand+i),icbt |
---|
442 | !pause |
---|
443 | else |
---|
444 | nrand=nrand+1 |
---|
445 | old_wp_buf=wp |
---|
446 | ath=-wp/tlw+sigw*dsigwdz+wp*wp/sigw*dsigwdz+sigw*sigw/rhoa*rhograd !1-note for inverse time should be -wp/tlw*ldirect+... calculated for wp=-wp |
---|
447 | !2-but since ldirect =-1 for inverse time and this must be calculated for (-wp) and |
---|
448 | !3-the gaussian pdf is symmetric (i.e. pdf(w)=pdf(-w) ldirect can be discarded |
---|
449 | bth=sigw*rannumb(nrand)*sqrt(2.*dtftlw) |
---|
450 | wp=(wp+ath*dtf+bth)*real(icbt) |
---|
451 | delz=wp*dtf |
---|
452 | del_test=(1.-wp)/wp !catch infinity value |
---|
453 | if (isnan(wp).or.isnan(del_test)) then |
---|
454 | nrand=nrand+1 |
---|
455 | wp=sigw*rannumb(nrand) |
---|
456 | delz=wp*dtf |
---|
457 | nan_count2=nan_count2+1 |
---|
458 | !print *,'NaN coutner equal to:', nan_count,'reduce ifine if this number became a non-negligible fraction of the particle number' |
---|
459 | end if |
---|
460 | end if |
---|
461 | !******************** END CBL option ******************************* |
---|
462 | !******************************************************************* |
---|
463 | else |
---|
464 | wp=((1.-dtftlw)*wp+rannumb(nrand+i)*sqrt(2.*dtftlw) & |
---|
465 | +dtf*(dsigwdz+rhoaux*sigw))*real(icbt) |
---|
466 | delz=wp*sigw*dtf |
---|
467 | end if |
---|
468 | else |
---|
469 | rw=exp(-dtftlw) |
---|
470 | wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2) & |
---|
471 | +tlw*(1.-rw)*(dsigwdz+rhoaux*sigw))*real(icbt) |
---|
472 | delz=wp*sigw*dtf |
---|
473 | endif |
---|
474 | |
---|
475 | else |
---|
476 | rw=exp(-dtftlw) |
---|
477 | wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2)*sigw & |
---|
478 | +tlw*(1.-rw)*(dsigw2dz+rhoaux*sigw**2))*real(icbt) |
---|
479 | delz=wp*dtf |
---|
480 | endif |
---|
481 | |
---|
482 | if (turboff) then |
---|
483 | !sec switch off turbulence |
---|
484 | up=0.0 |
---|
485 | vp=0.0 |
---|
486 | wp=0.0 |
---|
487 | delz=0. |
---|
488 | endif |
---|
489 | |
---|
490 | !**************************************************** |
---|
491 | ! Compute turbulent vertical displacement of particle |
---|
492 | !**************************************************** |
---|
493 | |
---|
494 | if (abs(delz).gt.h) delz=mod(delz,h) |
---|
495 | |
---|
496 | ! Determine if particle transfers to a "forbidden state" below the ground |
---|
497 | ! or above the mixing height |
---|
498 | !************************************************************************ |
---|
499 | |
---|
500 | if (delz.lt.-zt) then ! reflection at ground |
---|
501 | icbt=-1 |
---|
502 | zt=-zt-delz |
---|
503 | else if (delz.gt.(h-zt)) then ! reflection at h |
---|
504 | icbt=-1 |
---|
505 | zt=-zt-delz+2.*h |
---|
506 | else ! no reflection |
---|
507 | icbt=1 |
---|
508 | zt=zt+delz |
---|
509 | endif |
---|
510 | |
---|
511 | if (i.ne.ifine) then |
---|
512 | zeta=zt/h |
---|
513 | call hanna_short(zt) |
---|
514 | endif |
---|
515 | |
---|
516 | end do |
---|
517 | if (cblflag.ne.1) nrand=nrand+i |
---|
518 | |
---|
519 | ! Determine time step for next integration |
---|
520 | !***************************************** |
---|
521 | |
---|
522 | if (turbswitch) then |
---|
523 | ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), & |
---|
524 | 0.5/abs(dsigwdz))*ctl) |
---|
525 | else |
---|
526 | ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5))*ctl) |
---|
527 | endif |
---|
528 | ldt=max(ldt,mintime) |
---|
529 | |
---|
530 | |
---|
531 | ! If particle represents only a single species, add gravitational settling |
---|
532 | ! velocity. The settling velocity is zero for gases, or if particle |
---|
533 | ! represents more than one species |
---|
534 | !************************************************************************* |
---|
535 | |
---|
536 | if (mdomainfill.eq.0) then |
---|
537 | if (lsettling) then |
---|
538 | do nsp=1,nspec |
---|
539 | if (xmass(nrelpoint,nsp).gt.eps3) exit |
---|
540 | end do |
---|
541 | if (nsp.gt.nspec) then |
---|
542 | nsp=nspec |
---|
543 | end if |
---|
544 | if (density(nsp).gt.0.) then |
---|
545 | call get_settling(itime,real(xt),real(yt),zt,nsp,settling) !bugfix |
---|
546 | w=w+settling |
---|
547 | end if |
---|
548 | end if |
---|
549 | endif |
---|
550 | |
---|
551 | ! Horizontal displacements during time step dt are small real values compared |
---|
552 | ! to the position; adding the two, would result in large numerical errors. |
---|
553 | ! Thus, displacements are accumulated during lsynctime and are added to the |
---|
554 | ! position at the end |
---|
555 | !**************************************************************************** |
---|
556 | |
---|
557 | dxsave=dxsave+u*dt |
---|
558 | dysave=dysave+v*dt |
---|
559 | dawsave=dawsave+up*dt |
---|
560 | dcwsave=dcwsave+vp*dt |
---|
561 | zt=zt+w*dt*real(ldirect) |
---|
562 | |
---|
563 | ! HSO/AL: Particle managed to go over highest level -> interpolation error in goto 700 |
---|
564 | ! alias interpol_wind (division by zero) |
---|
565 | if (zt.ge.height(nz)) zt=height(nz)-100.*eps |
---|
566 | |
---|
567 | if (zt.gt.h) then |
---|
568 | if (itimec.eq.itime+lsynctime) goto 99 |
---|
569 | goto 700 ! complete the current interval above PBL |
---|
570 | endif |
---|
571 | |
---|
572 | |
---|
573 | !!! CHANGE: TEST OF THE WELL-MIXED CRITERION |
---|
574 | !!! These lines may be switched on to test the well-mixed criterion |
---|
575 | !if (zt.le.h) then |
---|
576 | ! zacc=zacc+zt/h*dt |
---|
577 | ! hsave=hsave+h*dt |
---|
578 | ! tacc=tacc+dt |
---|
579 | ! do 67 i=1,iclass |
---|
580 | ! if ((zt/h.gt.th(i-1)).and.(zt/h.le.th(i))) |
---|
581 | ! + t(i)=t(i)+dt |
---|
582 | !7 continue |
---|
583 | !endif |
---|
584 | !if ((mod(itime,10800).eq.0).and.dump) then |
---|
585 | ! dump=.false. |
---|
586 | ! write(550,'(i6,12f10.3)') itime,hsave/tacc,zacc/tacc, |
---|
587 | ! + (t(i)/tacc*real(iclass),i=1,iclass) |
---|
588 | ! zacc=0. |
---|
589 | ! tacc=0. |
---|
590 | ! do 68 i=1,iclass |
---|
591 | !8 t(i)=0. |
---|
592 | ! hsave=0. |
---|
593 | !endif |
---|
594 | !if (mod(itime,10800).ne.0) dump=.true. |
---|
595 | !!! CHANGE |
---|
596 | |
---|
597 | ! Determine probability of deposition |
---|
598 | !************************************ |
---|
599 | |
---|
600 | if ((DRYDEP).and.(zt.lt.2.*href)) then |
---|
601 | do ks=1,nspec |
---|
602 | if (DRYDEPSPEC(ks)) then |
---|
603 | if (depoindicator(ks)) then |
---|
604 | if (ngrid.le.0) then |
---|
605 | call interpol_vdep(ks,vdepo(ks)) |
---|
606 | else |
---|
607 | call interpol_vdep_nests(ks,vdepo(ks)) |
---|
608 | endif |
---|
609 | endif |
---|
610 | ! correction by Petra Seibert, 10 April 2001 |
---|
611 | ! this formulation means that prob(n) = 1 - f(0)*...*f(n) |
---|
612 | ! where f(n) is the exponential term |
---|
613 | prob(ks)=1.+(prob(ks)-1.)* & |
---|
614 | exp(-vdepo(ks)*abs(dt)/(2.*href)) |
---|
615 | endif |
---|
616 | end do |
---|
617 | endif |
---|
618 | |
---|
619 | if (zt.lt.0.) zt=min(h-eps2,-1.*zt) ! if particle below ground -> reflection |
---|
620 | |
---|
621 | if (itimec.eq.(itime+lsynctime)) then |
---|
622 | usig=0.5*(usigprof(indzp)+usigprof(indz)) |
---|
623 | vsig=0.5*(vsigprof(indzp)+vsigprof(indz)) |
---|
624 | wsig=0.5*(wsigprof(indzp)+wsigprof(indz)) |
---|
625 | goto 99 ! finished |
---|
626 | endif |
---|
627 | goto 100 |
---|
628 | |
---|
629 | ! END TIME LOOP |
---|
630 | !============== |
---|
631 | |
---|
632 | |
---|
633 | endif |
---|
634 | |
---|
635 | |
---|
636 | |
---|
637 | !********************************************************** |
---|
638 | ! For all particles that are outside the PBL, make a single |
---|
639 | ! time step. Only horizontal turbulent disturbances are |
---|
640 | ! calculated. Vertical disturbances are reset. |
---|
641 | !********************************************************** |
---|
642 | |
---|
643 | |
---|
644 | ! Interpolate the wind |
---|
645 | !********************* |
---|
646 | |
---|
647 | 700 continue |
---|
648 | if (ngrid.le.0) then |
---|
649 | xts=real(xt) |
---|
650 | yts=real(yt) |
---|
651 | call interpol_wind(itime,xts,yts,zt) |
---|
652 | else |
---|
653 | call interpol_wind_nests(itime,xtn,ytn,zt) |
---|
654 | endif |
---|
655 | |
---|
656 | |
---|
657 | ! Compute everything for above the PBL |
---|
658 | |
---|
659 | ! Assume constant, uncorrelated, turbulent perturbations |
---|
660 | ! In the stratosphere, use a small vertical diffusivity d_strat, |
---|
661 | ! in the troposphere, use a larger horizontal diffusivity d_trop. |
---|
662 | ! Turbulent velocity scales are determined based on sqrt(d_trop/dt) |
---|
663 | !****************************************************************** |
---|
664 | |
---|
665 | ldt=abs(lsynctime-itimec+itime) |
---|
666 | dt=real(ldt) |
---|
667 | |
---|
668 | if (zt.lt.tropop) then ! in the troposphere |
---|
669 | uxscale=sqrt(2.*d_trop/dt) |
---|
670 | if (nrand+1.gt.maxrand) nrand=1 |
---|
671 | ux=rannumb(nrand)*uxscale |
---|
672 | vy=rannumb(nrand+1)*uxscale |
---|
673 | nrand=nrand+2 |
---|
674 | wp=0. |
---|
675 | else if (zt.lt.tropop+1000.) then ! just above the tropopause: make transition |
---|
676 | weight=(zt-tropop)/1000. |
---|
677 | uxscale=sqrt(2.*d_trop/dt*(1.-weight)) |
---|
678 | if (nrand+2.gt.maxrand) nrand=1 |
---|
679 | ux=rannumb(nrand)*uxscale |
---|
680 | vy=rannumb(nrand+1)*uxscale |
---|
681 | wpscale=sqrt(2.*d_strat/dt*weight) |
---|
682 | wp=rannumb(nrand+2)*wpscale+d_strat/1000. |
---|
683 | nrand=nrand+3 |
---|
684 | else ! in the stratosphere |
---|
685 | if (nrand.gt.maxrand) nrand=1 |
---|
686 | ux=0. |
---|
687 | vy=0. |
---|
688 | wpscale=sqrt(2.*d_strat/dt) |
---|
689 | wp=rannumb(nrand)*wpscale |
---|
690 | nrand=nrand+1 |
---|
691 | endif |
---|
692 | |
---|
693 | if (turboff) then |
---|
694 | !sec switch off turbulence |
---|
695 | ux=0.0 |
---|
696 | vy=0.0 |
---|
697 | wp=0.0 |
---|
698 | endif |
---|
699 | |
---|
700 | ! If particle represents only a single species, add gravitational settling |
---|
701 | ! velocity. The settling velocity is zero for gases |
---|
702 | !************************************************************************* |
---|
703 | |
---|
704 | if (mdomainfill.eq.0) then |
---|
705 | if (lsettling) then |
---|
706 | do nsp=1,nspec |
---|
707 | if (xmass(nrelpoint,nsp).gt.eps3) exit |
---|
708 | end do |
---|
709 | if (nsp.gt.nspec) then |
---|
710 | nsp=nspec |
---|
711 | end if |
---|
712 | if (density(nsp).gt.0.) then |
---|
713 | call get_settling(itime,real(xt),real(yt),zt,nsp,settling) !bugfix |
---|
714 | w=w+settling |
---|
715 | end if |
---|
716 | endif |
---|
717 | end if |
---|
718 | |
---|
719 | |
---|
720 | ! Calculate position at time step itime+lsynctime |
---|
721 | !************************************************ |
---|
722 | |
---|
723 | dxsave=dxsave+(u+ux)*dt |
---|
724 | dysave=dysave+(v+vy)*dt |
---|
725 | zt=zt+(w+wp)*dt*real(ldirect) |
---|
726 | if (zt.lt.0.) zt=min(h-eps2,-1.*zt) ! if particle below ground -> reflection |
---|
727 | |
---|
728 | 99 continue |
---|
729 | |
---|
730 | |
---|
731 | |
---|
732 | !**************************************************************** |
---|
733 | ! Add mesoscale random disturbances |
---|
734 | ! This is done only once for the whole lsynctime interval to save |
---|
735 | ! computation time |
---|
736 | !**************************************************************** |
---|
737 | |
---|
738 | |
---|
739 | ! Mesoscale wind velocity fluctuations are obtained by scaling |
---|
740 | ! with the standard deviation of the grid-scale winds surrounding |
---|
741 | ! the particle location, multiplied by a factor turbmesoscale. |
---|
742 | ! The autocorrelation time constant is taken as half the |
---|
743 | ! time interval between wind fields |
---|
744 | !**************************************************************** |
---|
745 | |
---|
746 | r=exp(-2.*real(abs(lsynctime))/real(lwindinterv)) |
---|
747 | rs=sqrt(1.-r**2) |
---|
748 | if (nrand+2.gt.maxrand) nrand=1 |
---|
749 | usigold=r*usigold+rs*rannumb(nrand)*usig*turbmesoscale |
---|
750 | vsigold=r*vsigold+rs*rannumb(nrand+1)*vsig*turbmesoscale |
---|
751 | wsigold=r*wsigold+rs*rannumb(nrand+2)*wsig*turbmesoscale |
---|
752 | |
---|
753 | dxsave=dxsave+usigold*real(lsynctime) |
---|
754 | dysave=dysave+vsigold*real(lsynctime) |
---|
755 | |
---|
756 | zt=zt+wsigold*real(lsynctime) |
---|
757 | if (zt.lt.0.) zt=-1.*zt ! if particle below ground -> refletion |
---|
758 | |
---|
759 | !************************************************************* |
---|
760 | ! Transform along and cross wind components to xy coordinates, |
---|
761 | ! add them to u and v, transform u,v to grid units/second |
---|
762 | ! and calculate new position |
---|
763 | !************************************************************* |
---|
764 | |
---|
765 | call windalign(dxsave,dysave,dawsave,dcwsave,ux,vy) |
---|
766 | dxsave=dxsave+ux ! comment by mc: comment this line to stop the particles horizontally for test reasons |
---|
767 | dysave=dysave+vy |
---|
768 | if (ngrid.ge.0) then |
---|
769 | cosfact=dxconst/cos((yt*dy+ylat0)*pi180) |
---|
770 | xt=xt+real(dxsave*cosfact*real(ldirect),kind=dp) |
---|
771 | yt=yt+real(dysave*dyconst*real(ldirect),kind=dp) |
---|
772 | else if (ngrid.eq.-1) then ! around north pole |
---|
773 | xlon=xlon0+xt*dx !comment by mc: compute old particle position |
---|
774 | ylat=ylat0+yt*dy |
---|
775 | call cll2xy(northpolemap,ylat,xlon,xpol,ypol) !convert old particle position in polar stereographic |
---|
776 | gridsize=1000.*cgszll(northpolemap,ylat,xlon) !calculate size in m of grid element in polar stereographic coordinate |
---|
777 | dxsave=dxsave/gridsize !increment from meter to grdi unit |
---|
778 | dysave=dysave/gridsize |
---|
779 | xpol=xpol+dxsave*real(ldirect) !position in grid unit polar stereographic |
---|
780 | ypol=ypol+dysave*real(ldirect) |
---|
781 | call cxy2ll(northpolemap,xpol,ypol,ylat,xlon) !convert to lat long coordinate |
---|
782 | xt=(xlon-xlon0)/dx !convert to grid units in lat long coordinate, comment by mc |
---|
783 | yt=(ylat-ylat0)/dy |
---|
784 | else if (ngrid.eq.-2) then ! around south pole |
---|
785 | xlon=xlon0+xt*dx |
---|
786 | ylat=ylat0+yt*dy |
---|
787 | call cll2xy(southpolemap,ylat,xlon,xpol,ypol) |
---|
788 | gridsize=1000.*cgszll(southpolemap,ylat,xlon) |
---|
789 | dxsave=dxsave/gridsize |
---|
790 | dysave=dysave/gridsize |
---|
791 | xpol=xpol+dxsave*real(ldirect) |
---|
792 | ypol=ypol+dysave*real(ldirect) |
---|
793 | call cxy2ll(southpolemap,xpol,ypol,ylat,xlon) |
---|
794 | xt=(xlon-xlon0)/dx |
---|
795 | yt=(ylat-ylat0)/dy |
---|
796 | endif |
---|
797 | |
---|
798 | |
---|
799 | ! If global data are available, use cyclic boundary condition |
---|
800 | !************************************************************ |
---|
801 | |
---|
802 | if (xglobal) then |
---|
803 | if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1) |
---|
804 | if (xt.lt.0.) xt=xt+real(nxmin1) |
---|
805 | if (xt.le.eps) xt=eps |
---|
806 | if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps |
---|
807 | endif |
---|
808 | |
---|
809 | ! HSO/AL: Prevent particles from disappearing at the pole |
---|
810 | !****************************************************************** |
---|
811 | |
---|
812 | if ( yt.lt.0. ) then |
---|
813 | xt=mod(xt+180.,360.) |
---|
814 | yt=-yt |
---|
815 | else if ( yt.gt.real(nymin1) ) then |
---|
816 | xt=mod(xt+180.,360.) |
---|
817 | yt=2*real(nymin1)-yt |
---|
818 | endif |
---|
819 | |
---|
820 | ! Check position: If trajectory outside model domain, terminate it |
---|
821 | !***************************************************************** |
---|
822 | |
---|
823 | if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. & |
---|
824 | (yt.gt.real(nymin1))) then |
---|
825 | nstop=3 |
---|
826 | return |
---|
827 | endif |
---|
828 | |
---|
829 | ! If particle above highest model level, set it back into the domain |
---|
830 | !******************************************************************* |
---|
831 | |
---|
832 | if (zt.ge.height(nz)) zt=height(nz)-100.*eps |
---|
833 | |
---|
834 | |
---|
835 | !************************************************************************ |
---|
836 | ! Now we could finish, as this was done in FLEXPART versions up to 4.0. |
---|
837 | ! However, truncation errors of the advection can be significantly |
---|
838 | ! reduced by doing one iteration of the Petterssen scheme, if this is |
---|
839 | ! possible. |
---|
840 | ! Note that this is applied only to the grid-scale winds, not to |
---|
841 | ! the turbulent winds. |
---|
842 | !************************************************************************ |
---|
843 | |
---|
844 | ! The Petterssen scheme can only applied with long time steps (only then u |
---|
845 | ! is the "old" wind as required by the scheme); otherwise do nothing |
---|
846 | !************************************************************************* |
---|
847 | |
---|
848 | if (ldt.ne.abs(lsynctime)) return |
---|
849 | |
---|
850 | ! The Petterssen scheme can only be applied if the ending time of the time step |
---|
851 | ! (itime+ldt*ldirect) is still between the two wind fields held in memory; |
---|
852 | ! otherwise do nothing |
---|
853 | !****************************************************************************** |
---|
854 | |
---|
855 | if (abs(itime+ldt*ldirect).gt.abs(memtime(2))) return |
---|
856 | |
---|
857 | ! Apply it also only if starting and ending point of current time step are on |
---|
858 | ! the same grid; otherwise do nothing |
---|
859 | !***************************************************************************** |
---|
860 | if (nglobal.and.(yt.gt.switchnorthg)) then |
---|
861 | ngr=-1 |
---|
862 | else if (sglobal.and.(yt.lt.switchsouthg)) then |
---|
863 | ngr=-2 |
---|
864 | else |
---|
865 | ngr=0 |
---|
866 | do j=numbnests,1,-1 |
---|
867 | if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. & |
---|
868 | (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then |
---|
869 | ngr=j |
---|
870 | goto 43 |
---|
871 | endif |
---|
872 | end do |
---|
873 | 43 continue |
---|
874 | endif |
---|
875 | |
---|
876 | if (ngr.ne.ngrid) return |
---|
877 | |
---|
878 | ! Determine nested grid coordinates |
---|
879 | !********************************** |
---|
880 | |
---|
881 | if (ngrid.gt.0) then |
---|
882 | xtn=(xt-xln(ngrid))*xresoln(ngrid) |
---|
883 | ytn=(yt-yln(ngrid))*yresoln(ngrid) |
---|
884 | ix=int(xtn) |
---|
885 | jy=int(ytn) |
---|
886 | else |
---|
887 | ix=int(xt) |
---|
888 | jy=int(yt) |
---|
889 | endif |
---|
890 | ixp=ix+1 |
---|
891 | jyp=jy+1 |
---|
892 | |
---|
893 | |
---|
894 | ! Memorize the old wind |
---|
895 | !********************** |
---|
896 | |
---|
897 | uold=u |
---|
898 | vold=v |
---|
899 | wold=w |
---|
900 | |
---|
901 | ! Interpolate wind at new position and time |
---|
902 | !****************************************** |
---|
903 | |
---|
904 | if (ngrid.le.0) then |
---|
905 | xts=real(xt) |
---|
906 | yts=real(yt) |
---|
907 | call interpol_wind_short(itime+ldt*ldirect,xts,yts,zt) |
---|
908 | else |
---|
909 | call interpol_wind_short_nests(itime+ldt*ldirect,xtn,ytn,zt) |
---|
910 | endif |
---|
911 | |
---|
912 | if (mdomainfill.eq.0) then |
---|
913 | if (lsettling) then |
---|
914 | do nsp=1,nspec |
---|
915 | if (xmass(nrelpoint,nsp).gt.eps3) exit |
---|
916 | end do |
---|
917 | if (nsp.gt.nspec) then |
---|
918 | nsp=nspec |
---|
919 | end if |
---|
920 | if (density(nsp).gt.0.) then |
---|
921 | call get_settling(itime+ldt,real(xt),real(yt),zt,nsp,settling) !bugfix |
---|
922 | w=w+settling |
---|
923 | end if |
---|
924 | endif |
---|
925 | end if |
---|
926 | |
---|
927 | |
---|
928 | ! Determine the difference vector between new and old wind |
---|
929 | ! (use half of it to correct position according to Petterssen) |
---|
930 | !************************************************************* |
---|
931 | |
---|
932 | u=(u-uold)/2. |
---|
933 | v=(v-vold)/2. |
---|
934 | w=(w-wold)/2. |
---|
935 | |
---|
936 | |
---|
937 | ! Finally, correct the old position |
---|
938 | !********************************** |
---|
939 | |
---|
940 | zt=zt+w*real(ldt*ldirect) |
---|
941 | if (zt.lt.0.) zt=min(h-eps2,-1.*zt) ! if particle below ground -> reflection |
---|
942 | if (ngrid.ge.0) then |
---|
943 | cosfact=dxconst/cos((yt*dy+ylat0)*pi180) |
---|
944 | xt=xt+real(u*cosfact*real(ldt*ldirect),kind=dp) |
---|
945 | yt=yt+real(v*dyconst*real(ldt*ldirect),kind=dp) |
---|
946 | else if (ngrid.eq.-1) then ! around north pole |
---|
947 | xlon=xlon0+xt*dx |
---|
948 | ylat=ylat0+yt*dy |
---|
949 | call cll2xy(northpolemap,ylat,xlon,xpol,ypol) |
---|
950 | gridsize=1000.*cgszll(northpolemap,ylat,xlon) |
---|
951 | u=u/gridsize |
---|
952 | v=v/gridsize |
---|
953 | xpol=xpol+u*real(ldt*ldirect) |
---|
954 | ypol=ypol+v*real(ldt*ldirect) |
---|
955 | call cxy2ll(northpolemap,xpol,ypol,ylat,xlon) |
---|
956 | xt=(xlon-xlon0)/dx |
---|
957 | yt=(ylat-ylat0)/dy |
---|
958 | else if (ngrid.eq.-2) then ! around south pole |
---|
959 | xlon=xlon0+xt*dx |
---|
960 | ylat=ylat0+yt*dy |
---|
961 | call cll2xy(southpolemap,ylat,xlon,xpol,ypol) |
---|
962 | gridsize=1000.*cgszll(southpolemap,ylat,xlon) |
---|
963 | u=u/gridsize |
---|
964 | v=v/gridsize |
---|
965 | xpol=xpol+u*real(ldt*ldirect) |
---|
966 | ypol=ypol+v*real(ldt*ldirect) |
---|
967 | call cxy2ll(southpolemap,xpol,ypol,ylat,xlon) |
---|
968 | xt=(xlon-xlon0)/dx |
---|
969 | yt=(ylat-ylat0)/dy |
---|
970 | endif |
---|
971 | |
---|
972 | ! If global data are available, use cyclic boundary condition |
---|
973 | !************************************************************ |
---|
974 | |
---|
975 | if (xglobal) then |
---|
976 | if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1) |
---|
977 | if (xt.lt.0.) xt=xt+real(nxmin1) |
---|
978 | if (xt.le.eps) xt=eps |
---|
979 | if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps |
---|
980 | endif |
---|
981 | |
---|
982 | ! HSO/AL: Prevent particles from disappearing at the pole |
---|
983 | !****************************************************************** |
---|
984 | |
---|
985 | if ( yt.lt.0. ) then |
---|
986 | xt=mod(xt+180.,360.) |
---|
987 | yt=-yt |
---|
988 | else if ( yt.gt.real(nymin1) ) then |
---|
989 | xt=mod(xt+180.,360.) |
---|
990 | yt=2*real(nymin1)-yt |
---|
991 | endif |
---|
992 | |
---|
993 | ! Check position: If trajectory outside model domain, terminate it |
---|
994 | !***************************************************************** |
---|
995 | |
---|
996 | if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. & |
---|
997 | (yt.gt.real(nymin1))) then |
---|
998 | nstop=3 |
---|
999 | return |
---|
1000 | endif |
---|
1001 | |
---|
1002 | ! If particle above highest model level, set it back into the domain |
---|
1003 | !******************************************************************* |
---|
1004 | |
---|
1005 | if (zt.ge.height(nz)) zt=height(nz)-100.*eps |
---|
1006 | |
---|
1007 | |
---|
1008 | end subroutine advance |
---|
1009 | |
---|