source: flexpart.git/src/advance.f90 @ 027e844

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 027e844 was 5184a7c, checked in by Espen Sollum ATMOS <eso@…>, 7 years ago

Enable settling with multiple species if from separate releases

  • Property mode set to 100644
File size: 34.6 KB
Line 
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
22subroutine 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
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
19223   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
300100   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
3356       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.eps2) exit
540          end do
541          if (nsp.gt.nspec) then
542  ! This should never happen         
543            write(*,*) 'advance.f90: ERROR: could not find releasepoint'
544            stop
545          end if
546          if (density(nsp).gt.0.) then
547            call get_settling(itime,real(xt),real(yt),zt,nsp,settling)  !bugfix
548            w=w+settling
549          end if
550        end if
551      endif
552
553  ! Horizontal displacements during time step dt are small real values compared
554  ! to the position; adding the two, would result in large numerical errors.
555  ! Thus, displacements are accumulated during lsynctime and are added to the
556  ! position at the end
557  !****************************************************************************
558
559      dxsave=dxsave+u*dt
560      dysave=dysave+v*dt
561      dawsave=dawsave+up*dt
562      dcwsave=dcwsave+vp*dt
563      zt=zt+w*dt*real(ldirect)
564
565      ! HSO/AL: Particle managed to go over highest level -> interpolation error in goto 700
566      !          alias interpol_wind (division by zero)
567      if (zt.ge.height(nz)) zt=height(nz)-100.*eps
568
569      if (zt.gt.h) then
570        if (itimec.eq.itime+lsynctime) goto 99
571        goto 700    ! complete the current interval above PBL
572      endif
573
574
575  !!! CHANGE: TEST OF THE WELL-MIXED CRITERION
576  !!! These lines may be switched on to test the well-mixed criterion
577  !if (zt.le.h) then
578  !  zacc=zacc+zt/h*dt
579  !  hsave=hsave+h*dt
580  !  tacc=tacc+dt
581  !  do 67 i=1,iclass
582  !    if ((zt/h.gt.th(i-1)).and.(zt/h.le.th(i)))
583  !    +    t(i)=t(i)+dt
584  !7        continue
585  !endif
586  !if ((mod(itime,10800).eq.0).and.dump) then
587  ! dump=.false.
588  ! write(550,'(i6,12f10.3)') itime,hsave/tacc,zacc/tacc,
589  !    + (t(i)/tacc*real(iclass),i=1,iclass)
590  !  zacc=0.
591  !  tacc=0.
592  !  do 68 i=1,iclass
593  !8        t(i)=0.
594  !  hsave=0.
595  !endif
596  !if (mod(itime,10800).ne.0) dump=.true.
597  !!! CHANGE
598     
599  ! Determine probability of deposition
600  !************************************
601
602      if ((DRYDEP).and.(zt.lt.2.*href)) then
603        do ks=1,nspec
604          if (DRYDEPSPEC(ks)) then
605            if (depoindicator(ks)) then
606              if (ngrid.le.0) then
607                call interpol_vdep(ks,vdepo(ks))
608              else
609                call interpol_vdep_nests(ks,vdepo(ks))
610              endif
611            endif
612  ! correction by Petra Seibert, 10 April 2001
613  !   this formulation means that prob(n) = 1 - f(0)*...*f(n)
614  !   where f(n) is the exponential term
615               prob(ks)=1.+(prob(ks)-1.)* &
616                    exp(-vdepo(ks)*abs(dt)/(2.*href))
617          endif
618        end do
619      endif
620
621      if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
622
623      if (itimec.eq.(itime+lsynctime)) then
624        usig=0.5*(usigprof(indzp)+usigprof(indz))
625        vsig=0.5*(vsigprof(indzp)+vsigprof(indz))
626        wsig=0.5*(wsigprof(indzp)+wsigprof(indz))
627        goto 99  ! finished
628      endif
629      goto 100
630
631  ! END TIME LOOP
632  !==============
633
634
635  endif
636
637
638
639  !**********************************************************
640  ! For all particles that are outside the PBL, make a single
641  ! time step. Only horizontal turbulent disturbances are
642  ! calculated. Vertical disturbances are reset.
643  !**********************************************************
644
645
646  ! Interpolate the wind
647  !*********************
648
649700   continue
650  if (ngrid.le.0) then
651    xts=real(xt)
652    yts=real(yt)
653    call interpol_wind(itime,xts,yts,zt)
654  else
655    call interpol_wind_nests(itime,xtn,ytn,zt)
656  endif
657
658
659  ! Compute everything for above the PBL
660
661  ! Assume constant, uncorrelated, turbulent perturbations
662  ! In the stratosphere, use a small vertical diffusivity d_strat,
663  ! in the troposphere, use a larger horizontal diffusivity d_trop.
664  ! Turbulent velocity scales are determined based on sqrt(d_trop/dt)
665  !******************************************************************
666
667  ldt=abs(lsynctime-itimec+itime)
668  dt=real(ldt)
669
670  if (zt.lt.tropop) then  ! in the troposphere
671    uxscale=sqrt(2.*d_trop/dt)
672    if (nrand+1.gt.maxrand) nrand=1
673    ux=rannumb(nrand)*uxscale
674    vy=rannumb(nrand+1)*uxscale
675    nrand=nrand+2
676    wp=0.
677  else if (zt.lt.tropop+1000.) then     ! just above the tropopause: make transition
678    weight=(zt-tropop)/1000.
679    uxscale=sqrt(2.*d_trop/dt*(1.-weight))
680    if (nrand+2.gt.maxrand) nrand=1
681    ux=rannumb(nrand)*uxscale
682    vy=rannumb(nrand+1)*uxscale
683    wpscale=sqrt(2.*d_strat/dt*weight)
684    wp=rannumb(nrand+2)*wpscale+d_strat/1000.
685    nrand=nrand+3
686  else                 ! in the stratosphere
687    if (nrand.gt.maxrand) nrand=1
688    ux=0.
689    vy=0.
690    wpscale=sqrt(2.*d_strat/dt)
691    wp=rannumb(nrand)*wpscale
692    nrand=nrand+1
693  endif
694
695  if (turboff) then
696!sec switch off turbulence
697    ux=0.0
698    vy=0.0
699    wp=0.0
700  endif
701
702  ! If particle represents only a single species, add gravitational settling
703  ! velocity. The settling velocity is zero for gases
704  !*************************************************************************
705
706  if (mdomainfill.eq.0) then
707    if (lsettling) then
708      do nsp=1,nspec
709        if (xmass(nrelpoint,nsp).gt.eps2) exit
710      end do
711      if (nsp.gt.nspec) then
712  ! This should never happen         
713        write(*,*) 'advance.f90: ERROR: could not find releasepoint'
714        stop
715      end if
716      if (density(nsp).gt.0.) then
717        call get_settling(itime,real(xt),real(yt),zt,nsp,settling)  !bugfix
718        w=w+settling
719      end if
720    endif
721  end if
722
723
724  ! Calculate position at time step itime+lsynctime
725  !************************************************
726
727  dxsave=dxsave+(u+ux)*dt
728  dysave=dysave+(v+vy)*dt
729  zt=zt+(w+wp)*dt*real(ldirect)
730  if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
731
73299   continue
733
734
735
736  !****************************************************************
737  ! Add mesoscale random disturbances
738  ! This is done only once for the whole lsynctime interval to save
739  ! computation time
740  !****************************************************************
741
742
743  ! Mesoscale wind velocity fluctuations are obtained by scaling
744  ! with the standard deviation of the grid-scale winds surrounding
745  ! the particle location, multiplied by a factor turbmesoscale.
746  ! The autocorrelation time constant is taken as half the
747  ! time interval between wind fields
748  !****************************************************************
749
750  r=exp(-2.*real(abs(lsynctime))/real(lwindinterv))
751  rs=sqrt(1.-r**2)
752  if (nrand+2.gt.maxrand) nrand=1
753  usigold=r*usigold+rs*rannumb(nrand)*usig*turbmesoscale
754  vsigold=r*vsigold+rs*rannumb(nrand+1)*vsig*turbmesoscale
755  wsigold=r*wsigold+rs*rannumb(nrand+2)*wsig*turbmesoscale
756
757  dxsave=dxsave+usigold*real(lsynctime)
758  dysave=dysave+vsigold*real(lsynctime)
759
760  zt=zt+wsigold*real(lsynctime)
761  if (zt.lt.0.) zt=-1.*zt    ! if particle below ground -> refletion
762
763  !*************************************************************
764  ! Transform along and cross wind components to xy coordinates,
765  ! add them to u and v, transform u,v to grid units/second
766  ! and calculate new position
767  !*************************************************************
768
769  call windalign(dxsave,dysave,dawsave,dcwsave,ux,vy)
770  dxsave=dxsave+ux   ! comment by mc: comment this line to stop the particles horizontally for test reasons
771  dysave=dysave+vy
772  if (ngrid.ge.0) then
773    cosfact=dxconst/cos((yt*dy+ylat0)*pi180)
774    xt=xt+real(dxsave*cosfact*real(ldirect),kind=dp)
775    yt=yt+real(dysave*dyconst*real(ldirect),kind=dp)
776  else if (ngrid.eq.-1) then      ! around north pole
777    xlon=xlon0+xt*dx                                !comment by mc: compute old particle position
778    ylat=ylat0+yt*dy
779    call cll2xy(northpolemap,ylat,xlon,xpol,ypol)   !convert old particle position in polar stereographic
780    gridsize=1000.*cgszll(northpolemap,ylat,xlon)   !calculate size in m of grid element in polar stereographic coordinate
781    dxsave=dxsave/gridsize                          !increment from meter to grdi unit
782    dysave=dysave/gridsize
783    xpol=xpol+dxsave*real(ldirect)                  !position in grid unit polar stereographic
784    ypol=ypol+dysave*real(ldirect)
785    call cxy2ll(northpolemap,xpol,ypol,ylat,xlon)  !convert to lat long coordinate
786    xt=(xlon-xlon0)/dx                             !convert to grid units in lat long coordinate, comment by mc
787    yt=(ylat-ylat0)/dy
788  else if (ngrid.eq.-2) then    ! around south pole
789    xlon=xlon0+xt*dx
790    ylat=ylat0+yt*dy
791    call cll2xy(southpolemap,ylat,xlon,xpol,ypol)
792    gridsize=1000.*cgszll(southpolemap,ylat,xlon)
793    dxsave=dxsave/gridsize
794    dysave=dysave/gridsize
795    xpol=xpol+dxsave*real(ldirect)
796    ypol=ypol+dysave*real(ldirect)
797    call cxy2ll(southpolemap,xpol,ypol,ylat,xlon)
798    xt=(xlon-xlon0)/dx
799    yt=(ylat-ylat0)/dy
800  endif
801
802
803  ! If global data are available, use cyclic boundary condition
804  !************************************************************
805
806  if (xglobal) then
807    if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1)
808    if (xt.lt.0.) xt=xt+real(nxmin1)
809    if (xt.le.eps) xt=eps
810    if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps
811  endif
812
813  ! HSO/AL: Prevent particles from disappearing at the pole
814  !******************************************************************
815
816  if ( yt.lt.0. ) then
817    xt=mod(xt+180.,360.)
818    yt=-yt
819  else if ( yt.gt.real(nymin1) ) then
820    xt=mod(xt+180.,360.)
821    yt=2*real(nymin1)-yt
822  endif
823
824  ! Check position: If trajectory outside model domain, terminate it
825  !*****************************************************************
826
827  if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. &
828       (yt.gt.real(nymin1))) then
829    nstop=3
830    return
831  endif
832
833  ! If particle above highest model level, set it back into the domain
834  !*******************************************************************
835
836  if (zt.ge.height(nz)) zt=height(nz)-100.*eps
837
838
839  !************************************************************************
840  ! Now we could finish, as this was done in FLEXPART versions up to 4.0.
841  ! However, truncation errors of the advection can be significantly
842  ! reduced by doing one iteration of the Petterssen scheme, if this is
843  ! possible.
844  ! Note that this is applied only to the grid-scale winds, not to
845  ! the turbulent winds.
846  !************************************************************************
847
848  ! The Petterssen scheme can only applied with long time steps (only then u
849  ! is the "old" wind as required by the scheme); otherwise do nothing
850  !*************************************************************************
851
852  if (ldt.ne.abs(lsynctime)) return
853
854  ! The Petterssen scheme can only be applied if the ending time of the time step
855  ! (itime+ldt*ldirect) is still between the two wind fields held in memory;
856  ! otherwise do nothing
857  !******************************************************************************
858
859  if (abs(itime+ldt*ldirect).gt.abs(memtime(2))) return
860
861  ! Apply it also only if starting and ending point of current time step are on
862  ! the same grid; otherwise do nothing
863  !*****************************************************************************
864  if (nglobal.and.(yt.gt.switchnorthg)) then
865    ngr=-1
866  else if (sglobal.and.(yt.lt.switchsouthg)) then
867    ngr=-2
868  else
869    ngr=0
870    do j=numbnests,1,-1
871      if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. &
872           (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then
873        ngr=j
874        goto 43
875      endif
876    end do
87743   continue
878  endif
879
880  if (ngr.ne.ngrid) return
881
882  ! Determine nested grid coordinates
883  !**********************************
884
885  if (ngrid.gt.0) then
886    xtn=(xt-xln(ngrid))*xresoln(ngrid)
887    ytn=(yt-yln(ngrid))*yresoln(ngrid)
888    ix=int(xtn)
889    jy=int(ytn)
890  else
891    ix=int(xt)
892    jy=int(yt)
893  endif
894  ixp=ix+1
895  jyp=jy+1
896
897
898  ! Memorize the old wind
899  !**********************
900
901  uold=u
902  vold=v
903  wold=w
904
905  ! Interpolate wind at new position and time
906  !******************************************
907
908  if (ngrid.le.0) then
909    xts=real(xt)
910    yts=real(yt)
911    call interpol_wind_short(itime+ldt*ldirect,xts,yts,zt)
912  else
913    call interpol_wind_short_nests(itime+ldt*ldirect,xtn,ytn,zt)
914  endif
915
916  if (mdomainfill.eq.0) then
917    if (lsettling) then
918      do nsp=1,nspec
919        if (xmass(nrelpoint,nsp).gt.eps2) exit
920      end do
921      if (nsp.gt.nspec) then
922  ! This should never happen         
923        write(*,*) 'advance.f90: ERROR: could not find releasepoint'
924        stop
925      end if
926      if (density(nsp).gt.0.) then
927        call get_settling(itime+ldt,real(xt),real(yt),zt,nsp,settling) !bugfix
928        w=w+settling
929      end if
930    endif
931  end if
932
933
934  ! Determine the difference vector between new and old wind
935  ! (use half of it to correct position according to Petterssen)
936  !*************************************************************
937
938  u=(u-uold)/2.
939  v=(v-vold)/2.
940  w=(w-wold)/2.
941
942
943  ! Finally, correct the old position
944  !**********************************
945
946  zt=zt+w*real(ldt*ldirect)
947  if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
948  if (ngrid.ge.0) then
949    cosfact=dxconst/cos((yt*dy+ylat0)*pi180)
950    xt=xt+real(u*cosfact*real(ldt*ldirect),kind=dp)
951    yt=yt+real(v*dyconst*real(ldt*ldirect),kind=dp)
952  else if (ngrid.eq.-1) then      ! around north pole
953    xlon=xlon0+xt*dx
954    ylat=ylat0+yt*dy
955    call cll2xy(northpolemap,ylat,xlon,xpol,ypol)
956    gridsize=1000.*cgszll(northpolemap,ylat,xlon)
957    u=u/gridsize
958    v=v/gridsize
959    xpol=xpol+u*real(ldt*ldirect)
960    ypol=ypol+v*real(ldt*ldirect)
961    call cxy2ll(northpolemap,xpol,ypol,ylat,xlon)
962    xt=(xlon-xlon0)/dx
963    yt=(ylat-ylat0)/dy
964  else if (ngrid.eq.-2) then    ! around south pole
965    xlon=xlon0+xt*dx
966    ylat=ylat0+yt*dy
967    call cll2xy(southpolemap,ylat,xlon,xpol,ypol)
968    gridsize=1000.*cgszll(southpolemap,ylat,xlon)
969    u=u/gridsize
970    v=v/gridsize
971    xpol=xpol+u*real(ldt*ldirect)
972    ypol=ypol+v*real(ldt*ldirect)
973    call cxy2ll(southpolemap,xpol,ypol,ylat,xlon)
974    xt=(xlon-xlon0)/dx
975    yt=(ylat-ylat0)/dy
976  endif
977
978  ! If global data are available, use cyclic boundary condition
979  !************************************************************
980
981  if (xglobal) then
982    if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1)
983    if (xt.lt.0.) xt=xt+real(nxmin1)
984    if (xt.le.eps) xt=eps
985    if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps
986  endif
987
988  ! HSO/AL: Prevent particles from disappearing at the pole
989  !******************************************************************
990
991  if ( yt.lt.0. ) then
992    xt=mod(xt+180.,360.)
993    yt=-yt
994  else if ( yt.gt.real(nymin1) ) then
995    xt=mod(xt+180.,360.)
996    yt=2*real(nymin1)-yt
997  endif
998
999  ! Check position: If trajectory outside model domain, terminate it
1000  !*****************************************************************
1001
1002  if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. &
1003       (yt.gt.real(nymin1))) then
1004    nstop=3
1005    return
1006  endif
1007
1008  ! If particle above highest model level, set it back into the domain
1009  !*******************************************************************
1010
1011  if (zt.ge.height(nz)) zt=height(nz)-100.*eps
1012
1013
1014end subroutine advance
1015
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG