source: flexpart.git/src/timemanager.f90 @ 02095e3

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

Renamed lnokernel, corrected default setting.

  • Property mode set to 100644
File size: 29.4 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(metdata_format)
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  !  Unified ECMWF and GFS builds                                              *
53  !   Marian Harustak, 12.5.2017                                               *
54  !   - Added passing of metdata_format as it was needed by called routines    *
55  !*****************************************************************************
56  !                                                                            *
57  ! Variables:                                                                 *
58  ! DEP                .true. if either wet or dry deposition is switched on   *
59  ! decay(maxspec) [1/s] decay constant for radioactive decay                  *
60  ! DRYDEP             .true. if dry deposition is switched on                 *
61  ! ideltas [s]        modelling period                                        *
62  ! itime [s]          actual temporal position of calculation                 *
63  ! ldeltat [s]        time since computation of radioact. decay of depositions*
64  ! loutaver [s]       averaging period for concentration calculations         *
65  ! loutend [s]        end of averaging for concentration calculations         *
66  ! loutnext [s]       next time at which output fields shall be centered      *
67  ! loutsample [s]     sampling interval for averaging of concentrations       *
68  ! loutstart [s]      start of averaging for concentration calculations       *
69  ! loutstep [s]       time interval for which concentrations shall be         *
70  !                    calculated                                              *
71  ! npoint(maxpart)    index, which starting point the trajectory has          *
72  !                    starting positions of trajectories                      *
73  ! nstop              serves as indicator for fate of particles               *
74  !                    in the particle loop                                    *
75  ! nstop1             serves as indicator for wind fields (see getfields)     *
76  ! outnum             number of samples for each concentration calculation    *
77  ! outnum             number of samples for each concentration calculation    *
78  ! prob               probability of absorption at ground due to dry          *
79  !                    deposition                                              *
80  ! WETDEP             .true. if wet deposition is switched on                 *
81  ! weight             weight for each concentration sample (1/2 or 1)         *
82  ! uap(maxpart),ucp(maxpart),uzp(maxpart) = random velocities due to          *
83  !                    turbulence                                              *
84  ! us(maxpart),vs(maxpart),ws(maxpart) = random velocities due to inter-      *
85  !                    polation                                                *
86  ! xtra1(maxpart), ytra1(maxpart), ztra1(maxpart) =                           *
87  !                    spatial positions of trajectories                       *
88  ! metdata_format     format of metdata (ecmwf/gfs)                           *
89  !                                                                            *
90  ! Constants:                                                                 *
91  ! maxpart            maximum number of trajectories                          *
92  !                                                                            *
93  !*****************************************************************************
94
95  use unc_mod
96  use point_mod
97  use xmass_mod
98  use flux_mod
99  use outg_mod
100  use oh_mod
101  use par_mod
102  use com_mod
103  use netcdf_output_mod, only: concoutput_netcdf,concoutput_nest_netcdf,&
104       &concoutput_surf_netcdf,concoutput_surf_nest_netcdf
105
106  implicit none
107
108  integer :: metdata_format
109  integer :: j,ks,kp,l,n,itime=0,nstop,nstop1
110! integer :: ksp
111  integer :: loutnext,loutstart,loutend
112  integer :: ix,jy,ldeltat,itage,nage,idummy
113  integer :: i_nan=0,ii_nan,total_nan_intl=0  !added by mc to check instability in CBL scheme
114  real :: outnum,weight,prob_rec(maxspec),prob(maxspec),decfact,wetscav
115  ! real :: uap(maxpart),ucp(maxpart),uzp(maxpart)
116  ! real :: us(maxpart),vs(maxpart),ws(maxpart)
117  ! integer(kind=2) :: cbt(maxpart)
118  real(sp) :: gridtotalunc
119  real(dep_prec) :: drydeposit(maxspec),wetgridtotalunc,drygridtotalunc
120  real :: xold,yold,zold,xmassfract
121  real :: grfraction(3)
122  real, parameter :: e_inv = 1.0/exp(1.0)
123
124  !double precision xm(maxspec,maxpointspec_act),
125  !    +                 xm_depw(maxspec,maxpointspec_act),
126  !    +                 xm_depd(maxspec,maxpointspec_act)
127
128
129  !open(88,file='TEST.dat')
130
131  ! First output for time 0
132  !************************
133
134  loutnext=loutstep/2
135  outnum=0.
136  loutstart=loutnext-loutaver/2
137  loutend=loutnext+loutaver/2
138
139  !  open(127,file=path(2)(1:length(2))//'depostat.dat'
140  !    +  ,form='unformatted')
141  !write (*,*) 'writing deposition statistics depostat.dat!'
142
143  !**********************************************************************
144  ! Loop over the whole modelling period in time steps of mintime seconds
145  !**********************************************************************
146
147!ZHG 2015
148!CGZ-lifetime: set lifetime to 0
149  ! checklifetime(:,:)=0
150  ! species_lifetime(:,:)=0
151  ! print*, 'Initialized lifetime'
152!CGZ-lifetime: set lifetime to 0
153 
154  if (.not.lusekerneloutput) write(*,*) 'Not using the kernel'
155  if (turboff) write(*,*) 'Turbulence switched off'
156
157  write(*,46) float(itime)/3600,itime,numpart
158
159  if (verbosity.gt.0) then
160    write (*,*) 'timemanager> starting simulation'
161    if (verbosity.gt.1) then
162      CALL SYSTEM_CLOCK(count_clock)
163      WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
164    endif     
165  endif
166
167  do itime=0,ideltas,lsynctime
168
169  ! Computation of wet deposition, OH reaction and mass transfer
170  ! between two species every lsynctime seconds
171  ! maybe wet depo frequency can be relaxed later but better be on safe side
172  ! wetdepo must be called BEFORE new fields are read in but should not
173  ! be called in the very beginning before any fields are loaded, or
174  ! before particles are in the system
175  ! Code may not be correct for decay of deposition
176  ! changed by Petra Seibert 9/02
177  !********************************************************************
178
179    if (WETDEP .and. itime .ne. 0 .and. numpart .gt. 0) then
180        if (verbosity.gt.0) then
181           write (*,*) 'timemanager> call wetdepo'
182        endif     
183         call wetdepo(itime,lsynctime,loutnext)
184    endif
185
186    if (OHREA .and. itime .ne. 0 .and. numpart .gt. 0) &
187         call ohreaction(itime,lsynctime,loutnext)
188
189    if (ASSSPEC .and. itime .ne. 0 .and. numpart .gt. 0) then
190       stop 'associated species not yet implemented!'
191  !     call transferspec(itime,lsynctime,loutnext)
192    endif
193
194  ! compute convection for backward runs
195  !*************************************
196
197   if ((ldirect.eq.-1).and.(lconvection.eq.1).and.(itime.lt.0)) then
198        if (verbosity.gt.0) then
199           write (*,*) 'timemanager> call convmix -- backward'
200        endif         
201      call convmix(itime,metdata_format)
202        if (verbosity.gt.1) then
203          !CALL SYSTEM_CLOCK(count_clock, count_rate, count_max)
204          CALL SYSTEM_CLOCK(count_clock)
205          WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
206        endif
207   endif
208
209  ! Get necessary wind fields if not available
210  !*******************************************
211    if (verbosity.gt.0) then
212           write (*,*) 'timemanager> call getfields'
213    endif
214    call getfields(itime,nstop1,metdata_format)
215        if (verbosity.gt.1) then
216          CALL SYSTEM_CLOCK(count_clock)
217          WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
218        endif
219    if (nstop1.gt.1) stop 'NO METEO FIELDS AVAILABLE'
220
221  ! Get hourly OH fields if not available
222  !****************************************************
223    if (OHREA) then
224      if (verbosity.gt.0) then
225             write (*,*) 'timemanager> call gethourlyOH'
226      endif
227      call gethourlyOH(itime)
228          if (verbosity.gt.1) then
229            CALL SYSTEM_CLOCK(count_clock)
230            WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
231          endif
232    endif
233       
234  ! Release particles
235  !******************
236
237    if (verbosity.gt.0) then
238           write (*,*) 'timemanager>  Release particles'
239    endif
240
241    if (mdomainfill.ge.1) then
242      if (itime.eq.0) then
243        if (verbosity.gt.0) then
244          write (*,*) 'timemanager>  call init_domainfill'
245        endif       
246        call init_domainfill
247      else
248        if (verbosity.gt.0) then
249          write (*,*) 'timemanager>  call boundcond_domainfill'
250        endif   
251        call boundcond_domainfill(itime,loutend)
252      endif
253    else
254      if (verbosity.gt.0) then
255        print*,'call releaseparticles' 
256      endif
257      call releaseparticles(itime)
258      if (verbosity.gt.1) then
259        CALL SYSTEM_CLOCK(count_clock)
260        WRITE(*,*) 'timemanager> SYSTEM CLOCK',(count_clock - count_clock0)/real(count_rate)
261      endif
262    endif
263
264
265  ! Compute convective mixing for forward runs
266  ! for backward runs it is done before next windfield is read in
267  !**************************************************************
268
269   if ((ldirect.eq.1).and.(lconvection.eq.1)) then
270     if (verbosity.gt.0) then
271       write (*,*) 'timemanager> call convmix -- forward'
272     endif   
273     call convmix(itime,metdata_format)
274   endif
275
276  ! If middle of averaging period of output fields is reached, accumulated
277  ! deposited mass radioactively decays
278  !***********************************************************************
279
280    if (DEP.and.(itime.eq.loutnext).and.(ldirect.gt.0)) then
281      do ks=1,nspec
282      do kp=1,maxpointspec_act
283        if (decay(ks).gt.0.) then
284          do nage=1,nageclass
285            do l=1,nclassunc
286  ! Mother output grid
287              do jy=0,numygrid-1
288                do ix=0,numxgrid-1
289                  wetgridunc(ix,jy,ks,kp,l,nage)= &
290                       wetgridunc(ix,jy,ks,kp,l,nage)* &
291                       exp(-1.*outstep*decay(ks))
292                  drygridunc(ix,jy,ks,kp,l,nage)= &
293                       drygridunc(ix,jy,ks,kp,l,nage)* &
294                       exp(-1.*outstep*decay(ks))
295                end do
296              end do
297  ! Nested output grid
298              if (nested_output.eq.1) then
299                do jy=0,numygridn-1
300                  do ix=0,numxgridn-1
301                    wetgriduncn(ix,jy,ks,kp,l,nage)= &
302                         wetgriduncn(ix,jy,ks,kp,l,nage)* &
303                         exp(-1.*outstep*decay(ks))
304                    drygriduncn(ix,jy,ks,kp,l,nage)= &
305                         drygriduncn(ix,jy,ks,kp,l,nage)* &
306                         exp(-1.*outstep*decay(ks))
307                  end do
308                end do
309              endif
310            end do
311          end do
312        endif
313      end do
314      end do
315    endif
316
317  !!! CHANGE: These lines may be switched on to check the conservation
318  !!! of mass within FLEXPART
319  !   if (itime.eq.loutnext) then
320  !   do 247 ksp=1, nspec
321  !   do 247 kp=1, maxpointspec_act
322  !47         xm(ksp,kp)=0.
323
324  !   do 249 ksp=1, nspec
325  !     do 249 j=1,numpart
326  !          if (ioutputforeachrelease.eq.1) then
327  !            kp=npoint(j)
328  !          else
329  !            kp=1
330  !          endif
331  !       if (itra1(j).eq.itime) then
332  !          xm(ksp,kp)=xm(ksp,kp)+xmass1(j,ksp)
333  !         write(*,*) 'xmass: ',xmass1(j,ksp),j,ksp,nspec
334  !       endif
335  !49     continue
336  !  do 248 ksp=1,nspec
337  !  do 248 kp=1,maxpointspec_act
338  !  xm_depw(ksp,kp)=0.
339  !  xm_depd(ksp,kp)=0.
340  !     do 248 nage=1,nageclass
341  !       do 248 ix=0,numxgrid-1
342  !         do 248 jy=0,numygrid-1
343  !           do 248 l=1,nclassunc
344  !              xm_depw(ksp,kp)=xm_depw(ksp,kp)
345  !    +                  +wetgridunc(ix,jy,ksp,kp,l,nage)
346  !48                 xm_depd(ksp,kp)=xm_depd(ksp,kp)
347  !    +                  +drygridunc(ix,jy,ksp,kp,l,nage)
348  !             do 246 ksp=1,nspec
349  !46                    write(88,'(2i10,3e12.3)')
350  !    +              itime,ksp,(xm(ksp,kp),kp=1,maxpointspec_act),
351  !    +                (xm_depw(ksp,kp),kp=1,maxpointspec_act),
352  !    +                (xm_depd(ksp,kp),kp=1,maxpointspec_act)
353  !  endif
354  !!! CHANGE
355
356
357
358  ! Check whether concentrations are to be calculated
359  !**************************************************
360
361    if ((ldirect*itime.ge.ldirect*loutstart).and. &
362         (ldirect*itime.le.ldirect*loutend)) then ! add to grid
363      if (mod(itime-loutstart,loutsample).eq.0) then
364
365  ! If we are exactly at the start or end of the concentration averaging interval,
366  ! give only half the weight to this sample
367  !*****************************************************************************
368
369        if ((itime.eq.loutstart).or.(itime.eq.loutend)) then
370          weight=0.5
371        else
372          weight=1.0
373        endif
374        outnum=outnum+weight
375        call conccalc(itime,weight)
376      endif
377
378
379      if ((mquasilag.eq.1).and.(itime.eq.(loutstart+loutend)/2)) &
380           call partoutput_short(itime)    ! dump particle positions in extremely compressed format
381
382
383  ! Output and reinitialization of grid
384  ! If necessary, first sample of new grid is also taken
385  !*****************************************************
386
387      if ((itime.eq.loutend).and.(outnum.gt.0.)) then
388        if ((iout.le.3.).or.(iout.eq.5)) then
389          if (surf_only.ne.1) then
390            if (lnetcdfout.eq.1) then
391              call concoutput_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
392            else
393              call concoutput(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
394            endif
395          else 
396            if (verbosity.eq.1) then
397             print*,'call concoutput_surf '
398             call system_clock(count_clock)
399             write(*,*) 'system clock',count_clock - count_clock0   
400            endif
401            if (lnetcdfout.eq.1) then
402              call concoutput_surf_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
403            else
404              call concoutput_surf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
405              if (verbosity.eq.1) then
406                print*,'called concoutput_surf '
407                call system_clock(count_clock)
408                write(*,*) 'system clock',count_clock - count_clock0   
409              endif
410            endif
411          endif
412
413          if (nested_output .eq. 1) then
414            if (lnetcdfout.eq.0) then
415              if (surf_only.ne.1) then
416                call concoutput_nest(itime,outnum)
417              else
418                call concoutput_surf_nest(itime,outnum)
419              endif
420            else
421              if (surf_only.ne.1) then
422                call concoutput_nest_netcdf(itime,outnum)
423              else
424                call concoutput_surf_nest_netcdf(itime,outnum)
425              endif
426            endif
427          endif
428          outnum=0.
429        endif
430        if ((iout.eq.4).or.(iout.eq.5)) call plumetraj(itime)
431        if (iflux.eq.1) call fluxoutput(itime)
432        write(*,45) itime,numpart,gridtotalunc,wetgridtotalunc,drygridtotalunc
433 
434        !CGZ-lifetime: output species lifetime
435!ZHG
436        ! write(*,*) 'Overview species lifetime in days', &
437        !      real((species_lifetime(:,1)/species_lifetime(:,2))/real(3600.0*24.0))
438        ! write(*,*) 'all info:',species_lifetime
439!ZHG
440        !CGZ-lifetime: output species lifetime
441
442        !write(*,46) float(itime)/3600,itime,numpart
44345      format(i13,' SECONDS SIMULATED: ',i13, ' PARTICLES:    Uncertainty: ',3f7.3)
44446      format(' Simulated ',f7.1,' hours (',i13,' s), ',i13, ' particles')
445        if (ipout.ge.1) call partoutput(itime)    ! dump particle positions
446        loutnext=loutnext+loutstep
447        loutstart=loutnext-loutaver/2
448        loutend=loutnext+loutaver/2
449        if (itime.eq.loutstart) then
450          weight=0.5
451          outnum=outnum+weight
452          call conccalc(itime,weight)
453        endif
454
455
456  ! Check, whether particles are to be split:
457  ! If so, create new particles and attribute all information from the old
458  ! particles also to the new ones; old and new particles both get half the
459  ! mass of the old ones
460  !************************************************************************
461
462        if (ldirect*itime.ge.ldirect*itsplit) then
463          n=numpart
464          do j=1,numpart
465            if (ldirect*itime.ge.ldirect*itrasplit(j)) then
466              if (n.lt.maxpart) then
467                n=n+1
468                itrasplit(j)=2*(itrasplit(j)-itramem(j))+itramem(j)
469                itrasplit(n)=itrasplit(j)
470                itramem(n)=itramem(j)
471                itra1(n)=itra1(j)
472                idt(n)=idt(j)
473                npoint(n)=npoint(j)
474                nclass(n)=nclass(j)
475                xtra1(n)=xtra1(j)
476                ytra1(n)=ytra1(j)
477                ztra1(n)=ztra1(j)
478                uap(n)=uap(j)
479                ucp(n)=ucp(j)
480                uzp(n)=uzp(j)
481                us(n)=us(j)
482                vs(n)=vs(j)
483                ws(n)=ws(j)
484                cbt(n)=cbt(j)
485                do ks=1,nspec
486                  xmass1(j,ks)=xmass1(j,ks)/2.
487                  xmass1(n,ks)=xmass1(j,ks)
488                end do
489              endif
490            endif
491          end do
492          numpart=n
493        endif
494      endif
495    endif
496
497
498    if (itime.eq.ideltas) exit         ! almost finished
499
500  ! Compute interval since radioactive decay of deposited mass was computed
501  !************************************************************************
502
503    if (itime.lt.loutnext) then
504      ldeltat=itime-(loutnext-loutstep)
505    else                                  ! first half of next interval
506      ldeltat=itime-loutnext
507    endif
508
509
510  ! Loop over all particles
511  !************************
512  ! Various variables for testing reason of CBL scheme, by mc
513    well_mixed_vector=0. !erase vector to test well mixed condition: modified by mc
514    well_mixed_norm=0.   !erase normalization to test well mixed condition: modified by mc
515    avg_ol=0.
516    avg_wst=0.
517    avg_h=0.
518    avg_air_dens=0.  !erase vector to obtain air density at particle positions: modified by mc
519  !-----------------------------------------------------------------------------
520    do j=1,numpart
521
522
523  ! If integration step is due, do it
524  !**********************************
525
526      if (itra1(j).eq.itime) then
527
528        if (ioutputforeachrelease.eq.1) then
529            kp=npoint(j)
530        else
531            kp=1
532        endif
533  ! Determine age class of the particle
534        itage=abs(itra1(j)-itramem(j))
535        do nage=1,nageclass
536          if (itage.lt.lage(nage)) exit
537        end do
538
539  ! Initialize newly released particle
540  !***********************************
541
542        if ((itramem(j).eq.itime).or.(itime.eq.0)) &
543             call initialize(itime,idt(j),uap(j),ucp(j),uzp(j), &
544             us(j),vs(j),ws(j),xtra1(j),ytra1(j),ztra1(j),cbt(j))
545
546  ! Memorize particle positions
547  !****************************
548
549        xold=xtra1(j)
550        yold=ytra1(j)
551        zold=ztra1(j)
552
553   
554  ! RECEPTOR: dry/wet depovel
555  !****************************
556  ! Before the particle is moved
557  ! the calculation of the scavenged mass shall only be done once after release
558  ! xscav_frac1 was initialised with a negative value
559
560      if  (DRYBKDEP) then
561       do ks=1,nspec
562         if  ((xscav_frac1(j,ks).lt.0)) then
563            call get_vdep_prob(itime,xtra1(j),ytra1(j),ztra1(j),prob_rec)
564            if (DRYDEPSPEC(ks)) then        ! dry deposition
565               xscav_frac1(j,ks)=prob_rec(ks)
566             else
567                xmass1(j,ks)=0.
568                xscav_frac1(j,ks)=0.
569             endif
570         endif
571        enddo
572       endif
573
574       if (WETBKDEP) then
575       do ks=1,nspec
576         if  ((xscav_frac1(j,ks).lt.0)) then
577            call get_wetscav(itime,lsynctime,loutnext,j,ks,grfraction,idummy,idummy,wetscav)
578            if (wetscav.gt.0) then
579                xscav_frac1(j,ks)=wetscav* &
580                       (zpoint2(npoint(j))-zpoint1(npoint(j)))*grfraction(1)
581            else
582                xmass1(j,ks)=0.
583                xscav_frac1(j,ks)=0.
584            endif
585         endif
586        enddo
587       endif
588
589  ! Integrate Lagevin equation for lsynctime seconds
590  !*************************************************
591
592        if (verbosity.gt.0) then
593           if (j.eq.1) then
594             write (*,*) 'timemanager> call advance'
595           endif     
596        endif
597     
598        call advance(itime,npoint(j),idt(j),uap(j),ucp(j),uzp(j), &
599             us(j),vs(j),ws(j),nstop,xtra1(j),ytra1(j),ztra1(j),prob, &
600             cbt(j))
601!        write (*,*) 'advance: ',prob(1),xmass1(j,1),ztra1(j)
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