advance.f90 34.6 KB
Newer Older
Matthias Langer's avatar
 
Matthias Langer committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
!                                                                     *
! This file is part of FLEXPART.                                      *
!                                                                     *
! FLEXPART is free software: you can redistribute it and/or modify    *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or   *
! (at your option) any later version.                                 *
!                                                                     *
! FLEXPART is distributed in the hope that it will be useful,         *
! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
! GNU General Public License for more details.                        *
!                                                                     *
! You should have received a copy of the GNU General Public License   *
! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
!**********************************************************************

subroutine advance(itime,nrelpoint,ldt,up,vp,wp, &
       usigold,vsigold,wsigold,nstop,xt,yt,zt,prob,icbt)
  !                     i    i  i/oi/oi/o
  !  i/o     i/o     i/o     o  i/oi/oi/o i/o  i/o
  !*****************************************************************************
  !                                                                            *
  !  Calculation of turbulent particle trajectories utilizing a                *
  !  zero-acceleration scheme, which is corrected by a numerically more        *
  !  accurate Petterssen scheme whenever possible.                             *
  !                                                                            *
  !  Particle positions are read in, incremented, and returned to the calling  *
  !  program.                                                                  *
  !                                                                            *
  !  In different regions of the atmosphere (PBL vs. free troposphere),        *
  !  different parameters are needed for advection, parameterizing turbulent   *
  !  velocities, etc. For efficiency, different interpolation routines have    *
  !  been written for these different cases, with the disadvantage that there  *
  !  exist several routines doing almost the same. They all share the          *
  !  included file 'interpol_mod'. The following                               *
  !  interpolation routines are used:                                          *
  !                                                                            *
  !  interpol_all(_nests)     interpolates everything (called inside the PBL)  *
  !  interpol_misslev(_nests) if a particle moves vertically in the PBL,       *
  !                           additional parameters are interpolated if it     *
  !                           crosses a model level                            *
  !  interpol_wind(_nests)    interpolates the wind and determines the         *
  !                           standard deviation of the wind (called outside   *
  !                           PBL) also interpolates potential vorticity       *
  !  interpol_wind_short(_nests) only interpolates the wind (needed for the    *
  !                           Petterssen scheme)                               *
  !  interpol_vdep(_nests)    interpolates deposition velocities               *
  !                                                                            *
  !                                                                            *
  !     Author: A. Stohl                                                       *
  !                                                                            *
  !     16 December 1997                                                       *
  !                                                                            *
  !  Changes:                                                                  *
  !                                                                            *
  !  8 April 2000: Deep convection parameterization                            *
  !                                                                            *
  !  May 2002: Petterssen scheme introduced                                    *
  !                                                                            *
  !*****************************************************************************
  !                                                                            *
  ! Variables:                                                                 *
  ! icbt               1 if particle not transferred to forbidden state,       *
  !                    else -1                                                 *
  ! dawsave            accumulated displacement in along-wind direction        *
  ! dcwsave            accumulated displacement in cross-wind direction        *
  ! dxsave             accumulated displacement in longitude                   *
  ! dysave             accumulated displacement in latitude                    *
  ! h [m]              Mixing height                                           *
  ! lwindinterv [s]    time interval between two wind fields                   *
  ! itime [s]          time at which this subroutine is entered                *
  ! itimec [s]         actual time, which is incremented in this subroutine    *
  ! href [m]           height for which dry deposition velocity is calculated  *
  ! ladvance [s]       Total integration time period                           *
  ! ldirect            1 forward, -1 backward                                  *
  ! ldt [s]            Time step for the next integration                      *
  ! lsynctime [s]      Synchronisation interval of FLEXPART                    *
  ! ngrid              index which grid is to be used                          *
  ! nrand              index for a variable to be picked from rannumb          *
  ! nstop              if > 1 particle has left domain and must be stopped     *
  ! prob               probability of absorption due to dry deposition         *
  ! rannumb(maxrand)   normally distributed random variables                   *
  ! rhoa               air density                                             *
  ! rhograd            vertical gradient of the air density                    *
  ! up,vp,wp           random velocities due to turbulence (along wind, cross  *
  !                    wind, vertical wind                                     *
  ! usig,vsig,wsig     mesoscale wind fluctuations                             *
  ! usigold,vsigold,wsigold  like usig, etc., but for the last time step       *
  ! vdepo              Deposition velocities for all species                   *
  ! xt,yt,zt           Particle position                                       *
  !                                                                            *
  !*****************************************************************************

  use point_mod
  use par_mod
  use com_mod
  use interpol_mod
  use hanna_mod
  use cmapf_mod
105
  use random_mod, only: ran3
Matthias Langer's avatar
 
Matthias Langer committed
106
107
108
109
110

  implicit none

  real(kind=dp) :: xt,yt
  real :: zt,xts,yts,weight
111
  integer :: itime,itimec,nstop,ldt,i,j,k,nrand,loop,memindnext,mind
Matthias Langer's avatar
 
Matthias Langer committed
112
113
114
115
116
117
118
  integer :: ngr,nix,njy,ks,nsp,nrelpoint
  real :: dz,dz1,dz2,xlon,ylat,xpol,ypol,gridsize
  real :: ru,rv,rw,dt,ux,vy,cosfact,xtn,ytn,tropop
  real :: prob(maxspec),up,vp,wp,dxsave,dysave,dawsave
  real :: dcwsave
  real :: usigold,vsigold,wsigold,r,rs
  real :: uold,vold,wold,vdepo(maxspec)
Sabine's avatar
Sabine committed
119
  real :: h1(2)
Matthias Langer's avatar
 
Matthias Langer committed
120
121
122
  !real uprof(nzmax),vprof(nzmax),wprof(nzmax)
  !real usigprof(nzmax),vsigprof(nzmax),wsigprof(nzmax)
  !real rhoprof(nzmax),rhogradprof(nzmax)
123
  real :: rhoa,rhograd,delz,dtf,rhoaux,dtftlw,uxscale,wpscale
Matthias Langer's avatar
 
Matthias Langer committed
124
125
  integer(kind=2) :: icbt
  real,parameter :: eps=nxmax/3.e5,eps2=1.e-9
126
127
128
  real :: ptot_lhh,Q_lhh,phi_lhh,ath,bth !modified by mc 
  real :: old_wp_buf,dcas,dcas1,del_test !added by mc
  integer :: i_well,jj,flagrein !test well mixed: modified by mc
Matthias Langer's avatar
 
Matthias Langer committed
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225


  !!! CHANGE: TEST OF THE WELL-MIXED CRITERION
  !  integer,parameter :: iclass=10
  !  real(kind=dp) :: zacc,tacc,t(iclass),th(0:iclass),hsave
  !  logical dump
  !  save zacc,tacc,t,th,hsave,dump
  !!! CHANGE

  integer :: idummy = -7
  real    :: settling = 0.


  !!! CHANGE: TEST OF THE WELL-MIXED CRITERION
  !if (idummy.eq.-7) then
  !open(550,file='WELLMIXEDTEST')
  !do 17 i=0,iclass
  !7      th(i)=real(i)/real(iclass)
  !endif
  !!! CHANGE


  nstop=0
  do i=1,nmixz
    indzindicator(i)=.true.
  end do


  if (DRYDEP) then    ! reset probability for deposition
    do ks=1,nspec
      depoindicator(ks)=.true.
      prob(ks)=0.
    end do
  endif

  dxsave=0.           ! reset position displacements
  dysave=0.           ! due to mean wind
  dawsave=0.          ! and turbulent wind
  dcwsave=0.

  itimec=itime

  nrand=int(ran3(idummy)*real(maxrand-1))+1


  ! Determine whether lat/long grid or polarstereographic projection
  ! is to be used
  ! Furthermore, determine which nesting level to be used
  !*****************************************************************

  if (nglobal.and.(yt.gt.switchnorthg)) then
    ngrid=-1
  else if (sglobal.and.(yt.lt.switchsouthg)) then
    ngrid=-2
  else
    ngrid=0
    do j=numbnests,1,-1
      if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. &
           (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then
        ngrid=j
        goto 23
      endif
    end do
23   continue
  endif


  !***************************
  ! Interpolate necessary data
  !***************************

  if (abs(itime-memtime(1)).lt.abs(itime-memtime(2))) then
    memindnext=1
  else
    memindnext=2
  endif

  ! Determine nested grid coordinates
  !**********************************

  if (ngrid.gt.0) then
    xtn=(xt-xln(ngrid))*xresoln(ngrid)
    ytn=(yt-yln(ngrid))*yresoln(ngrid)
    ix=int(xtn)
    jy=int(ytn)
    nix=nint(xtn)
    njy=nint(ytn)
  else
    ix=int(xt)
    jy=int(yt)
    nix=nint(xt)
    njy=nint(yt)
  endif
  ixp=ix+1
  jyp=jy+1


Sabine's avatar
Sabine committed
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
 ! Determine the lower left corner and its distance to the current position
  !*************************************************************************

  ddx=xt-real(ix)
  ddy=yt-real(jy)
  rddx=1.-ddx
  rddy=1.-ddy
  p1=rddx*rddy
  p2=ddx*rddy
  p3=rddx*ddy
  p4=ddx*ddy

 ! Calculate variables for time interpolation
  !*******************************************

  dt1=real(itime-memtime(1))
  dt2=real(memtime(2)-itime)
  dtt=1./(dt1+dt2)

245
246
247
248
249
250
! eso: Temporary fix for particle exactly at north pole
  if (jyp >= nymax) then
    ! write(*,*) 'WARNING: advance.f90 jyp >= nymax. xt,yt:',xt,yt
    jyp=jyp-1
  end if

Matthias Langer's avatar
 
Matthias Langer committed
251
252
253
254
255
256
  ! Compute maximum mixing height around particle position
  !*******************************************************

  h=0.
  if (ngrid.le.0) then
    do k=1,2
257
       mind=memind(k) ! eso: compatibility with 3-field version
258
       if (interpolhmix) then
259
260
261
262
263
264
265
266
267
268
269
             h1(k)=p1*hmix(ix ,jy ,1,mind) &
                 + p2*hmix(ixp,jy ,1,mind) &
                 + p3*hmix(ix ,jyp,1,mind) &
                 + p4*hmix(ixp,jyp,1,mind)
        else
          do j=jy,jyp
            do i=ix,ixp
               if (hmix(i,j,1,mind).gt.h) h=hmix(i,j,1,mind)
            end do
          end do
        endif
Matthias Langer's avatar
 
Matthias Langer committed
270
271
272
273
    end do
    tropop=tropopause(nix,njy,1,1)
  else
    do k=1,2
274
      mind=memind(k)
Matthias Langer's avatar
 
Matthias Langer committed
275
276
      do j=jy,jyp
        do i=ix,ixp
277
          if (hmixn(i,j,1,mind,ngrid).gt.h) h=hmixn(i,j,1,mind,ngrid)
Matthias Langer's avatar
 
Matthias Langer committed
278
279
280
281
282
283
        end do
      end do
    end do
    tropop=tropopausen(nix,njy,1,1,ngrid)
  endif

284
  if (interpolhmix) h=(h1(1)*dt2+h1(2)*dt1)*dtt 
Matthias Langer's avatar
 
Matthias Langer committed
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
  zeta=zt/h



  !*************************************************************
  ! If particle is in the PBL, interpolate once and then make a
  ! time loop until end of interval is reached
  !*************************************************************

  if (zeta.le.1.) then

  ! BEGIN TIME LOOP
  !================

    loop=0
100   loop=loop+1
      if (method.eq.1) then
        ldt=min(ldt,abs(lsynctime-itimec+itime))
        itimec=itimec+ldt*ldirect
      else
        ldt=abs(lsynctime)
        itimec=itime+lsynctime
      endif
      dt=real(ldt)

      zeta=zt/h


      if (loop.eq.1) then
        if (ngrid.le.0) then
          xts=real(xt)
          yts=real(yt)
          call interpol_all(itime,xts,yts,zt)
        else
          call interpol_all_nests(itime,xtn,ytn,zt)
        endif

      else


  ! Determine the level below the current position for u,v,rho
  !***********************************************************

        do i=2,nz
          if (height(i).gt.zt) then
            indz=i-1
            indzp=i
            goto 6
          endif
        end do
6       continue

  ! If one of the levels necessary is not yet available,
  ! calculate it
  !*****************************************************

        do i=indz,indzp
          if (indzindicator(i)) then
            if (ngrid.le.0) then
              call interpol_misslev(i)
            else
              call interpol_misslev_nests(i)
            endif
          endif
        end do
      endif


  ! Vertical interpolation of u,v,w,rho and drhodz
  !***********************************************

  ! Vertical distance to the level below and above current position
  ! both in terms of (u,v) and (w) fields
  !****************************************************************

      dz=1./(height(indzp)-height(indz))
      dz1=(zt-height(indz))*dz
      dz2=(height(indzp)-zt)*dz

      u=dz1*uprof(indzp)+dz2*uprof(indz)
      v=dz1*vprof(indzp)+dz2*vprof(indz)
      w=dz1*wprof(indzp)+dz2*wprof(indz)
      rhoa=dz1*rhoprof(indzp)+dz2*rhoprof(indz)
      rhograd=dz1*rhogradprof(indzp)+dz2*rhogradprof(indz)


  ! Compute the turbulent disturbances
  ! Determine the sigmas and the timescales
  !****************************************

      if (turbswitch) then
        call hanna(zt)
      else
        call hanna1(zt)
      endif


  !*****************************************
  ! Determine the new diffusivity velocities
  !*****************************************

  ! Horizontal components
  !**********************

      if (nrand+1.gt.maxrand) nrand=1
      if (dt/tlu.lt..5) then
        up=(1.-dt/tlu)*up+rannumb(nrand)*sigu*sqrt(2.*dt/tlu)
      else
        ru=exp(-dt/tlu)
        up=ru*up+rannumb(nrand)*sigu*sqrt(1.-ru**2)
      endif
      if (dt/tlv.lt..5) then
        vp=(1.-dt/tlv)*vp+rannumb(nrand+1)*sigv*sqrt(2.*dt/tlv)
      else
        rv=exp(-dt/tlv)
        vp=rv*vp+rannumb(nrand+1)*sigv*sqrt(1.-rv**2)
      endif
      nrand=nrand+2


      if (nrand+ifine.gt.maxrand) nrand=1
      rhoaux=rhograd/rhoa
      dtf=dt*fine

      dtftlw=dtf/tlw

  ! Loop over ifine short time steps for vertical component
  !********************************************************

      do i=1,ifine

  ! Determine the drift velocity and density correction velocity
  !*************************************************************

        if (turbswitch) then
          if (dtftlw.lt..5) then
421
422
  !*************************************************************
  !************** CBL options added by mc see routine cblf90 ***
pesei's avatar
pesei committed
423
            if (iflagcbl.eq.1) then  !modified by mc
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
              if (-h/ol.gt.5) then  !modified by mc
              !if (ol.lt.0.) then   !modified by mc  
              !if (ol.gt.0.) then   !modified by mc : for test
                  !print  *,zt,wp,ath,bth,tlw,dtf,'prima'
                  flagrein=0
                  nrand=nrand+1
                  old_wp_buf=wp
                  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
                  wp=(wp+ath*dtf+bth*rannumb(nrand)*sqrt(dtf))*real(icbt) 
                  ! wp=(wp+ath*dtf+bth*gasdev2(mydum)*sqrt(dtf))*real(icbt) 
                  delz=wp*dtf
                  if (flagrein.eq.1) then
                      call re_initialize_particle(zt,ust,wst,h,sigw,old_wp_buf,nrand,ol)
                      wp=old_wp_buf
                      delz=wp*dtf
                      nan_count=nan_count+1
                  end if
                  !print  *,zt,wp,ath,bth,tlw,dtf,rannumb(nrand+i),icbt
                  !pause                  
              else 
                  nrand=nrand+1
                  old_wp_buf=wp
                  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
                                                                                      !2-but since ldirect =-1 for inverse time and this must be calculated for (-wp) and
                                                                                      !3-the gaussian pdf is symmetric (i.e. pdf(w)=pdf(-w) ldirect can be discarded
                  bth=sigw*rannumb(nrand)*sqrt(2.*dtftlw)
                  wp=(wp+ath*dtf+bth)*real(icbt)  
                  delz=wp*dtf
                  del_test=(1.-wp)/wp !catch infinity value
                  if (isnan(wp).or.isnan(del_test)) then 
                      nrand=nrand+1                      
                      wp=sigw*rannumb(nrand)
                      delz=wp*dtf
                      nan_count2=nan_count2+1
                      !print *,'NaN coutner equal to:', nan_count,'reduce ifine if this number became a non-negligible fraction of the particle number'
                  end if  
              end if
  !******************** END CBL option *******************************            
  !*******************************************************************            
            else
                 wp=((1.-dtftlw)*wp+rannumb(nrand+i)*sqrt(2.*dtftlw) &
                 +dtf*(dsigwdz+rhoaux*sigw))*real(icbt) 
                 delz=wp*sigw*dtf
            end if
Matthias Langer's avatar
 
Matthias Langer committed
468
469
470
471
          else
            rw=exp(-dtftlw)
            wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2) &
                 +tlw*(1.-rw)*(dsigwdz+rhoaux*sigw))*real(icbt)
472
            delz=wp*sigw*dtf
Matthias Langer's avatar
 
Matthias Langer committed
473
          endif
474
          
Matthias Langer's avatar
 
Matthias Langer committed
475
476
477
478
479
480
481
        else
          rw=exp(-dtftlw)
          wp=(rw*wp+rannumb(nrand+i)*sqrt(1.-rw**2)*sigw &
               +tlw*(1.-rw)*(dsigw2dz+rhoaux*sigw**2))*real(icbt)
          delz=wp*dtf
        endif

482
483
484
485
486
487
488
489
        if (turboff) then
!sec switch off turbulence
          up=0.0
          vp=0.0
          wp=0.0
          delz=0.
        endif

Matthias Langer's avatar
 
Matthias Langer committed
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
  !****************************************************
  ! Compute turbulent vertical displacement of particle
  !****************************************************

        if (abs(delz).gt.h) delz=mod(delz,h)

  ! Determine if particle transfers to a "forbidden state" below the ground
  ! or above the mixing height
  !************************************************************************

        if (delz.lt.-zt) then         ! reflection at ground
          icbt=-1
          zt=-zt-delz
        else if (delz.gt.(h-zt)) then ! reflection at h
          icbt=-1
          zt=-zt-delz+2.*h
        else                         ! no reflection
          icbt=1
          zt=zt+delz
        endif

        if (i.ne.ifine) then
          zeta=zt/h
          call hanna_short(zt)
        endif

      end do
pesei's avatar
pesei committed
517
      if (iflagcbl.ne.1) nrand=nrand+i
Matthias Langer's avatar
 
Matthias Langer committed
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536

  ! Determine time step for next integration
  !*****************************************

      if (turbswitch) then
        ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), &
             0.5/abs(dsigwdz))*ctl)
      else
        ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5))*ctl)
      endif
      ldt=max(ldt,mintime)


  ! If particle represents only a single species, add gravitational settling
  ! velocity. The settling velocity is zero for gases, or if particle
  ! represents more than one species
  !*************************************************************************

      if (mdomainfill.eq.0) then
