source: flexpart.git/src/advance.f90

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

add SPDX-License-Identifier to all .f90 files

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