source: flexpart.git/src/timemanager.f90 @ 2bec33e

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

resolved kernel merge: removed all uskernel

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