537
538
539
540
541
542
543
544
545
546
547
548
549
        if (lsettling) then
          do nsp=1,nspec
            if (xmass(nrelpoint,nsp).gt.eps2) exit
          end do
          if (nsp.gt.nspec) then
  ! This should never happen          
            write(*,*) 'advance.f90: ERROR: could not find releasepoint'
            stop
          end if
          if (density(nsp).gt.0.) then
            call get_settling(itime,real(xt),real(yt),zt,nsp,settling)  !bugfix
            w=w+settling
          end if
550
        end if
Matthias Langer's avatar
 
Matthias Langer committed
551
552
553
554
555
556
557
558
559
560
561
562
563
564
      endif

  ! Horizontal displacements during time step dt are small real values compared
  ! to the position; adding the two, would result in large numerical errors.
  ! Thus, displacements are accumulated during lsynctime and are added to the
  ! position at the end
  !****************************************************************************

      dxsave=dxsave+u*dt
      dysave=dysave+v*dt
      dawsave=dawsave+up*dt
      dcwsave=dcwsave+vp*dt
      zt=zt+w*dt*real(ldirect)

565
566
567
568
      ! HSO/AL: Particle managed to go over highest level -> interpolation error in goto 700
      !          alias interpol_wind (division by zero)
      if (zt.ge.height(nz)) zt=height(nz)-100.*eps

Matthias Langer's avatar
 
Matthias Langer committed
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
      if (zt.gt.h) then
        if (itimec.eq.itime+lsynctime) goto 99
        goto 700    ! complete the current interval above PBL
      endif


  !!! CHANGE: TEST OF THE WELL-MIXED CRITERION
  !!! These lines may be switched on to test the well-mixed criterion
  !if (zt.le.h) then
  !  zacc=zacc+zt/h*dt
  !  hsave=hsave+h*dt
  !  tacc=tacc+dt
  !  do 67 i=1,iclass
  !    if ((zt/h.gt.th(i-1)).and.(zt/h.le.th(i)))
  !    +    t(i)=t(i)+dt
  !7        continue
  !endif
  !if ((mod(itime,10800).eq.0).and.dump) then
  ! dump=.false.
  ! write(550,'(i6,12f10.3)') itime,hsave/tacc,zacc/tacc,
  !    + (t(i)/tacc*real(iclass),i=1,iclass)
  !  zacc=0.
  !  tacc=0.
  !  do 68 i=1,iclass
  !8        t(i)=0.
  !  hsave=0.
  !endif
  !if (mod(itime,10800).ne.0) dump=.true.
  !!! CHANGE
598
      
Matthias Langer's avatar
 
Matthias Langer committed
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
  ! Determine probability of deposition
  !************************************

      if ((DRYDEP).and.(zt.lt.2.*href)) then
        do ks=1,nspec
          if (DRYDEPSPEC(ks)) then
            if (depoindicator(ks)) then
              if (ngrid.le.0) then
                call interpol_vdep(ks,vdepo(ks))
              else
                call interpol_vdep_nests(ks,vdepo(ks))
              endif
            endif
  ! correction by Petra Seibert, 10 April 2001
  !   this formulation means that prob(n) = 1 - f(0)*...*f(n)
  !   where f(n) is the exponential term
               prob(ks)=1.+(prob(ks)-1.)* &
                    exp(-vdepo(ks)*abs(dt)/(2.*href))
          endif
        end do
      endif

      if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection

      if (itimec.eq.(itime+lsynctime)) then
        usig=0.5*(usigprof(indzp)+usigprof(indz))
        vsig=0.5*(vsigprof(indzp)+vsigprof(indz))
        wsig=0.5*(wsigprof(indzp)+wsigprof(indz))
        goto 99  ! finished
      endif
      goto 100

  ! END TIME LOOP
  !==============


  endif



  !**********************************************************
  ! For all particles that are outside the PBL, make a single
  ! time step. Only horizontal turbulent disturbances are
  ! calculated. Vertical disturbances are reset.
  !**********************************************************


  ! Interpolate the wind
  !*********************

