source: flexpart.git/flexpart_code/timemanager.F90 @ b398fb6

FPv9.3.1FPv9.3.1b_testingFPv9.3.2fp9.3.1-20161214-nc4grib2nc4_repair
Last change on this file since b398fb6 was 496c607, checked in by Don Morton <Don.Morton@…>, 8 years ago

Initial commit of FPv9.3.1

Currently, this is a clone of snapshot FPv9.3.0

  • Property mode set to 100644
File size: 21.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(metdata_format)
23
24  !*****************************************************************************
25  !                                                                            *
26  ! Handles the computation of trajectories, i.e. determines which             *
27  ! trajectories have to be computed at what time.                             *
28  ! Manages dry+wet deposition routines, radioactive decay and the computation *
29  ! of concentrations.                                                         *
30  !                                                                            *
31  !     Author: A. Stohl                                                       *
32  !                                                                            *
33  !     20 May 1996                                                            *
34  !                                                                            *
35  !*****************************************************************************
36  !  Changes, Bernd C. Krueger, Feb. 2001:                                     *
37  !        Call of convmix when new windfield is read                          *
38  !------------------------------------                                        *
39  !  Changes Petra Seibert, Sept 2002                                          *
40  !     fix wet scavenging problem                                             *
41  !     Code may not be correct for decay of deposition!                       *
42  !  Changes Petra Seibert, Nov 2002                                           *
43  !     call convection BEFORE new fields are read in BWD mode                 *
44  !  Changes Caroline Forster, Feb 2005                                        *
45  !new interface between flexpart and convection scheme                        *
46  !Emanuel's latest subroutine convect43c.f is used                            *
47  !*****************************************************************************
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,metdata_format
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  !write (*,*) 'starting simulation'
131  do itime=0,ideltas,lsynctime
132
133
134  ! Computation of wet deposition, OH reaction and mass transfer
135  ! between two species every lsynctime seconds
136  ! maybe wet depo frequency can be relaxed later but better be on safe side
137  ! wetdepo must be called BEFORE new fields are read in but should not
138  ! be called in the very beginning before any fields are loaded, or
139  ! before particles are in the system
140  ! Code may not be correct for decay of deposition
141  ! changed by Petra Seibert 9/02
142  !********************************************************************
143
144    if (WETDEP .and. itime .ne. 0 .and. numpart .gt. 0) &
145         call wetdepo(itime,lsynctime,loutnext)
146
147    if (OHREA .and. itime .ne. 0 .and. numpart .gt. 0) &
148         call ohreaction(itime,lsynctime,loutnext)
149
150    if (ASSSPEC .and. itime .ne. 0 .and. numpart .gt. 0) then
151       stop 'associated species not yet implemented!'
152  !     call transferspec(itime,lsynctime,loutnext)
153    endif
154
155  ! compute convection for backward runs
156  !*************************************
157
158    if (metdata_format.eq.ecmwf_metdata) then
159      if ((ldirect.eq.-1).and.(lconvection.eq.1).and.(itime.lt.0)) &
160           call convmix_ecmwf(itime)
161    endif
162    if (metdata_format.eq.gfs_metdata) then
163      if ((ldirect.eq.-1).and.(lconvection.eq.1).and.(itime.lt.0)) &
164           call convmix_gfs(itime)
165    endif
166
167
168  ! Get necessary wind fields if not available
169  !*******************************************
170
171#ifdef TESTUSEGETFPFIELDS
172    ! Set up
173    call getfpfields(itime,nstop1,metdata_format)
174#else
175    if ( preprocessed_metdata.eq.1 ) then
176      call getfpfields(itime,nstop1,metdata_format)
177    else
178      call getfields(itime,nstop1,metdata_format)
179    endif
180#endif
181    if (nstop1.gt.1) stop 'NO METEO FIELDS AVAILABLE'
182  ! Release particles
183  !******************
184
185    if (mdomainfill.ge.1) then
186      if (itime.eq.0) then
187        call init_domainfill
188      else
189        call boundcond_domainfill(itime,loutend)
190      endif
191    else
192      call releaseparticles(itime)
193    endif
194
195
196  ! Compute convective mixing for forward runs
197  ! for backward runs it is done before next windfield is read in
198  !**************************************************************
199    if (metdata_format.eq.ecmwf_metdata) then
200      if ((ldirect.eq.1).and.(lconvection.eq.1)) &
201           call convmix_ecmwf(itime)
202    endif
203    if (metdata_format.eq.gfs_metdata) then
204      if ((ldirect.eq.1).and.(lconvection.eq.1)) &
205           call convmix_gfs(itime)
206    endif
207
208
209  ! If middle of averaging period of output fields is reached, accumulated
210  ! deposited mass radioactively decays
211  !***********************************************************************
212
213    if (DEP.and.(itime.eq.loutnext).and.(ldirect.gt.0)) then
214      do ks=1,nspec
215      do kp=1,maxpointspec_act
216        if (decay(ks).gt.0.) then
217          do nage=1,nageclass
218            do l=1,nclassunc
219  ! Mother output grid
220              do jy=0,numygrid-1
221                do ix=0,numxgrid-1
222                  wetgridunc(ix,jy,ks,kp,l,nage)= &
223                       wetgridunc(ix,jy,ks,kp,l,nage)* &
224                       exp(-1.*outstep*decay(ks))
225                  drygridunc(ix,jy,ks,kp,l,nage)= &
226                       drygridunc(ix,jy,ks,kp,l,nage)* &
227                       exp(-1.*outstep*decay(ks))
228                end do
229              end do
230  ! Nested output grid
231              if (nested_output.eq.1) then
232                do jy=0,numygridn-1
233                  do ix=0,numxgridn-1
234                    wetgriduncn(ix,jy,ks,kp,l,nage)= &
235                         wetgriduncn(ix,jy,ks,kp,l,nage)* &
236                         exp(-1.*outstep*decay(ks))
237                    drygriduncn(ix,jy,ks,kp,l,nage)= &
238                         drygriduncn(ix,jy,ks,kp,l,nage)* &
239                         exp(-1.*outstep*decay(ks))
240                  end do
241                end do
242              endif
243            end do
244          end do
245        endif
246      end do
247      end do
248    endif
249
250  !!! CHANGE: These lines may be switched on to check the conservation
251  !!! of mass within FLEXPART
252  !   if (itime.eq.loutnext) then
253  !   do 247 ksp=1, nspec
254  !   do 247 kp=1, maxpointspec_act
255  !47         xm(ksp,kp)=0.
256
257  !   do 249 ksp=1, nspec
258  !     do 249 j=1,numpart
259  !          if (ioutputforeachrelease.eq.1) then
260  !            kp=npoint(j)
261  !          else
262  !            kp=1
263  !          endif
264  !       if (itra1(j).eq.itime) then
265  !          xm(ksp,kp)=xm(ksp,kp)+xmass1(j,ksp)
266  !         write(*,*) 'xmass: ',xmass1(j,ksp),j,ksp,nspec
267  !       endif
268  !49     continue
269  !  do 248 ksp=1,nspec
270  !  do 248 kp=1,maxpointspec_act
271  !  xm_depw(ksp,kp)=0.
272  !  xm_depd(ksp,kp)=0.
273  !     do 248 nage=1,nageclass
274  !       do 248 ix=0,numxgrid-1
275  !         do 248 jy=0,numygrid-1
276  !           do 248 l=1,nclassunc
277  !              xm_depw(ksp,kp)=xm_depw(ksp,kp)
278  !    +                  +wetgridunc(ix,jy,ksp,kp,l,nage)
279  !48                 xm_depd(ksp,kp)=xm_depd(ksp,kp)
280  !    +                  +drygridunc(ix,jy,ksp,kp,l,nage)
281  !             do 246 ksp=1,nspec
282  !46                    write(88,'(2i10,3e12.3)')
283  !    +              itime,ksp,(xm(ksp,kp),kp=1,maxpointspec_act),
284  !    +                (xm_depw(ksp,kp),kp=1,maxpointspec_act),
285  !    +                (xm_depd(ksp,kp),kp=1,maxpointspec_act)
286  !  endif
287  !!! CHANGE
288
289
290
291  ! Check whether concentrations are to be calculated
292  !**************************************************
293
294    if ((ldirect*itime.ge.ldirect*loutstart).and. &
295         (ldirect*itime.le.ldirect*loutend)) then ! add to grid
296      if (mod(itime-loutstart,loutsample).eq.0) then
297
298  ! If we are exactly at the start or end of the concentration averaging interval,
299  ! give only half the weight to this sample
300  !*****************************************************************************
301
302        if ((itime.eq.loutstart).or.(itime.eq.loutend)) then
303          weight=0.5
304        else
305          weight=1.0
306        endif
307        outnum=outnum+weight
308        call conccalc(itime,weight)
309      endif
310
311
312      if ((mquasilag.eq.1).and.(itime.eq.(loutstart+loutend)/2)) &
313           call partoutput_short(itime)    ! dump particle positions in extremely compressed format
314
315
316  ! Output and reinitialization of grid
317  ! If necessary, first sample of new grid is also taken
318  !*****************************************************
319
320      if ((itime.eq.loutend).and.(outnum.gt.0.)) then
321        if ((iout.le.3.).or.(iout.eq.5)) then
322          call concoutput(itime,outnum,gridtotalunc, &
323               wetgridtotalunc,drygridtotalunc)
324          if (nested_output.eq.1) call concoutput_nest(itime,outnum)
325          outnum=0.
326        endif
327        if ((iout.eq.4).or.(iout.eq.5)) call plumetraj(itime)
328        if (iflux.eq.1) call fluxoutput(itime)
329        write(*,45) itime,numpart,gridtotalunc,wetgridtotalunc, &
330             drygridtotalunc
33145      format(i9,' SECONDS SIMULATED: ',i8, &
332             ' PARTICLES:    Uncertainty: ',3f7.3)
333        if (ipout.ge.1) call partoutput(itime)    ! dump particle positions
334        loutnext=loutnext+loutstep
335        loutstart=loutnext-loutaver/2
336        loutend=loutnext+loutaver/2
337        if (itime.eq.loutstart) then
338          weight=0.5
339          outnum=outnum+weight
340          call conccalc(itime,weight)
341        endif
342
343
344  ! Check, whether particles are to be split:
345  ! If so, create new particles and attribute all information from the old
346  ! particles also to the new ones; old and new particles both get half the
347  ! mass of the old ones
348  !************************************************************************
349
350        if (ldirect*itime.ge.ldirect*itsplit) then
351          n=numpart
352          do j=1,numpart
353            if (ldirect*itime.ge.ldirect*itrasplit(j)) then
354              if (n.lt.maxpart) then
355                n=n+1
356                itrasplit(j)=2*(itrasplit(j)-itramem(j))+itramem(j)
357                itrasplit(n)=itrasplit(j)
358                itramem(n)=itramem(j)
359                itra1(n)=itra1(j)
360                idt(n)=idt(j)
361                npoint(n)=npoint(j)
362                nclass(n)=nclass(j)
363                xtra1(n)=xtra1(j)
364                ytra1(n)=ytra1(j)
365                ztra1(n)=ztra1(j)
366                uap(n)=uap(j)
367                ucp(n)=ucp(j)
368                uzp(n)=uzp(j)
369                us(n)=us(j)
370                vs(n)=vs(j)
371                ws(n)=ws(j)
372                cbt(n)=cbt(j)
373                do ks=1,nspec
374                  xmass1(j,ks)=xmass1(j,ks)/2.
375                  xmass1(n,ks)=xmass1(j,ks)
376                end do
377              endif
378            endif
379          end do
380          numpart=n
381        endif
382      endif
383    endif
384
385
386    if (itime.eq.ideltas) exit         ! almost finished
387
388  ! Compute interval since radioactive decay of deposited mass was computed
389  !************************************************************************
390
391    if (itime.lt.loutnext) then
392      ldeltat=itime-(loutnext-loutstep)
393    else                                  ! first half of next interval
394      ldeltat=itime-loutnext
395    endif
396
397
398  ! Loop over all particles
399  !************************
400
401    do j=1,numpart
402
403
404  ! If integration step is due, do it
405  !**********************************
406
407      if (itra1(j).eq.itime) then
408
409        if (ioutputforeachrelease.eq.1) then
410            kp=npoint(j)
411        else
412            kp=1
413        endif
414  ! Determine age class of the particle
415        itage=abs(itra1(j)-itramem(j))
416        do nage=1,nageclass
417          if (itage.lt.lage(nage)) exit
418        end do
419
420  ! Initialize newly released particle
421  !***********************************
422
423        if ((itramem(j).eq.itime).or.(itime.eq.0)) &
424             call initialize(itime,idt(j),uap(j),ucp(j),uzp(j), &
425             us(j),vs(j),ws(j),xtra1(j),ytra1(j),ztra1(j),cbt(j))
426
427  ! Memorize particle positions
428  !****************************
429
430        xold=xtra1(j)
431        yold=ytra1(j)
432        zold=ztra1(j)
433
434  ! Integrate Lagevin equation for lsynctime seconds
435  !*************************************************
436
437        call advance(itime,npoint(j),idt(j),uap(j),ucp(j),uzp(j), &
438             us(j),vs(j),ws(j),nstop,xtra1(j),ytra1(j),ztra1(j),prob, &
439             cbt(j))
440
441  ! Calculate the gross fluxes across layer interfaces
442  !***************************************************
443
444        if (iflux.eq.1) call calcfluxes(nage,j,xold,yold,zold)
445
446
447  ! Determine, when next time step is due
448  ! If trajectory is terminated, mark it
449  !**************************************
450
451        if (nstop.gt.1) then
452          if (linit_cond.ge.1) call initial_cond_calc(itime,j)
453          itra1(j)=-999999999
454        else
455          itra1(j)=itime+lsynctime
456
457
458  ! Dry deposition and radioactive decay for each species
459  ! Also check maximum (of all species) of initial mass remaining on the particle;
460  ! if it is below a threshold value, terminate particle
461  !*****************************************************************************
462
463          xmassfract=0.
464          do ks=1,nspec
465            if (decay(ks).gt.0.) then             ! radioactive decay
466              decfact=exp(-real(abs(lsynctime))*decay(ks))
467            else
468              decfact=1.
469            endif
470
471            if (DRYDEPSPEC(ks)) then        ! dry deposition
472              drydeposit(ks)=xmass1(j,ks)*prob(ks)*decfact
473              xmass1(j,ks)=xmass1(j,ks)*(1.-prob(ks))*decfact
474              if (decay(ks).gt.0.) then   ! correct for decay (see wetdepo)
475                drydeposit(ks)=drydeposit(ks)* &
476                     exp(real(abs(ldeltat))*decay(ks))
477              endif
478            else                           ! no dry deposition
479              xmass1(j,ks)=xmass1(j,ks)*decfact
480            endif
481
482
483            if (mdomainfill.eq.0) then
484              if (xmass(npoint(j),ks).gt.0.) &
485                   xmassfract=max(xmassfract,real(npart(npoint(j)))* &
486                   xmass1(j,ks)/xmass(npoint(j),ks))
487            else
488              xmassfract=1.
489            endif
490          end do
491
492          if (xmassfract.lt.0.0001) then   ! terminate all particles carrying less mass
493            itra1(j)=-999999999
494          endif
495
496  !        Sabine Eckhardt, June 2008
497  !        don't create depofield for backward runs
498          if (DRYDEP.AND.(ldirect.eq.1)) then
499            call drydepokernel(nclass(j),drydeposit,real(xtra1(j)), &
500                 real(ytra1(j)),nage,kp)
501            if (nested_output.eq.1) call drydepokernel_nest( &
502                 nclass(j),drydeposit,real(xtra1(j)),real(ytra1(j)), &
503                 nage,kp)
504          endif
505
506  ! Terminate trajectories that are older than maximum allowed age
507  !***************************************************************
508
509          if (abs(itra1(j)-itramem(j)).ge.lage(nageclass)) then
510            if (linit_cond.ge.1) &
511                 call initial_cond_calc(itime+lsynctime,j)
512            itra1(j)=-999999999
513          endif
514        endif
515
516      endif
517
518    end do
519
520  end do
521
522
523  ! Complete the calculation of initial conditions for particles not yet terminated
524  !*****************************************************************************
525
526  do j=1,numpart
527    if (linit_cond.ge.1) call initial_cond_calc(itime,j)
528  end do
529
530  if (ipout.eq.2) call partoutput(itime)     ! dump particle positions
531
532  if (linit_cond.ge.1) call initial_cond_output(itime)   ! dump initial cond. field
533
534  close(104)
535
536  ! De-allocate memory and end
537  !***************************
538
539  if (iflux.eq.1) then
540      deallocate(flux)
541  endif
542  if (OHREA.eqv..TRUE.) then
543      deallocate(OH_field,OH_field_height)
544  endif
545  if (ldirect.gt.0) then
546  deallocate(drygridunc,wetgridunc)
547  endif
548  deallocate(gridunc)
549  deallocate(xpoint1,xpoint2,ypoint1,ypoint2,zpoint1,zpoint2,xmass)
550  deallocate(ireleasestart,ireleaseend,npart,kindz)
551  deallocate(xmasssave)
552  if (nested_output.eq.1) then
553     deallocate(orooutn, arean, volumen)
554     if (ldirect.gt.0) then
555     deallocate(griduncn,drygriduncn,wetgriduncn)
556     endif
557  endif
558  deallocate(outheight,outheighthalf)
559  deallocate(oroout, area, volume)
560
561end subroutine timemanager
562
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG