source: flexpart.git/src/timemanager.f90 @ b7ae015

10.4.1_peseiFPv9.3.1FPv9.3.1b_testingFPv9.3.2GFS_025bugfixes+enhancementsdevfp9.3.1-20161214-nc4grib2nc4_repairrelease-10release-10.4.1scaling-bugunivie
Last change on this file since b7ae015 was b7ae015, checked in by Ignacio Pisso <Ignacio.Pisso@…>, 9 years ago

clean up runtime messages and adjust verbosity level

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