700   continue
  if (ngrid.le.0) then
    xts=real(xt)
    yts=real(yt)
    call interpol_wind(itime,xts,yts,zt)
  else
    call interpol_wind_nests(itime,xtn,ytn,zt)
  endif


  ! Compute everything for above the PBL

  ! Assume constant, uncorrelated, turbulent perturbations
  ! In the stratosphere, use a small vertical diffusivity d_strat,
  ! in the troposphere, use a larger horizontal diffusivity d_trop.
  ! Turbulent velocity scales are determined based on sqrt(d_trop/dt)
  !******************************************************************

  ldt=abs(lsynctime-itimec+itime)
  dt=real(ldt)

  if (zt.lt.tropop) then  ! in the troposphere
    uxscale=sqrt(2.*d_trop/dt)
    if (nrand+1.gt.maxrand) nrand=1
    ux=rannumb(nrand)*uxscale
    vy=rannumb(nrand+1)*uxscale
    nrand=nrand+2
    wp=0.
  else if (zt.lt.tropop+1000.) then     ! just above the tropopause: make transition
    weight=(zt-tropop)/1000.
    uxscale=sqrt(2.*d_trop/dt*(1.-weight))
    if (nrand+2.gt.maxrand) nrand=1
    ux=rannumb(nrand)*uxscale
    vy=rannumb(nrand+1)*uxscale
    wpscale=sqrt(2.*d_strat/dt*weight)
    wp=rannumb(nrand+2)*wpscale+d_strat/1000.
    nrand=nrand+3
  else                 ! in the stratosphere
    if (nrand.gt.maxrand) nrand=1
    ux=0.
    vy=0.
    wpscale=sqrt(2.*d_strat/dt)
    wp=rannumb(nrand)*wpscale
    nrand=nrand+1
  endif

