source: flexpart.git/src/advance.f90 @ f9ce123

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since f9ce123 was 5f9d14a, checked in by Espen Sollum ATMOS <eso@…>, 9 years ago

Updated wet depo scheme

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