source: branches/flexpart91_hasod/src/timemanager.f90 @ 7

Last change on this file since 7 was 7, checked in by hasod, 11 years ago

Initial import

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