695
696
697
698
699
700
  if (turboff) then
!sec switch off turbulence
    ux=0.0
    vy=0.0
    wp=0.0
  endif
Matthias Langer's avatar
 
Matthias Langer committed
701
702
703
704
705

  ! If particle represents only a single species, add gravitational settling
  ! velocity. The settling velocity is zero for gases
  !*************************************************************************

706
707
708
709
710
711
712
713
714
715
716
717
718
  if (mdomainfill.eq.0) then
    if (lsettling) then
      do nsp=1,nspec
        if (xmass(nrelpoint,nsp).gt.eps2) exit
      end do
      if (nsp.gt.nspec) then
  ! This should never happen          
        write(*,*) 'advance.f90: ERROR: could not find releasepoint'
        stop
      end if
      if (density(nsp).gt.0.) then
        call get_settling(itime,real(xt),real(yt),zt,nsp,settling)  !bugfix
        w=w+settling
719
      end if
Matthias Langer's avatar
 
Matthias Langer committed
720
    endif
721
722
  end if

Matthias Langer's avatar
 
Matthias Langer committed
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769

  ! Calculate position at time step itime+lsynctime
  !************************************************

  dxsave=dxsave+(u+ux)*dt
  dysave=dysave+(v+vy)*dt
  zt=zt+(w+wp)*dt*real(ldirect)
  if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection

99   continue



  !****************************************************************
  ! Add mesoscale random disturbances
  ! This is done only once for the whole lsynctime interval to save
  ! computation time
  !****************************************************************


  ! Mesoscale wind velocity fluctuations are obtained by scaling
  ! with the standard deviation of the grid-scale winds surrounding
  ! the particle location, multiplied by a factor turbmesoscale.
  ! The autocorrelation time constant is taken as half the
  ! time interval between wind fields
  !****************************************************************

  r=exp(-2.*real(abs(lsynctime))/real(lwindinterv))
  rs=sqrt(1.-r**2)
  if (nrand+2.gt.maxrand) nrand=1
  usigold=r*usigold+rs*rannumb(nrand)*usig*turbmesoscale
  vsigold=r*vsigold+rs*rannumb(nrand+1)*vsig*turbmesoscale
  wsigold=r*wsigold+rs*rannumb(nrand+2)*wsig*turbmesoscale

  dxsave=dxsave+usigold*real(lsynctime)
  dysave=dysave+vsigold*real(lsynctime)

  zt=zt+wsigold*real(lsynctime)
  if (zt.lt.0.) zt=-1.*zt    ! if particle below ground -> refletion

  !*************************************************************
  ! Transform along and cross wind components to xy coordinates,
  ! add them to u and v, transform u,v to grid units/second
  ! and calculate new position
  !*************************************************************

  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 
Matthias Langer's avatar
 
Matthias Langer committed
771
772
773
774
775
776
  dysave=dysave+vy
  if (ngrid.ge.0) then
    cosfact=dxconst/cos((yt*dy+ylat0)*pi180)
    xt=xt+real(dxsave*cosfact*real(ldirect),kind=dp)
    yt=yt+real(dysave*dyconst*real(ldirect),kind=dp)
  else if (ngrid.eq.-1) then      ! around north pole
777
    xlon=xlon0+xt*dx                                !comment by mc: compute old particle position
Matthias Langer's avatar
 
Matthias Langer committed
778
    ylat=ylat0+yt*dy
779
780
781
    call cll2xy(northpolemap,ylat,xlon,xpol,ypol)   !convert old particle position in polar stereographic
    gridsize=1000.*cgszll(northpolemap,ylat,xlon)   !calculate size in m of grid element in polar stereographic coordinate
    dxsave=dxsave/gridsize                          !increment from meter to grdi unit
Matthias Langer's avatar
 
Matthias Langer committed
782
    dysave=dysave/gridsize
783
    xpol=xpol+dxsave*real(ldirect)                  !position in grid unit polar stereographic
Matthias Langer's avatar
 
Matthias Langer committed
784
    ypol=ypol+dysave*real(ldirect)
785
786
    call cxy2ll(northpolemap,xpol,ypol,ylat,xlon)  !convert to lat long coordinate
    xt=(xlon-xlon0)/dx                             !convert to grid units in lat long coordinate, comment by mc
Matthias Langer's avatar
 
Matthias Langer committed
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
    yt=(ylat-ylat0)/dy
  else if (ngrid.eq.-2) then    ! around south pole
    xlon=xlon0+xt*dx
    ylat=ylat0+yt*dy
    call cll2xy(southpolemap,ylat,xlon,xpol,ypol)
    gridsize=1000.*cgszll(southpolemap,ylat,xlon)
    dxsave=dxsave/gridsize
    dysave=dysave/gridsize
    xpol=xpol+dxsave*real(ldirect)
    ypol=ypol+dysave*real(ldirect)
    call cxy2ll(southpolemap,xpol,ypol,ylat,xlon)
    xt=(xlon-xlon0)/dx
    yt=(ylat-ylat0)/dy
  endif


  ! If global data are available, use cyclic boundary condition
  !************************************************************

  if (xglobal) then
    if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1)
    if (xt.lt.0.) xt=xt+real(nxmin1)
    if (xt.le.eps) xt=eps
    if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps
  endif

813
814
815
816
817
818
819
820
821
822
  ! HSO/AL: Prevent particles from disappearing at the pole
  !******************************************************************

  if ( yt.lt.0. ) then
    xt=mod(xt+180.,360.)
    yt=-yt
  else if ( yt.gt.real(nymin1) ) then
    xt=mod(xt+180.,360.)
    yt=2*real(nymin1)-yt
  endif
Matthias Langer's avatar
 
Matthias Langer committed
823
824
825
826
827

  ! Check position: If trajectory outside model domain, terminate it
  !*****************************************************************

  if ((xt.lt.0.).or.(xt.ge.real(nxmin1)).or.(yt.lt.0.).or. &
828
       (yt.gt.real(nymin1))) then
Matthias Langer's avatar
 
Matthias Langer committed
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
    nstop=3
    return
  endif

  ! If particle above highest model level, set it back into the domain
  !*******************************************************************

  if (zt.ge.height(nz)) zt=height(nz)-100.*eps


  !************************************************************************
  ! Now we could finish, as this was done in FLEXPART versions up to 4.0.
  ! However, truncation errors of the advection can be significantly
  ! reduced by doing one iteration of the Petterssen scheme, if this is
  ! possible.
  ! Note that this is applied only to the grid-scale winds, not to
  ! the turbulent winds.
  !************************************************************************

  ! The Petterssen scheme can only applied with long time steps (only then u
  ! is the "old" wind as required by the scheme); otherwise do nothing
  !*************************************************************************

  if (ldt.ne.abs(lsynctime)) return

  ! The Petterssen scheme can only be applied if the ending time of the time step
  ! (itime+ldt*ldirect) is still between the two wind fields held in memory;
  ! otherwise do nothing
  !******************************************************************************

  if (abs(itime+ldt*ldirect).gt.abs(memtime(2))) return

  ! Apply it also only if starting and ending point of current time step are on
  ! the same grid; otherwise do nothing
  !*****************************************************************************
  if (nglobal.and.(yt.gt.switchnorthg)) then
    ngr=-1
  else if (sglobal.and.(yt.lt.switchsouthg)) then
    ngr=-2
  else
    ngr=0
    do j=numbnests,1,-1
      if ((xt.gt.xln(j)+eps).and.(xt.lt.xrn(j)-eps).and. &
           (yt.gt.yln(j)+eps).and.(yt.lt.yrn(j)-eps)) then
        ngr=j
        goto 43
      endif
    end do
