source: flexpart.git/src/timemanager.f90 @ 2363b24

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 2363b24 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: 29.1 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 timemanager
23
24  !*****************************************************************************
25  !                                                                            *
26  ! Handles the computation of trajectories, i.e. determines which             *
27  ! trajectories have to be computed at what time.                             *
28  ! Manages dry+wet deposition routines, radioactive decay and the computation *
29  ! of concentrations.                                                         *
30  !                                                                            *
31  !     Author: A. Stohl                                                       *
32  !                                                                            *
33  !     20 May 1996                                                            *
34  !                                                                            *
35  !*****************************************************************************
36  !  Changes, Bernd C. Krueger, Feb. 2001:                                     *
37  !        Call of convmix when new windfield is read                          *
38  !------------------------------------                                        *
39  !  Changes Petra Seibert, Sept 2002                                          *
40  !     fix wet scavenging problem                                             *
41  !     Code may not be correct for decay of deposition!                       *
42  !  Changes Petra Seibert, Nov 2002                                           *
43  !     call convection BEFORE new fields are read in BWD mode                 *
44  !  Changes Caroline Forster, Feb 2005                                        *
45  !   new interface between flexpart and convection scheme                     *
46  !   Emanuel's latest subroutine convect43c.f is used                         *
47  !  Changes Stefan Henne, Harald Sodemann, 2013-2014                          *
48  !   added netcdf output code                                                 *
49  !  Changes Espen Sollum 2014                                                 *
50  !   For compatibility with MPI version,                                      *
51  !   variables uap,ucp,uzp,us,vs,ws,cbt now in module com_mod                 *
52  !*****************************************************************************
53  !                                                                            *
54  ! Variables:                                                                 *
55  ! DEP                .true. if either wet or dry deposition is switched on   *
56  ! decay(maxspec) [1/s] decay constant for radioactive decay                  *
57  ! DRYDEP             .true. if dry deposition is switched on                 *
58  ! ideltas [s]        modelling period                                        *
59  ! itime [s]          actual temporal position of calculation                 *
60  ! ldeltat [s]        time since computation of radioact. decay of depositions*
61  ! loutaver [s]       averaging period for concentration calculations         *
62  ! loutend [s]        end of averaging for concentration calculations         *
63  ! loutnext [s]       next time at which output fields shall be centered      *
64  ! loutsample [s]     sampling interval for averaging of concentrations       *
65  ! loutstart [s]      start of averaging for concentration calculations       *
66  ! loutstep [s]       time interval for which concentrations shall be         *
67  !                    calculated                                              *
68  ! npoint(maxpart)    index, which starting point the trajectory has          *
69  !                    starting positions of trajectories                      *
70  ! nstop              serves as indicator for fate of particles               *
71  !                    in the particle loop                                    *
72  ! nstop1             serves as indicator for wind fields (see getfields)     *
73  ! outnum             number of samples for each concentration calculation    *
74  ! outnum             number of samples for each concentration calculation    *
75  ! prob               probability of absorption at ground due to dry          *
76  !                    deposition                                              *
77  ! WETDEP             .true. if wet deposition is switched on                 *
78  ! weight             weight for each concentration sample (1/2 or 1)         *
79  ! uap(maxpart),ucp(maxpart),uzp(maxpart) = random velocities due to          *
80  !                    turbulence                                              *
81  ! us(maxpart),vs(maxpart),ws(maxpart) = random velocities due to inter-      *
82  !                    polation                                                *
83  ! xtra1(maxpart), ytra1(maxpart), ztra1(maxpart) =                           *
84  !                    spatial positions of trajectories                       *
85  !                                                                            *
86  ! Constants:                                                                 *
87  ! maxpart            maximum number of trajectories                          *
88  !                                                                            *
89  !*****************************************************************************
90
91  use unc_mod
92  use point_mod
93  use xmass_mod
94  use flux_mod
95  use outg_mod
96  use oh_mod
97  use par_mod
98  use com_mod
99  use netcdf_output_mod, only: concoutput_netcdf,concoutput_nest_netcdf,&
100       &concoutput_surf_netcdf,concoutput_surf_nest_netcdf
101
102  implicit none
103
104  integer :: j,ks,kp,l,n,itime=0,nstop,nstop1
105! integer :: ksp
106  integer :: loutnext,loutstart,loutend
107  integer :: ix,jy,ldeltat,itage,nage,idummy
108  integer :: i_nan=0,ii_nan,total_nan_intl=0  !added by mc to check instability in CBL scheme
109  real :: outnum,weight,prob_rec(maxspec),prob(maxspec),decfact,wetscav(maxspec)
110  ! real :: uap(maxpart),ucp(maxpart),uzp(maxpart)
111  ! real :: us(maxpart),vs(maxpart),ws(maxpart)
112  ! integer(kind=2) :: cbt(maxpart)
113  real(sp) :: gridtotalunc
114  real(dep_prec) :: drydeposit(maxspec),wetgridtotalunc,drygridtotalunc
115  real :: xold,yold,zold,xmassfract
116  real :: grfraction(3)
117  real, parameter :: e_inv = 1.0/exp(1.0)
118
119  !double precision xm(maxspec,maxpointspec_act),
120  !    +                 xm_depw(maxspec,maxpointspec_act),
121  !    +                 xm_depd(maxspec,maxpointspec_act)
122
123
124  !open(88,file='TEST.dat')
125
126  ! First output for time 0
127  !************************
128
129  loutnext=loutstep/2
130  outnum=0.
131  loutstart=loutnext-loutaver/2
132  loutend=loutnext+loutaver/2
133
134
135  !**********************************************************************
136  ! Loop over the whole modelling period in time steps of mintime seconds
137  !**********************************************************************
138
139!ZHG 2015
140!CGZ-lifetime: set lifetime to 0
141  ! checklifetime(:,:)=0
142  ! species_lifetime(:,:)=0
143  ! print*, 'Initialized lifetime'
144!CGZ-lifetime: set lifetime to 0
145 
146  if (interpolhmix) write(*,*) 'Interpolating hmix'
147
148  write(*,46) float(itime)/3600,itime,numpart
149
150  if (verbosity.gt.0) then
151    write (*,*) 'timemanager> starting simulation'
152    if (verbosity.gt.1) then
153      CALL SYSTEM_CLOCK(count_clock)
154      WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
155    endif     
156  endif
157
158  do itime=0,ideltas,lsynctime
159
160  ! Computation of wet deposition, OH reaction and mass transfer
161  ! between two species every lsynctime seconds
162  ! maybe wet depo frequency can be relaxed later but better be on safe side
163  ! wetdepo must be called BEFORE new fields are read in but should not
164  ! be called in the very beginning before any fields are loaded, or
165  ! before particles are in the system
166  ! Code may not be correct for decay of deposition
167  ! changed by Petra Seibert 9/02
168  !********************************************************************
169
170    if (WETDEP .and. itime .ne. 0 .and. numpart .gt. 0) then
171        if (verbosity.gt.0) then
172           write (*,*) 'timemanager> call wetdepo'
173        endif     
174         call wetdepo(itime,lsynctime,loutnext)
175    endif
176
177    if (OHREA .and. itime .ne. 0 .and. numpart .gt. 0) &
178         call ohreaction(itime,lsynctime,loutnext)
179
180    if (ASSSPEC .and. itime .ne. 0 .and. numpart .gt. 0) then
181       stop 'associated species not yet implemented!'
182  !     call transferspec(itime,lsynctime,loutnext)
183    endif
184
185  ! compute convection for backward runs
186  !*************************************
187
188   if ((ldirect.eq.-1).and.(lconvection.eq.1).and.(itime.lt.0)) then
189        if (verbosity.gt.0) then
190           write (*,*) 'timemanager> call convmix -- backward'
191        endif         
192      call convmix(itime)
193        if (verbosity.gt.1) then
194          !CALL SYSTEM_CLOCK(count_clock, count_rate, count_max)
195          CALL SYSTEM_CLOCK(count_clock)
196          WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
197        endif
198   endif
199
200  ! Get necessary wind fields if not available
201  !*******************************************
202    if (verbosity.gt.0) then
203           write (*,*) 'timemanager> call getfields'
204    endif
205    call getfields(itime,nstop1)
206        if (verbosity.gt.1) then
207          CALL SYSTEM_CLOCK(count_clock)
208          WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
209        endif
210    if (nstop1.gt.1) stop 'NO METEO FIELDS AVAILABLE'
211
212  ! Get hourly OH fields if not available
213  !****************************************************
214    if (OHREA) then
215      if (verbosity.gt.0) then
216             write (*,*) 'timemanager> call gethourlyOH'
217      endif
218      call gethourlyOH(itime)
219          if (verbosity.gt.1) then
220            CALL SYSTEM_CLOCK(count_clock)
221            WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
222          endif
223    endif
224       
225  ! Release particles
226  !******************
227
228    if (verbosity.gt.0) then
229           write (*,*) 'timemanager>  Release particles'
230    endif
231
232    if (mdomainfill.ge.1) then
233      if (itime.eq.0) then
234        if (verbosity.gt.0) then
235          write (*,*) 'timemanager>  call init_domainfill'
236        endif       
237        call init_domainfill
238      else
239        if (verbosity.gt.0) then
240          write (*,*) 'timemanager>  call boundcond_domainfill'
241        endif   
242        call boundcond_domainfill(itime,loutend)
243      endif
244    else
245      if (verbosity.gt.0) then
246        print*,'call releaseparticles' 
247      endif
248      call releaseparticles(itime)
249      if (verbosity.gt.1) then
250        CALL SYSTEM_CLOCK(count_clock)
251        WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
252      endif
253    endif
254
255
256  ! Compute convective mixing for forward runs
257  ! for backward runs it is done before next windfield is read in
258  !**************************************************************
259
260   if ((ldirect.eq.1).and.(lconvection.eq.1)) then
261     if (verbosity.gt.0) then
262       write (*,*) 'timemanager> call convmix -- forward'
263     endif   
264     call convmix(itime)
265   endif
266
267  ! If middle of averaging period of output fields is reached, accumulated
268  ! deposited mass radioactively decays
269  !***********************************************************************
270
271    if (DEP.and.(itime.eq.loutnext).and.(ldirect.gt.0)) then
272      do ks=1,nspec
273      do kp=1,maxpointspec_act
274        if (decay(ks).gt.0.) then
275          do nage=1,nageclass
276            do l=1,nclassunc
277  ! Mother output grid
278              do jy=0,numygrid-1
279                do ix=0,numxgrid-1
280                  wetgridunc(ix,jy,ks,kp,l,nage)= &
281                       wetgridunc(ix,jy,ks,kp,l,nage)* &
282                       exp(-1.*outstep*decay(ks))
283                  drygridunc(ix,jy,ks,kp,l,nage)= &
284                       drygridunc(ix,jy,ks,kp,l,nage)* &
285                       exp(-1.*outstep*decay(ks))
286                end do
287              end do
288  ! Nested output grid
289              if (nested_output.eq.1) then
290                do jy=0,numygridn-1
291                  do ix=0,numxgridn-1
292                    wetgriduncn(ix,jy,ks,kp,l,nage)= &
293                         wetgriduncn(ix,jy,ks,kp,l,nage)* &
294                         exp(-1.*outstep*decay(ks))
295                    drygriduncn(ix,jy,ks,kp,l,nage)= &
296                         drygriduncn(ix,jy,ks,kp,l,nage)* &
297                         exp(-1.*outstep*decay(ks))
298                  end do
299                end do
300              endif
301            end do
302          end do
303        endif
304      end do
305      end do
306    endif
307
308  !!! CHANGE: These lines may be switched on to check the conservation
309  !!! of mass within FLEXPART
310  !   if (itime.eq.loutnext) then
311  !   do 247 ksp=1, nspec
312  !    247 kp=1, maxpointspec_act
313  !47         xm(ksp,kp)=0.
314
315  !   do 249 ksp=1, nspec
316  !     do 249 j=1,numpart
317  !          if (ioutputforeachrelease.eq.1) then
318  !            kp=npoint(j)
319  !          else
320  !            kp=1
321  !          endif
322  !       if (itra1(j).eq.itime) then
323  !          xm(ksp,kp)=xm(ksp,kp)+xmass1(j,ksp)
324  !         write(*,*) 'xmass: ',xmass1(j,ksp),j,ksp,nspec
325  !       endif
326  !49     continue
327  !  do 248 ksp=1,nspec
328  !  do 248 kp=1,maxpointspec_act
329  !  xm_depw(ksp,kp)=0.
330  !  xm_depd(ksp,kp)=0.
331  !     do 248 nage=1,nageclass
332  !       do 248 ix=0,numxgrid-1
333  !         do 248 jy=0,numygrid-1
334  !           do 248 l=1,nclassunc
335  !              xm_depw(ksp,kp)=xm_depw(ksp,kp)
336  !    +                  +wetgridunc(ix,jy,ksp,kp,l,nage)
337  !48                 xm_depd(ksp,kp)=xm_depd(ksp,kp)
338  !    +                  +drygridunc(ix,jy,ksp,kp,l,nage)
339  !             do 246 ksp=1,nspec
340  !46                    write(88,'(2i10,3e12.3)')
341  !    +              itime,ksp,(xm(ksp,kp),kp=1,maxpointspec_act),
342  !    +                (xm_depw(ksp,kp),kp=1,maxpointspec_act),
343  !    +                (xm_depd(ksp,kp),kp=1,maxpointspec_act)
344  !  endif
345  !!! CHANGE
346
347
348
349  ! Check whether concentrations are to be calculated
350  !**************************************************
351
352    if ((ldirect*itime.ge.ldirect*loutstart).and. &
353         (ldirect*itime.le.ldirect*loutend)) then ! add to grid
354      if (mod(itime-loutstart,loutsample).eq.0) then
355
356   
357    do j=1,numpart
358  ! RECEPTOR: dry/wet depovel
359  !****************************
360  ! Before the particle is moved
361  ! the calculation of the scavenged mass shall only be done once after release
362  ! xscav_frac1 was initialised with a negative value
363
364!      write (*,*) 'out:',j,xtra1(j),ytra1(j),prob_rec(1),xscav_frac1(j,1)
365      if  (DRYBKDEP) then
366       do ks=1,nspec
367         if  ((xscav_frac1(j,ks).lt.0)) then
368            call get_vdep_prob(itime,xtra1(j),ytra1(j),ztra1(j),prob_rec)
369            if (DRYDEPSPEC(ks)) then        ! dry deposition
370!        write (*,*) 'j,..:',j,xtra1(j),ytra1(j),prob_rec(1),xscav_frac1(j,1)
371               xscav_frac1(j,ks)=prob_rec(ks)
372             else
373                xmass1(j,ks)=0
374                xscav_frac1(j,ks)=0.
375             endif
376         endif
377        enddo
378       endif
379
380       if (WETBKDEP) then
381       do ks=1,nspec
382         if  ((xscav_frac1(j,ks).lt.0)) then
383            call get_wetscav(itime,lsynctime,loutnext,j,ks,grfraction,idummy,idummy,wetscav)
384            if (wetscav(ks).gt.0) then
385                xscav_frac1(j,ks)=wetscav(ks)* &
386                       (zpoint2(npoint(j))-zpoint1(npoint(j)))*grfraction(1)
387            else
388                xmass1(j,ks)=0.
389                xscav_frac1(j,ks)=0.
390            endif
391         endif
392        enddo
393       endif
394
395  enddo
396
397  ! If we are exactly at the start or end of the concentration averaging interval,
398  ! give only half the weight to this sample
399  !*****************************************************************************
400
401        if ((itime.eq.loutstart).or.(itime.eq.loutend)) then
402          weight=0.5
403        else
404          weight=1.0
405        endif
406        outnum=outnum+weight
407        call conccalc(itime,weight)
408      endif
409
410
411      if ((mquasilag.eq.1).and.(itime.eq.(loutstart+loutend)/2)) &
412           call partoutput_short(itime)    ! dump particle positions in extremely compressed format
413
414
415  ! Output and reinitialization of grid
416  ! If necessary, first sample of new grid is also taken
417  !*****************************************************
418
419      if ((itime.eq.loutend).and.(outnum.gt.0.)) then
420        if ((iout.le.3.).or.(iout.eq.5)) then
421          if (surf_only.ne.1) then
422            if (lnetcdfout.eq.1) then
423              call concoutput_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
424            else
425              call concoutput(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
426            endif
427          else 
428            if (verbosity.eq.1) then
429             print*,'call concoutput_surf '
430             call system_clock(count_clock)
431             write(*,*) 'system clock',count_clock - count_clock0   
432            endif
433            if (lnetcdfout.eq.1) then
434              call concoutput_surf_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
435            else
436              call concoutput_surf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
437              if (verbosity.eq.1) then
438                print*,'called concoutput_surf '
439                call system_clock(count_clock)
440                write(*,*) 'system clock',count_clock - count_clock0   
441              endif
442            endif
443          endif
444
445          if (nested_output .eq. 1) then
446            if (lnetcdfout.eq.0) then
447              if (surf_only.ne.1) then
448                call concoutput_nest(itime,outnum)
449              else
450                call concoutput_surf_nest(itime,outnum)
451              endif
452            else
453              if (surf_only.ne.1) then
454                call concoutput_nest_netcdf(itime,outnum)
455              else
456                call concoutput_surf_nest_netcdf(itime,outnum)
457              endif
458            endif
459          endif
460          outnum=0.
461        endif
462        if ((iout.eq.4).or.(iout.eq.5)) call plumetraj(itime)
463        if (iflux.eq.1) call fluxoutput(itime)
464        write(*,45) itime,numpart,gridtotalunc,wetgridtotalunc,drygridtotalunc
465 
466        !CGZ-lifetime: output species lifetime
467!ZHG
468        ! write(*,*) 'Overview species lifetime in days', &
469        !      real((species_lifetime(:,1)/species_lifetime(:,2))/real(3600.0*24.0))
470        ! write(*,*) 'all info:',species_lifetime
471!ZHG
472        !CGZ-lifetime: output species lifetime
473
474        !write(*,46) float(itime)/3600,itime,numpart
47545      format(i13,' SECONDS SIMULATED: ',i13, ' PARTICLES:    Uncertainty: ',3f7.3)
47646      format(' Simulated ',f7.1,' hours (',i13,' s), ',i13, ' particles')
477        if (ipout.ge.1) call partoutput(itime)    ! dump particle positions
478        loutnext=loutnext+loutstep
479        loutstart=loutnext-loutaver/2
480        loutend=loutnext+loutaver/2
481        if (itime.eq.loutstart) then
482          weight=0.5
483          outnum=outnum+weight
484          call conccalc(itime,weight)
485        endif
486
487
488  ! Check, whether particles are to be split:
489  ! If so, create new particles and attribute all information from the old
490  ! particles also to the new ones; old and new particles both get half the
491  ! mass of the old ones
492  !************************************************************************
493
494        if (ldirect*itime.ge.ldirect*itsplit) then
495          n=numpart
496          do j=1,numpart
497            if (ldirect*itime.ge.ldirect*itrasplit(j)) then
498              if (n.lt.maxpart) then
499                n=n+1
500                itrasplit(j)=2*(itrasplit(j)-itramem(j))+itramem(j)
501                itrasplit(n)=itrasplit(j)
502                itramem(n)=itramem(j)
503                itra1(n)=itra1(j)
504                idt(n)=idt(j)
505                npoint(n)=npoint(j)
506                nclass(n)=nclass(j)
507                xtra1(n)=xtra1(j)
508                ytra1(n)=ytra1(j)
509                ztra1(n)=ztra1(j)
510                uap(n)=uap(j)
511                ucp(n)=ucp(j)
512                uzp(n)=uzp(j)
513                us(n)=us(j)
514                vs(n)=vs(j)
515                ws(n)=ws(j)
516                cbt(n)=cbt(j)
517                do ks=1,nspec
518                  xmass1(j,ks)=xmass1(j,ks)/2.
519                  xmass1(n,ks)=xmass1(j,ks)
520                end do
521              endif
522            endif
523          end do
524          numpart=n
525        endif
526      endif
527    endif
528
529
530    if (itime.eq.ideltas) exit         ! almost finished
531
532  ! Compute interval since radioactive decay of deposited mass was computed
533  !************************************************************************
534
535    if (itime.lt.loutnext) then
536      ldeltat=itime-(loutnext-loutstep)
537    else                                  ! first half of next interval
538      ldeltat=itime-loutnext
539    endif
540
541
542  ! Loop over all particles
543  !************************
544  ! Various variables for testing reason of CBL scheme, by mc
545    well_mixed_vector=0. !erase vector to test well mixed condition: modified by mc
546    well_mixed_norm=0.   !erase normalization to test well mixed condition: modified by mc
547    avg_ol=0.
548    avg_wst=0.
549    avg_h=0.
550    avg_air_dens=0.  !erase vector to obtain air density at particle positions: modified by mc
551  !-----------------------------------------------------------------------------
552    do j=1,numpart
553
554
555  ! If integration step is due, do it
556  !**********************************
557
558      if (itra1(j).eq.itime) then
559
560        if (ioutputforeachrelease.eq.1) then
561            kp=npoint(j)
562        else
563            kp=1
564        endif
565  ! Determine age class of the particle
566        itage=abs(itra1(j)-itramem(j))
567        do nage=1,nageclass
568          if (itage.lt.lage(nage)) exit
569        end do
570
571  ! Initialize newly released particle
572  !***********************************
573
574        if ((itramem(j).eq.itime).or.(itime.eq.0)) &
575             call initialize(itime,idt(j),uap(j),ucp(j),uzp(j), &
576             us(j),vs(j),ws(j),xtra1(j),ytra1(j),ztra1(j),cbt(j))
577
578  ! Memorize particle positions
579  !****************************
580
581        xold=xtra1(j)
582        yold=ytra1(j)
583        zold=ztra1(j)
584
585
586  ! Integrate Lagevin equation for lsynctime seconds
587  !*************************************************
588
589        if (verbosity.gt.0) then
590           if (j.eq.1) then
591             write (*,*) 'timemanager> call advance'
592           endif     
593        endif
594     
595!        write (*,*) ' before advanced: ',j,prob(1),xmass1(j,1),ztra1(j)
596        call advance(itime,npoint(j),idt(j),uap(j),ucp(j),uzp(j), &
597             us(j),vs(j),ws(j),nstop,xtra1(j),ytra1(j),ztra1(j),prob, &
598             cbt(j))
599        if (ztra1(j).lt.30) then
600!        write (*,*) 'advanced: ',ztra1(j),j,xtra1(j),ytra1(j)
601        endif
602
603  ! Calculate the gross fluxes across layer interfaces
604  !***************************************************
605
606        if (iflux.eq.1) call calcfluxes(nage,j,xold,yold,zold)
607
608
609  ! Determine, when next time step is due
610  ! If trajectory is terminated, mark it
611  !**************************************
612
613        if (nstop.gt.1) then
614          if (linit_cond.ge.1) call initial_cond_calc(itime,j)
615          itra1(j)=-999999999
616        else
617          itra1(j)=itime+lsynctime
618
619
620  ! Dry deposition and radioactive decay for each species
621  ! Also check maximum (of all species) of initial mass remaining on the particle;
622  ! if it is below a threshold value, terminate particle
623  !*****************************************************************************
624
625          xmassfract=0.
626          do ks=1,nspec
627            if (decay(ks).gt.0.) then             ! radioactive decay
628              decfact=exp(-real(abs(lsynctime))*decay(ks))
629            else
630              decfact=1.
631            endif
632
633            if (DRYDEPSPEC(ks)) then        ! dry deposition
634              drydeposit(ks)=xmass1(j,ks)*prob(ks)*decfact
635              xmass1(j,ks)=xmass1(j,ks)*(1.-prob(ks))*decfact
636              if (decay(ks).gt.0.) then   ! correct for decay (see wetdepo)
637                drydeposit(ks)=drydeposit(ks)* &
638                     exp(real(abs(ldeltat))*decay(ks))
639              endif
640            else                           ! no dry deposition
641              xmass1(j,ks)=xmass1(j,ks)*decfact
642            endif
643
644! Skip check on mass fraction when npoint represents particle number
645            if (mdomainfill.eq.0.and.mquasilag.eq.0) then
646              if (xmass(npoint(j),ks).gt.0.) &
647                   xmassfract=max(xmassfract,real(npart(npoint(j)))* &
648                   xmass1(j,ks)/xmass(npoint(j),ks))
649!ZHG 2015
650                  !CGZ-lifetime: Check mass fraction left/save lifetime
651                   ! if(real(npart(npoint(j)))*xmass1(j,ks)/xmass(npoint(j),ks).lt.e_inv.and.checklifetime(j,ks).eq.0.)then
652                       !Mass below 1% of initial >register lifetime
653                       ! checklifetime(j,ks)=abs(itra1(j)-itramem(j))
654                       ! species_lifetime(ks,1)=species_lifetime(ks,1)+abs(itra1(j)-itramem(j))
655                       ! species_lifetime(ks,2)= species_lifetime(ks,2)+1
656                   ! endif
657                   !CGZ-lifetime: Check mass fraction left/save lifetime
658!ZHG 2015
659            else
660              xmassfract=1.0
661            end if
662          end do
663
664          if (xmassfract.lt.minmass) then   ! terminate all particles carrying less mass
665            itra1(j)=-999999999
666            if (verbosity.gt.0) then
667              print*,'terminated particle ',j,' for small mass'
668            endif
669          endif
670
671  !        Sabine Eckhardt, June 2008
672  !        don't create depofield for backward runs
673          if (DRYDEP.AND.(ldirect.eq.1)) then
674            call drydepokernel(nclass(j),drydeposit,real(xtra1(j)), &
675                 real(ytra1(j)),nage,kp)
676            if (nested_output.eq.1) call drydepokernel_nest( &
677                 nclass(j),drydeposit,real(xtra1(j)),real(ytra1(j)), &
678                 nage,kp)
679          endif
680
681  ! Terminate trajectories that are older than maximum allowed age
682  !***************************************************************
683
684          if (abs(itra1(j)-itramem(j)).ge.lage(nageclass)) then
685            if (linit_cond.ge.1) call initial_cond_calc(itime+lsynctime,j)
686            itra1(j)=-999999999
687            if (verbosity.gt.0) then
688              print*,'terminated particle ',j,' for age'
689            endif
690          endif
691        endif
692
693      endif
694
695    end do !loop over particles
696   
697  ! Counter of "unstable" particle velocity during a time scale of
698  ! maximumtl=20 minutes (defined in com_mod)
699  !***************************************************************
700   
701    total_nan_intl=0
702    i_nan=i_nan+1 ! added by mc to count nan during a time of maxtl (i.e. maximum tl fixed here to 20 minutes, see com_mod)
703    sum_nan_count(i_nan)=nan_count
704    if (i_nan > maxtl/lsynctime) i_nan=1 !lsynctime must be <= maxtl
705    do ii_nan=1, (maxtl/lsynctime)
706      total_nan_intl=total_nan_intl+sum_nan_count(ii_nan)
707    end do
708  ! Output to keep track of the numerical instabilities in CBL simulation and if
709  ! they are compromising the final result (or not)
710    if (cblflag.eq.1) print *,j,itime,'nan_synctime',nan_count,'nan_tl',total_nan_intl 
711         
712  end do
713
714
715  ! Complete the calculation of initial conditions for particles not yet terminated
716  !*****************************************************************************
717
718  do j=1,numpart
719    if (linit_cond.ge.1) call initial_cond_calc(itime,j)
720  end do
721
722  if (ipout.eq.2) call partoutput(itime)     ! dump particle positions
723
724  if (linit_cond.ge.1) call initial_cond_output(itime)   ! dump initial cond. field
725
726  !close(104)
727
728  ! De-allocate memory and end
729  !***************************
730
731  if (iflux.eq.1) then
732      deallocate(flux)
733  endif
734  if (OHREA) then
735      deallocate(OH_field,OH_hourly,lonOH,latOH,altOH)
736  endif
737  if (ldirect.gt.0) then
738  deallocate(drygridunc,wetgridunc)
739  endif
740  deallocate(gridunc)
741  deallocate(xpoint1,xpoint2,ypoint1,ypoint2,zpoint1,zpoint2,xmass)
742  deallocate(ireleasestart,ireleaseend,npart,kindz)
743  deallocate(xmasssave)
744  if (nested_output.eq.1) then
745     deallocate(orooutn, arean, volumen)
746     if (ldirect.gt.0) then
747     deallocate(griduncn,drygriduncn,wetgriduncn)
748     endif
749  endif
750  deallocate(outheight,outheighthalf)
751  deallocate(oroout, area, volume)
752
753end subroutine timemanager
754
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG