source: flexpart.git/src/timemanager.f90 @ 54cbd6c

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

all f90 files for dry/wet backward mode

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