43   continue
  endif

  if (ngr.ne.ngrid) return

  ! Determine nested grid coordinates
  !**********************************

  if (ngrid.gt.0) then
    xtn=(xt-xln(ngrid))*xresoln(ngrid)
    ytn=(yt-yln(ngrid))*yresoln(ngrid)
    ix=int(xtn)
    jy=int(ytn)
  else
    ix=int(xt)
    jy=int(yt)
  endif
  ixp=ix+1
  jyp=jy+1


  ! Memorize the old wind
  !**********************

  uold=u
  vold=v
  wold=w

  ! Interpolate wind at new position and time
  !******************************************

  if (ngrid.le.0) then
    xts=real(xt)
    yts=real(yt)
    call interpol_wind_short(itime+ldt*ldirect,xts,yts,zt)
  else
    call interpol_wind_short_nests(itime+ldt*ldirect,xtn,ytn,zt)
  endif

  if (mdomainfill.eq.0) then
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
    if (lsettling) then
      do nsp=1,nspec
        if (xmass(nrelpoint,nsp).gt.eps2) exit
      end do
      if (nsp.gt.nspec) then
  ! This should never happen          
        write(*,*) 'advance.f90: ERROR: could not find releasepoint'
        stop
      end if
      if (density(nsp).gt.0.) then
        call get_settling(itime+ldt,real(xt),real(yt),zt,nsp,settling) !bugfix
        w=w+settling
      end if
    endif
  end if
Matthias Langer's avatar
 
Matthias Langer committed
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987


  ! Determine the difference vector between new and old wind
  ! (use half of it to correct position according to Petterssen)
  !*************************************************************

  u=(u-uold)/2.
  v=(v-vold)/2.
  w=(w-wold)/2.


  ! Finally, correct the old position
  !**********************************

  zt=zt+w*real(ldt*ldirect)
  if (zt.lt.0.) zt=min(h-eps2,-1.*zt)    ! if particle below ground -> reflection
  if (ngrid.ge.0) then
    cosfact=dxconst/cos((yt*dy+ylat0)*pi180)
    xt=xt+real(u*cosfact*real(ldt*ldirect),kind=dp)
    yt=yt+real(v*dyconst*real(ldt*ldirect),kind=dp)
  else if (ngrid.eq.-1) then      ! around north pole
    xlon=xlon0+xt*dx
    ylat=ylat0+yt*dy
    call cll2xy(northpolemap,ylat,xlon,xpol,ypol)
    gridsize=1000.*cgszll(northpolemap,ylat,xlon)
    u=u/gridsize
    v=v/gridsize
    xpol=xpol+u*real(ldt*ldirect)
    ypol=ypol+v*real(ldt*ldirect)
    call cxy2ll(northpolemap,xpol,ypol,ylat,xlon)
    xt=(xlon-xlon0)/dx
    yt=(ylat-ylat0)/dy
  else if (ngrid.eq.-2) then    ! around south pole
    xlon=xlon0+xt*dx
    ylat=ylat0+yt*dy
    call cll2xy(southpolemap,ylat,xlon,xpol,ypol)
    gridsize=1000.*cgszll(southpolemap,ylat,xlon)
    u=u/gridsize
    v=v/gridsize
    xpol=xpol+u*real(ldt*ldirect)
    ypol=ypol+v*real(ldt*ldirect)
    call cxy2ll(southpolemap,xpol,ypol,ylat,xlon)
    xt=(xlon-xlon0)/dx
    yt=(ylat-ylat0)/dy
  endif

  ! If global data are available, use cyclic boundary condition
  !************************************************************

  if (xglobal) then
    if (xt.ge.real(nxmin1)) xt=xt-real(nxmin1)
    if (xt.lt.0.) xt=xt+real(nxmin1)
    if (xt.le.eps) xt=eps
    if (abs(xt-real(nxmin1)).le.eps) xt=real(nxmin1)-eps
  endif

988
989
990
991
992
993
994
995
996
997
998
  ! HSO/AL: Prevent particles from disappearing at the pole
  !******************************************************************

  if ( yt.lt.0. ) then
    xt=mod(xt+180.,360.)
    yt=-yt
  else if ( yt.gt.real(nymin1) ) then
    xt=mod(xt+180.,360.)
    yt=2*real(nymin1)-yt
  endif

Matthias Langer's avatar
 
Matthias Langer committed
999
1000
  ! Check position: If trajectory outside model domain, terminate it
  !*****************************************************************