source: flexpart.git/src/timemanager.f90 @ 5f53d0e

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

writing information if kernel is used

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