source: flexpart.git/src/advance.f90 @ 75a4ded

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 75a4ded was 2363b24, checked in by Sabine <sabine.eckhardt@…>, 7 years ago

added hmix interpolation option to com_mod

  • Property mode set to 100644
File size: 34.3 KB
RevLine 
[e200b7a]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
[8a65cb0]105  use random_mod, only: ran3
[e200b7a]106
107  implicit none
108
109  real(kind=dp) :: xt,yt
110  real :: zt,xts,yts,weight
[8a65cb0]111  integer :: itime,itimec,nstop,ldt,i,j,k,nrand,loop,memindnext,mind
[e200b7a]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)
[7725d96]119  real :: h1(2)
[e200b7a]120  !real uprof(nzmax),vprof(nzmax),wprof(nzmax)
121  !real usigprof(nzmax),vsigprof(nzmax),wsigprof(nzmax)
122  !real rhoprof(nzmax),rhogradprof(nzmax)
[8a65cb0]123  real :: rhoa,rhograd,delz,dtf,rhoaux,dtftlw,uxscale,wpscale
[e200b7a]124  integer(kind=2) :: icbt
125  real,parameter :: eps=nxmax/3.e5,eps2=1.e-9
[8a65cb0]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
[e200b7a]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
[7725d96]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
[e200b7a]245  ! Compute maximum mixing height around particle position
246  !*******************************************************
247
248  h=0.
249  if (ngrid.le.0) then
250    do k=1,2
[cea125d]251       mind=memind(k) ! eso: compatibility with 3-field version
[2363b24]252       if (interpolhmix) then
[cea125d]253             h1(k)=p1*hmix(ix ,jy ,1,mind) &
254                 + p2*hmix(ixp,jy ,1,mind) &
255                 + p3*hmix(ix ,jyp,1,mind) &
256                 + p4*hmix(ixp,jyp,1,mind)
257        else
258          do j=jy,jyp
259            do i=ix,ixp
260               if (hmix(i,j,1,mind).gt.h) h=hmix(i,j,1,mind)
261            end do
262          end do
263        endif
[e200b7a]264    end do
265    tropop=tropopause(nix,njy,1,1)
266  else
267    do k=1,2
[8a65cb0]268      mind=memind(k)
[e200b7a]269      do j=jy,jyp
270        do i=ix,ixp
[8a65cb0]271          if (hmixn(i,j,1,mind,ngrid).gt.h) h=hmixn(i,j,1,mind,ngrid)
[e200b7a]272        end do
273      end do
274    end do
275    tropop=tropopausen(nix,njy,1,1,ngrid)
276  endif
277
[cea125d]278  if (interpolhmix) h=(h1(1)*dt2+h1(2)*dt1)*dtt
[e200b7a]279  zeta=zt/h
280
281
282
283  !*************************************************************
284  ! If particle is in the PBL, interpolate once and then make a
285  ! time loop until end of interval is reached
286  !*************************************************************
287
288  if (zeta.le.1.) then
289
290  ! BEGIN TIME LOOP
291  !================
292
293    loop=0
294100   loop=loop+1
295      if (method.eq.1) then
296        ldt=min(ldt,abs(lsynctime-itimec+itime))
297        itimec=itimec+ldt*ldirect
298      else
299        ldt=abs(lsynctime)
300        itimec=itime+lsynctime
301      endif
302      dt=real(ldt)
303
304      zeta=zt/h
305
306
307      if (loop.eq.1) then
308        if (ngrid.le.0) then
309          xts=real(xt)
310          yts=real(yt)
311          call interpol_all(itime,xts,yts,zt)
312        else
313          call interpol_all_nests(itime,xtn,ytn,zt)
314        endif
315
316      else
317
318
319  ! Determine the level below the current position for u,v,rho
320  !***********************************************************
321
322        do i=2,nz
323          if (height(i).gt.zt) then
324            indz=i-1
325            indzp=i
326            goto 6
327          endif
328        end do
3296       continue
330
331  ! If one of the levels necessary is not yet available,
332  ! calculate it
333  !*****************************************************
334
335        do i=indz,indzp
336          if (indzindicator(i)) then
337            if (ngrid.le.0) then
338              call interpol_misslev(i)
339            else
340              call interpol_misslev_nests(i)
341            endif
342          endif
343        end do
344      endif
345
346
347  ! Vertical interpolation of u,v,w,rho and drhodz
348  !***********************************************
349
350  ! Vertical distance to the level below and above current position
351  ! both in terms of (u,v) and (w) fields
352  !****************************************************************
353
354      dz=1./(height(indzp)-height(indz))
355      dz1=(zt-height(indz))*dz
356      dz2=(height(indzp)-zt)*dz
357
358      u=dz1*uprof(indzp)+dz2*uprof(indz)
359      v=dz1*vprof(indzp)+dz2*vprof(indz)
360      w=dz1*wprof(indzp)+dz2*wprof(indz)
361      rhoa=dz1*rhoprof(indzp)+dz2*rhoprof(indz)
362      rhograd=dz1*rhogradprof(indzp)+dz2*rhogradprof(indz)
363
364
365  ! Compute the turbulent disturbances
366  ! Determine the sigmas and the timescales
367  !****************************************
368
369      if (turbswitch) then
370        call hanna(zt)
371      else
372        call hanna1(zt)
373      endif
374
375
376  !*****************************************
377  ! Determine the new diffusivity velocities
378  !*****************************************
379
380  ! Horizontal components
381  !**********************
382
383      if (nrand+1.gt.maxrand) nrand=1
384      if (dt/tlu.lt..5) then
385        up=(1.-dt/tlu)*up+rannumb(nrand)*sigu*sqrt(2.*dt/tlu)
386      else
387        ru=exp(-dt/tlu)
388        up=ru*up+rannumb(nrand)*sigu*sqrt(1.-ru**2)
389      endif
390      if (dt/tlv.lt..5) then
391        vp=(1.-dt/tlv)*vp+rannumb(nrand+1)*sigv*sqrt(2.*dt/tlv)
392      else
393        rv=exp(-dt/tlv)
394        vp=rv*vp+rannumb(nrand+1)*sigv*sqrt(1.-rv**2)
395      endif
396      nrand=nrand+2
397
398
399      if (nrand+ifine.gt.maxrand) nrand=1
400      rhoaux=rhograd/rhoa
401      dtf=dt*fine
402
403      dtftlw=dtf/tlw
404
405  ! Loop over ifine short time steps for vertical component
406  !********************************************************
407
408      do i=1,ifine
409
410  ! Determine the drift velocity and density correction velocity
411  !*************************************************************
412
413        if (turbswitch) then
414          if (dtftlw.lt..5) then
[8a65cb0]415  !*************************************************************
416  !************** CBL options added by mc see routine cblf90 ***
417            if (cblflag.eq.1) then  !modified by mc
418              if (-h/ol.gt.5) then  !modified by mc
419              !if (ol.lt.0.) then   !modified by mc 
420              !if (ol.gt.0.) then   !modified by mc : for test
421                  !print  *,zt,wp,ath,bth,tlw,dtf,'prima'
422                  flagrein=0
423                  nrand=nrand+1
424                  old_wp_buf=wp
425                  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
426                  wp=(wp+ath*dtf+bth*rannumb(nrand)*sqrt(dtf))*real(icbt)
427                  ! wp=(wp+ath*dtf+bth*gasdev2(mydum)*sqrt(dtf))*real(icbt)
428                  delz=wp*dtf
429                  if (flagrein.eq.1) then
430                      call re_initialize_particle(zt,ust,wst,h,sigw,old_wp_buf,nrand,ol)
431                      wp=old_wp_buf
432                      delz=wp*dtf
433                      nan_count=nan_count+1
434                  end if
435                  !print  *,zt,wp,ath,bth,tlw,dtf,rannumb(nrand+i),icbt
436                  !pause                 
437              else
438                  nrand=nrand+1
439                  old_wp_buf=wp
440                  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
441                                                                                      !2-but since ldirect =-1 for inverse time and this must be calculated for (-wp) and
442                                                                                      !3-the gaussian pdf is symmetric (i.e. pdf(w)=pdf(-w) ldirect can be discarded
443                  bth=sigw*rannumb(nrand)*sqrt(2.*dtftlw)
444                  wp=(wp+ath*dtf+bth)*real(icbt) 
445                  delz=wp*dtf
446                  del_test=(1.-wp)/wp !catch infinity value
447                  if (isnan(wp).or.isnan(del_test)) then
448                      nrand=nrand+1                     
449                      wp=sigw*rannumb(nrand)
450                      delz=wp*dtf
451                      nan_count2=nan_count2+1
452                      !print *,'NaN coutner equal to:', nan_count,'reduce ifine if this number became a non-negligible fraction of the particle number'
453                  end if 
454              end if
455  !******************** END CBL option *******************************           
456  !*******************************************************************           
457            else
458                 wp=((1.-dtftlw)*wp+rannumb(nrand+i)*sqrt(2.*dtftlw) &
459                 +dtf*(dsigwdz+rhoaux*sigw))*real(icbt)
460                 delz=wp*sigw*dtf
461            end if
[e200b7a]462          else
463            rw=exp(-dtftlw)
464            wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2) &
465                 +tlw*(1.-rw)*(dsigwdz+rhoaux*sigw))*real(icbt)
[8a65cb0]466            delz=wp*sigw*dtf
[e200b7a]467          endif
[8a65cb0]468         
[e200b7a]469        else
470          rw=exp(-dtftlw)
471          wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2)*sigw &
472               +tlw*(1.-rw)*(dsigw2dz+rhoaux*sigw**2))*real(icbt)
473          delz=wp*dtf
474        endif
475
476  !****************************************************
477  ! Compute turbulent vertical displacement of particle
478  !****************************************************
479
480        if (abs(delz).gt.h) delz=mod(delz,h)
481
482  ! Determine if particle transfers to a "forbidden state" below the ground
483  ! or above the mixing height
484  !************************************************************************
485
486        if (delz.lt.-zt) then         ! reflection at ground
487          icbt=-1
488          zt=-zt-delz
489        else if (delz.gt.(h-zt)) then ! reflection at h
490          icbt=-1
491          zt=-zt-delz+2.*h
492        else                         ! no reflection
493          icbt=1
494          zt=zt+delz
495        endif
496
497        if (i.ne.ifine) then
498          zeta=zt/h
499          call hanna_short(zt)
500        endif
501
502      end do
[8a65cb0]503      if (cblflag.ne.1) nrand=nrand+i
[e200b7a]504
505  ! Determine time step for next integration
506  !*****************************************
507
508      if (turbswitch) then
509        ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), &
510             0.5/abs(dsigwdz))*ctl)
511      else
512        ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5))*ctl)
513      endif
514      ldt=max(ldt,mintime)
515
516
517  ! If particle represents only a single species, add gravitational settling
518  ! velocity. The settling velocity is zero for gases, or if particle
519  ! represents more than one species
520  !*************************************************************************
521
522      if (mdomainfill.eq.0) then
[18adf60]523! ESO 05.2015  Changed this to fix MQUASILAG option, where nrelpoint is
524!              particle number and thus xmass array goes out of bounds
525!        do nsp=1,nspec
526!           if (xmass(nrelpoint,nsp).gt.eps2) goto 887
527!         end do
528! 887     nsp=min(nsp,nspec)
529        if (nspec.eq.1.and.density(1).gt.0.) then
[a652cd5]530          call get_settling(itime,real(xt),real(yt),zt,nspec,settling)  !bugfix
[18adf60]531        end if
[e200b7a]532        w=w+settling
533      endif
534
535  ! Horizontal displacements during time step dt are small real values compared
536  ! to the position; adding the two, would result in large numerical errors.
537  ! Thus, displacements are accumulated during lsynctime and are added to the
538  ! position at the end
539  !****************************************************************************
540
541      dxsave=dxsave+u*dt
542      dysave=dysave+v*dt
543      dawsave=dawsave+up*dt
544      dcwsave=dcwsave+vp*dt
545      zt=zt+w*dt*real(ldirect)
546
[8a65cb0]547      ! HSO/AL: Particle managed to go over highest level -> interpolation error in goto 700
548      !          alias interpol_wind (division by zero)
549      if (zt.ge.height(nz)) zt=height(nz)-100.*eps
550
[e200b7a]551      if (zt.gt.h) then
552        if (itimec.eq.itime+lsynctime) goto 99
553        goto 700    ! complete the current interval above PBL
554      endif
555
556
557  !!! CHANGE: TEST OF THE WELL-MIXED CRITERION
558  !!! These lines may be switched on to test the well-mixed criterion
559  !if (zt.le.h) then
560  !  zacc=zacc+zt/h*dt
561  !  hsave=hsave+h*dt
562  !  tacc=tacc+dt
563  !  do 67 i=1,iclass
564  !    if ((zt/h.gt.th(i-1)).and.(zt/h.le.th(i)))
565  !    +    t(i)=t(i)+dt
566  !7        continue
567  !endif
568  !if ((mod(itime,10800).eq.0).and.dump) then
569  ! dump=.false.
570  ! write(550,'(i6,12f10.3)') itime,hsave/tacc,zacc/tacc,
571  !    + (t(i)/tacc*real(iclass),i=1,iclass)
572  !  zacc=0.
573  !  tacc=0.
574  !  do 68 i=1,iclass
575  !8        t(i)=0.
576  !  hsave=0.
577  !endif
578  !if (mod(itime,10800).ne.0) dump=.true.
579  !!! CHANGE
[8a65cb0]580     
[e200b7a]581  ! Determine probability of deposition
582  !************************************
583
584      if ((DRYDEP).and.(zt.lt.2.*href)) then
585        do ks=1,nspec
586          if (DRYDEPSPEC(ks)) then
587            if (depoindicator(ks)) then
588              if (ngrid.le.0) then
589                call interpol_vdep(ks,vdepo(ks))
590              else
591                call interpol_vdep_nests(ks,vdepo(ks))
592              endif
593            endif
594  ! correction by Petra Seibert, 10 April 2001
595  !   this formulation means that prob(n) = 1 - f(0)*...*f(n)
596  !   where f(n) is the exponential term
597               prob(ks)=1.+(prob(ks)-1.)* &
598                    exp(-vdepo(ks)*abs(dt)/(2.*href))
599          endif
600        end do
601      endif
602
603      if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
604
605      if (itimec.eq.(itime+lsynctime)) then
606        usig=0.5*(usigprof(indzp)+usigprof(indz))
607        vsig=0.5*(vsigprof(indzp)+vsigprof(indz))
608        wsig=0.5*(wsigprof(indzp)+wsigprof(indz))
609        goto 99  ! finished
610      endif
611      goto 100
612
613  ! END TIME LOOP
614  !==============
615
616
617  endif
618
619
620
621  !**********************************************************
622  ! For all particles that are outside the PBL, make a single
623  ! time step. Only horizontal turbulent disturbances are
624  ! calculated. Vertical disturbances are reset.
625  !**********************************************************
626
627
628  ! Interpolate the wind
629  !*********************
630
631700   continue
632  if (ngrid.le.0) then
633    xts=real(xt)
634    yts=real(yt)
635    call interpol_wind(itime,xts,yts,zt)
636  else
637    call interpol_wind_nests(itime,xtn,ytn,zt)
638  endif
639
640
641  ! Compute everything for above the PBL
642
643  ! Assume constant, uncorrelated, turbulent perturbations
644  ! In the stratosphere, use a small vertical diffusivity d_strat,
645  ! in the troposphere, use a larger horizontal diffusivity d_trop.
646  ! Turbulent velocity scales are determined based on sqrt(d_trop/dt)
647  !******************************************************************
648
649  ldt=abs(lsynctime-itimec+itime)
650  dt=real(ldt)
651
652  if (zt.lt.tropop) then  ! in the troposphere
653    uxscale=sqrt(2.*d_trop/dt)
654    if (nrand+1.gt.maxrand) nrand=1
655    ux=rannumb(nrand)*uxscale
656    vy=rannumb(nrand+1)*uxscale
657    nrand=nrand+2
658    wp=0.
659  else if (zt.lt.tropop+1000.) then     ! just above the tropopause: make transition
660    weight=(zt-tropop)/1000.
661    uxscale=sqrt(2.*d_trop/dt*(1.-weight))
662    if (nrand+2.gt.maxrand) nrand=1
663    ux=rannumb(nrand)*uxscale
664    vy=rannumb(nrand+1)*uxscale
665    wpscale=sqrt(2.*d_strat/dt*weight)
666    wp=rannumb(nrand+2)*wpscale+d_strat/1000.
667    nrand=nrand+3
668  else                 ! in the stratosphere
669    if (nrand.gt.maxrand) nrand=1
670    ux=0.
671    vy=0.
672    wpscale=sqrt(2.*d_strat/dt)
673    wp=rannumb(nrand)*wpscale
674    nrand=nrand+1
675  endif
676
677
678  ! If particle represents only a single species, add gravitational settling
679  ! velocity. The settling velocity is zero for gases
680  !*************************************************************************
681
682
683
684    if (mdomainfill.eq.0) then
[18adf60]685! ESO 05.2015  Changed this to fix MQUASILAG option, where nrelpoint is
686!              particle number and thus xmass array goes out of bounds
687
688!      do nsp=1,nspec
689!         if (xmass(nrelpoint,nsp).gt.eps2) goto 888
690!       end do
691! 888   nsp=min(nsp,nspec)
692!        if (density(nsp).gt.0.) then
693      if (nspec.eq.1.and.density(1).gt.0.) then
[a652cd5]694        call get_settling(itime,real(xt),real(yt),zt,nspec,settling)  !bugfix
[18adf60]695      end if
[e200b7a]696      w=w+settling
697    endif
698
699  ! Calculate position at time step itime+lsynctime
700  !************************************************
701
702  dxsave=dxsave+(u+ux)*dt
703  dysave=dysave+(v+vy)*dt
704  zt=zt+(w+wp)*dt*real(ldirect)
705  if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
706
70799   continue
708
709
710
711  !****************************************************************
712  ! Add mesoscale random disturbances
713  ! This is done only once for the whole lsynctime interval to save
714  ! computation time
715  !****************************************************************
716
717
718  ! Mesoscale wind velocity fluctuations are obtained by scaling
719  ! with the standard deviation of the grid-scale winds surrounding
720  ! the particle location, multiplied by a factor turbmesoscale.
721  ! The autocorrelation time constant is taken as half the
722  ! time interval between wind fields
723  !****************************************************************
724
725  r=exp(-2.*real(abs(lsynctime))/real(lwindinterv))
726  rs=sqrt(1.-r**2)
727  if (nrand+2.gt.maxrand) nrand=1
728  usigold=r*usigold+rs*rannumb(nrand)*usig*turbmesoscale
729  vsigold=r*vsigold+rs*rannumb(nrand+1)*vsig*turbmesoscale
730  wsigold=r*wsigold+rs*rannumb(nrand+2)*wsig*turbmesoscale
731
732  dxsave=dxsave+usigold*real(lsynctime)
733  dysave=dysave+vsigold*real(lsynctime)
734
735  zt=zt+wsigold*real(lsynctime)
736  if (zt.lt.0.) zt=-1.*zt    ! if particle below ground -> refletion
737
738  !*************************************************************
739  ! Transform along and cross wind components to xy coordinates,
740  ! add them to u and v, transform u,v to grid units/second
741  ! and calculate new position
742  !*************************************************************
743
744  call windalign(dxsave,dysave,dawsave,dcwsave,ux,vy)
[8a65cb0]745  dxsave=dxsave+ux   ! comment by mc: comment this line to stop the particles horizontally for test reasons
[e200b7a]746  dysave=dysave+vy
747  if (ngrid.ge.0) then
748    cosfact=dxconst/cos((yt*dy+ylat0)*pi180)
749    xt=xt+real(dxsave*cosfact*real(ldirect),kind=dp)
750    yt=yt+real(dysave*dyconst*real(ldirect),kind=dp)
751  else if (ngrid.eq.-1) then      ! around north pole
[8a65cb0]752    xlon=xlon0+xt*dx                                !comment by mc: compute old particle position
[e200b7a]753    ylat=ylat0+yt*dy
[8a65cb0]754    call cll2xy(northpolemap,ylat,xlon,xpol,ypol)   !convert old particle position in polar stereographic
755    gridsize=1000.*cgszll(northpolemap,ylat,xlon)   !calculate size in m of grid element in polar stereographic coordinate
756    dxsave=dxsave/gridsize                          !increment from meter to grdi unit
[e200b7a]757    dysave=dysave/gridsize
[8a65cb0]758    xpol=xpol+dxsave*real(ldirect)                  !position in grid unit polar stereographic
[e200b7a]759    ypol=ypol+dysave*real(ldirect)
[8a65cb0]760    call cxy2ll(northpolemap,xpol,ypol,ylat,xlon)  !convert to lat long coordinate
761    xt=(xlon-xlon0)/dx                             !convert to grid units in lat long coordinate, comment by mc
[e200b7a]762    yt=(ylat-ylat0)/dy
763  else if (ngrid.eq.-2) then    ! around south pole
764    xlon=xlon0+xt*dx
765    ylat=ylat0+yt*dy
766    call cll2xy(southpolemap,ylat,xlon,xpol,ypol)
767    gridsize=1000.*cgszll(southpolemap,ylat,xlon)
768    dxsave=dxsave/gridsize
769    dysave=dysave/gridsize
770    xpol=xpol+dxsave*real(ldirect)
771    ypol=ypol+dysave*real(ldirect)
772    call cxy2ll(southpolemap,xpol,ypol,ylat,xlon)
773    xt=(xlon-xlon0)/dx
774    yt=(ylat-ylat0)/dy
775  endif
776
777
778  ! If global data are available, use cyclic boundary condition
779  !************************************************************
780
781  if (xglobal) then
782    if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1)
783    if (xt.lt.0.) xt=xt+real(nxmin1)
784    if (xt.le.eps) xt=eps
785    if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps
786  endif
787
[8a65cb0]788  ! HSO/AL: Prevent particles from disappearing at the pole
789  !******************************************************************
790
791  if ( yt.lt.0. ) then
792    xt=mod(xt+180.,360.)
793    yt=-yt
794  else if ( yt.gt.real(nymin1) ) then
795    xt=mod(xt+180.,360.)
796    yt=2*real(nymin1)-yt
797  endif
[e200b7a]798
799  ! Check position: If trajectory outside model domain, terminate it
800  !*****************************************************************
801
802  if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. &
[8a65cb0]803       (yt.gt.real(nymin1))) then
[e200b7a]804    nstop=3
805    return
806  endif
807
808  ! If particle above highest model level, set it back into the domain
809  !*******************************************************************
810
811  if (zt.ge.height(nz)) zt=height(nz)-100.*eps
812
813
814  !************************************************************************
815  ! Now we could finish, as this was done in FLEXPART versions up to 4.0.
816  ! However, truncation errors of the advection can be significantly
817  ! reduced by doing one iteration of the Petterssen scheme, if this is
818  ! possible.
819  ! Note that this is applied only to the grid-scale winds, not to
820  ! the turbulent winds.
821  !************************************************************************
822
823  ! The Petterssen scheme can only applied with long time steps (only then u
824  ! is the "old" wind as required by the scheme); otherwise do nothing
825  !*************************************************************************
826
827  if (ldt.ne.abs(lsynctime)) return
828
829  ! The Petterssen scheme can only be applied if the ending time of the time step
830  ! (itime+ldt*ldirect) is still between the two wind fields held in memory;
831  ! otherwise do nothing
832  !******************************************************************************
833
834  if (abs(itime+ldt*ldirect).gt.abs(memtime(2))) return
835
836  ! Apply it also only if starting and ending point of current time step are on
837  ! the same grid; otherwise do nothing
838  !*****************************************************************************
839  if (nglobal.and.(yt.gt.switchnorthg)) then
840    ngr=-1
841  else if (sglobal.and.(yt.lt.switchsouthg)) then
842    ngr=-2
843  else
844    ngr=0
845    do j=numbnests,1,-1
846      if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. &
847           (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then
848        ngr=j
849        goto 43
850      endif
851    end do
85243   continue
853  endif
854
855  if (ngr.ne.ngrid) return
856
857  ! Determine nested grid coordinates
858  !**********************************
859
860  if (ngrid.gt.0) then
861    xtn=(xt-xln(ngrid))*xresoln(ngrid)
862    ytn=(yt-yln(ngrid))*yresoln(ngrid)
863    ix=int(xtn)
864    jy=int(ytn)
865  else
866    ix=int(xt)
867    jy=int(yt)
868  endif
869  ixp=ix+1
870  jyp=jy+1
871
872
873  ! Memorize the old wind
874  !**********************
875
876  uold=u
877  vold=v
878  wold=w
879
880  ! Interpolate wind at new position and time
881  !******************************************
882
883  if (ngrid.le.0) then
884    xts=real(xt)
885    yts=real(yt)
886    call interpol_wind_short(itime+ldt*ldirect,xts,yts,zt)
887  else
888    call interpol_wind_short_nests(itime+ldt*ldirect,xtn,ytn,zt)
889  endif
890
891  if (mdomainfill.eq.0) then
[18adf60]892! ESO 05.2015  Changed this to fix MQUASILAG option, where nrelpoint is
893!              particle number and thus xmass array goes out of bounds
894!    do nsp=1,nspec
895!       if (xmass(nrelpoint,nsp).gt.eps2) goto 889
896!     end do
897! 889   nsp=min(nsp,nspec)
898!      if (density(nsp).gt.0.) then
899    if (nspec.eq.1.and.density(1).gt.0.) then
[a652cd5]900      call get_settling(itime+ldt,real(xt),real(yt),zt,nspec,settling)  !bugfix
[18adf60]901    end if
[e200b7a]902    w=w+settling
903  endif
904
905
906  ! Determine the difference vector between new and old wind
907  ! (use half of it to correct position according to Petterssen)
908  !*************************************************************
909
910  u=(u-uold)/2.
911  v=(v-vold)/2.
912  w=(w-wold)/2.
913
914
915  ! Finally, correct the old position
916  !**********************************
917
918  zt=zt+w*real(ldt*ldirect)
919  if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
920  if (ngrid.ge.0) then
921    cosfact=dxconst/cos((yt*dy+ylat0)*pi180)
922    xt=xt+real(u*cosfact*real(ldt*ldirect),kind=dp)
923    yt=yt+real(v*dyconst*real(ldt*ldirect),kind=dp)
924  else if (ngrid.eq.-1) then      ! around north pole
925    xlon=xlon0+xt*dx
926    ylat=ylat0+yt*dy
927    call cll2xy(northpolemap,ylat,xlon,xpol,ypol)
928    gridsize=1000.*cgszll(northpolemap,ylat,xlon)
929    u=u/gridsize
930    v=v/gridsize
931    xpol=xpol+u*real(ldt*ldirect)
932    ypol=ypol+v*real(ldt*ldirect)
933    call cxy2ll(northpolemap,xpol,ypol,ylat,xlon)
934    xt=(xlon-xlon0)/dx
935    yt=(ylat-ylat0)/dy
936  else if (ngrid.eq.-2) then    ! around south pole
937    xlon=xlon0+xt*dx
938    ylat=ylat0+yt*dy
939    call cll2xy(southpolemap,ylat,xlon,xpol,ypol)
940    gridsize=1000.*cgszll(southpolemap,ylat,xlon)
941    u=u/gridsize
942    v=v/gridsize
943    xpol=xpol+u*real(ldt*ldirect)
944    ypol=ypol+v*real(ldt*ldirect)
945    call cxy2ll(southpolemap,xpol,ypol,ylat,xlon)
946    xt=(xlon-xlon0)/dx
947    yt=(ylat-ylat0)/dy
948  endif
949
950  ! If global data are available, use cyclic boundary condition
951  !************************************************************
952
953  if (xglobal) then
954    if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1)
955    if (xt.lt.0.) xt=xt+real(nxmin1)
956    if (xt.le.eps) xt=eps
957    if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps
958  endif
959
[8a65cb0]960  ! HSO/AL: Prevent particles from disappearing at the pole
961  !******************************************************************
962
963  if ( yt.lt.0. ) then
964    xt=mod(xt+180.,360.)
965    yt=-yt
966  else if ( yt.gt.real(nymin1) ) then
967    xt=mod(xt+180.,360.)
968    yt=2*real(nymin1)-yt
969  endif
970
[e200b7a]971  ! Check position: If trajectory outside model domain, terminate it
972  !*****************************************************************
973
974  if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. &
[8a65cb0]975       (yt.gt.real(nymin1))) then
[e200b7a]976    nstop=3
977    return
978  endif
979
980  ! If particle above highest model level, set it back into the domain
981  !*******************************************************************
982
983  if (zt.ge.height(nz)) zt=height(nz)-100.*eps
984
985
986end subroutine advance
987
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG