source: flexpart.git/src/advance.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

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