source: flexpart.git/src/timemanager_mpi.f90 @ fdc0f03

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since fdc0f03 was fdc0f03, checked in by Espen Sollum ATMOS <eso@…>, 8 years ago

Code for cloud water should be correct if the total cw is stored in field clwc (old scheme) or in field qc (new scheme). Minor edits in some files.

  • Property mode set to 100644
File size: 32.7 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!   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! memstat            additional indicator for wind fields (see getfields)    *
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 mpi_mod
100  use netcdf_output_mod, only: concoutput_netcdf,concoutput_nest_netcdf,&
101       &concoutput_surf_netcdf,concoutput_surf_nest_netcdf
102
103  implicit none
104
105  logical :: reqv_state=.false. ! .true. if waiting for a MPI_Irecv to complete
106  integer :: j,ks,kp,l,n,itime=0,nstop,nstop1,memstat=0,mind
107! integer :: ksp
108  integer :: ip
109  integer :: loutnext,loutstart,loutend
110  integer :: ix,jy,ldeltat,itage,nage
111  integer :: i_nan=0,ii_nan,total_nan_intl=0  !added by mc to check instability in CBL scheme
112  integer :: numpart_tot_mpi ! for summing particles on all processes
113  real :: outnum,weight,prob(maxspec)
114  real :: decfact
115
116  real :: drydeposit(maxspec),gridtotalunc,wetgridtotalunc
117  real :: drygridtotalunc,xold,yold,zold,xmassfract
118  real, parameter :: e_inv = 1.0/exp(1.0)
119!double precision xm(maxspec,maxpointspec_act),
120!    +                 xm_depw(maxspec,maxpointspec_act),
121!    +                 xm_depd(maxspec,maxpointspec_act)
122
123!open(88,file='TEST.dat')
124
125! Measure time spent in timemanager
126  if (mp_measure_time) call mpif_mtime('timemanager',0)
127
128! First output for time 0
129!************************
130
131  loutnext=loutstep/2
132  outnum=0.
133  loutstart=loutnext-loutaver/2
134  loutend=loutnext+loutaver/2
135
136
137!**********************************************************************
138! Loop over the whole modelling period in time steps of mintime seconds
139!**********************************************************************
140
141
142!  itime=0
143  if (lroot) then
144  !  write(*,45) itime,numpart*mp_partgroup_np,gridtotalunc,wetgridtotalunc,drygridtotalunc
145    write(*,46) float(itime)/3600,itime,numpart*mp_partgroup_np
146   
147    if (verbosity.gt.0) then
148      write (*,*) 'timemanager> starting simulation'
149    end if
150  end if ! (lroot)
151
152!CGZ-lifetime: set lifetime to 0
153  ! checklifetime(:,:)=0
154  ! species_lifetime(:,:)=0
155  ! print*, 'Initialized lifetime'
156!CGZ-lifetime: set lifetime to 0
157
158
159
160  do itime=0,ideltas,lsynctime
161
162! Computation of wet deposition, OH reaction and mass transfer
163! between two species every lsynctime seconds
164! maybe wet depo frequency can be relaxed later but better be on safe side
165! wetdepo must be called BEFORE new fields are read in but should not
166! be called in the very beginning before any fields are loaded, or
167! before particles are in the system
168! Code may not be correct for decay of deposition
169! changed by Petra Seibert 9/02
170!********************************************************************
171
172    if (mp_dev_mode) write(*,*) 'myid, itime: ',mp_pid,itime
173   
174    if (wetdep .and. itime .ne. 0 .and. numpart .gt. 0) then
175      if (verbosity.gt.0) then
176        write (*,*) 'timemanager> call wetdepo'
177      endif
178
179      if (mp_measure_time) call mpif_mtime('wetdepo',0)
180       
181! readwind process skips this step
182      if (.not.(lmpreader.and.lmp_use_reader)) then
183        call wetdepo(itime,lsynctime,loutnext)
184      end if
185
186      if (mp_measure_time) call mpif_mtime('wetdepo',1)
187    end if
188
189    if (OHREA .and. itime .ne. 0 .and. numpart .gt. 0) then
190! readwind process skips this step
191      if (.not.(lmpreader.and.lmp_use_reader)) then
192        call ohreaction(itime,lsynctime,loutnext)
193      endif
194    end if
195
196    if (assspec .and. itime .ne. 0 .and. numpart .gt. 0) then
197      stop 'associated species not yet implemented!'
198!     call transferspec(itime,lsynctime,loutnext)
199    endif
200
201
202! compute convection for backward runs
203!*************************************
204
205    if ((ldirect.eq.-1).and.(lconvection.eq.1).and.(itime.lt.0)) then
206      if (verbosity.gt.0) then
207        write (*,*) 'timemanager> call convmix -- backward'
208      endif
209! readwind process skips this step
210      if (.not.(lmpreader.and.lmp_use_reader)) call convmix(itime)
211    endif
212
213! Get necessary wind fields if not available
214!*******************************************
215    if (verbosity.gt.0 .and. lmpreader) then
216      write (*,*) 'timemanager> call getfields'
217    endif
218
219! This time measure includes reading/MPI communication (for the reader process),
220! or MPI communication time only (for other processes)
221    if (mp_measure_time) call mpif_mtime('getfields',0)
222
223    call getfields(itime,nstop1,memstat)
224
225    if (mp_measure_time) call mpif_mtime('getfields',1)
226
227
228! Broadcast fields to all MPI processes
229! Skip if all processes have called getfields or if no new fields
230!*****************************************************************
231
232    if (mp_measure_time.and..not.(lmpreader.and.lmp_use_reader)) call mpif_mtime('getfields',0)
233
234! Version 1 (lmp_sync=.true.) uses a read-ahead process where send/recv is done
235! in sync at start of each new field time interval
236    if (lmp_sync.and.lmp_use_reader.and.memstat.gt.0) then
237      call mpif_gf_send_vars(memstat)
238      call mpif_gf_send_vars_nest(memstat)
239! Version 2  (lmp_sync=.false., see below) is also used whenever 2 new fields are
240! read (as at first time step), in which case async send/recv is impossible.
241    else if (.not.lmp_sync.and.lmp_use_reader.and.memstat.ge.32) then
242      call mpif_gf_send_vars(memstat)
243      call mpif_gf_send_vars_nest(memstat)
244    end if
245
246! Version 2 (lmp_sync=.false.) is for holding three fields in memory. Uses a
247! read-ahead process where sending/receiving of the 3rd fields is done in
248! the background in parallel with performing computations with fields 1&2
249!********************************************************************************
250    if (.not.lmp_sync) then
251   
252! READER PROCESS:
253      if (memstat.gt.0..and.memstat.lt.32.and.lmp_use_reader.and.lmpreader) then
254        if (mp_dev_mode) write(*,*) 'Reader process: calling mpif_gf_send_vars_async'
255        call mpif_gf_send_vars_async(memstat)
256      end if
257
258! COMPLETION CHECK:
259! Issued at start of each new field period.
260      if (memstat.ne.0.and.memstat.lt.32.and.lmp_use_reader) then
261! TODO: z0(7) changes with time, so should be dimension (numclass,2) to
262! allow transfer of the future value in the background
263        call MPI_Bcast(z0,numclass,mp_pp,id_read,MPI_COMM_WORLD,mp_ierr)
264        call mpif_gf_request
265      end if
266
267! RECVEIVING PROCESS(ES):
268      ! eso TODO: at this point we do not know if clwc/ciwc will be available
269      ! at next time step. Issue receive request anyway, cancel at mpif_gf_request
270      if (memstat.gt.0.and.lmp_use_reader.and..not.lmpreader) then
271        if (mp_dev_mode) write(*,*) 'Receiving process: calling mpif_gf_send_vars_async. PID: ', mp_pid
272        call mpif_gf_recv_vars_async(memstat)
273      end if
274
275    end if
276
277    if (mp_measure_time.and..not.(lmpreader.and.lmp_use_reader)) call mpif_mtime('getfields',1)
278
279! For validation and tests: call the function below to set all fields to simple
280! homogeneous values
281    if (mp_dev_mode) call set_fields_synthetic
282
283!*******************************************************************************
284
285    if (lmpreader.and.nstop1.gt.1) stop 'NO METEO FIELDS AVAILABLE'
286
287! Reader process goes back to top of time loop (unless simulation end)
288!*********************************************************************
289
290    if (lmpreader.and.lmp_use_reader) then
291      if (itime.lt.ideltas) then
292        cycle
293      else
294        goto 999
295      end if
296    end if
297
298
299! Get hourly OH fields if not available
300!****************************************************
301    if (OHREA) then
302      if (verbosity.gt.0) then
303        write (*,*) 'timemanager> call gethourlyOH'
304      endif
305      call gethourlyOH(itime)
306    endif
307
308
309! Release particles
310!******************
311
312    if (verbosity.gt.0.and.lroot) then
313      write (*,*) 'timemanager>  Release particles'
314    endif
315
316    if (mdomainfill.ge.1) then
317      if (itime.eq.0) then
318        if (verbosity.gt.0) then
319          write (*,*) 'timemanager>  call init_domainfill'
320        endif
321        call init_domainfill
322      else
323        if (verbosity.gt.0.and.lroot) then
324          write (*,*) 'timemanager>  call boundcond_domainfill'
325        endif
326        call boundcond_domainfill(itime,loutend)
327      endif
328    else
329      if (verbosity.gt.0.and.lroot) then
330        print*,'call releaseparticles' 
331      endif
332      call releaseparticles(itime)
333    endif
334
335
336! Compute convective mixing for forward runs
337! for backward runs it is done before next windfield is read in
338!**************************************************************
339
340    if ((ldirect.eq.1).and.(lconvection.eq.1)) then
341      if (verbosity.gt.0) then
342        write (*,*) 'timemanager> call convmix -- forward'
343      endif
344      call convmix(itime)
345    endif
346
347! If middle of averaging period of output fields is reached, accumulated
348! deposited mass radioactively decays
349!***********************************************************************
350
351    if (DEP.and.(itime.eq.loutnext).and.(ldirect.gt.0)) then
352      do ks=1,nspec
353        do kp=1,maxpointspec_act
354          if (decay(ks).gt.0.) then ! TODO move this statement up 2 levels
355            do nage=1,nageclass
356              do l=1,nclassunc
357! Mother output grid
358                do jy=0,numygrid-1
359                  do ix=0,numxgrid-1
360                    wetgridunc(ix,jy,ks,kp,l,nage)= &
361                         wetgridunc(ix,jy,ks,kp,l,nage)* &
362                         exp(-1.*outstep*decay(ks))
363                    drygridunc(ix,jy,ks,kp,l,nage)= &
364                         drygridunc(ix,jy,ks,kp,l,nage)* &
365                         exp(-1.*outstep*decay(ks))
366                  end do
367                end do
368! Nested output grid
369                if (nested_output.eq.1) then
370                  do jy=0,numygridn-1
371                    do ix=0,numxgridn-1
372                      wetgriduncn(ix,jy,ks,kp,l,nage)= &
373                           wetgriduncn(ix,jy,ks,kp,l,nage)* &
374                           exp(-1.*outstep*decay(ks))
375                      drygriduncn(ix,jy,ks,kp,l,nage)= &
376                           drygriduncn(ix,jy,ks,kp,l,nage)* &
377                           exp(-1.*outstep*decay(ks))
378                    end do
379                  end do
380                endif
381              end do
382            end do
383          endif
384        end do
385      end do
386    endif
387
388
389!!! CHANGE: These lines may be switched on to check the conservation
390!!! of mass within FLEXPART
391!   if (itime.eq.loutnext) then
392!   do 247 ksp=1, nspec
393!   do 247 kp=1, maxpointspec_act
394!47         xm(ksp,kp)=0.
395
396!   do 249 ksp=1, nspec
397!     do 249 j=1,numpart
398!          if (ioutputforeachrelease.eq.1) then
399!            kp=npoint(j)
400!          else
401!            kp=1
402!          endif
403!       if (itra1(j).eq.itime) then
404!          xm(ksp,kp)=xm(ksp,kp)+xmass1(j,ksp)
405!         write(*,*) 'xmass: ',xmass1(j,ksp),j,ksp,nspec
406!       endif
407!49     continue
408!  do 248 ksp=1,nspec
409!  do 248 kp=1,maxpointspec_act
410!  xm_depw(ksp,kp)=0.
411!  xm_depd(ksp,kp)=0.
412!     do 248 nage=1,nageclass
413!       do 248 ix=0,numxgrid-1
414!         do 248 jy=0,numygrid-1
415!           do 248 l=1,nclassunc
416!              xm_depw(ksp,kp)=xm_depw(ksp,kp)
417!    +                  +wetgridunc(ix,jy,ksp,kp,l,nage)
418!48                 xm_depd(ksp,kp)=xm_depd(ksp,kp)
419!    +                  +drygridunc(ix,jy,ksp,kp,l,nage)
420!             do 246 ksp=1,nspec
421!46                    write(88,'(2i10,3e12.3)')
422!    +              itime,ksp,(xm(ksp,kp),kp=1,maxpointspec_act),
423!    +                (xm_depw(ksp,kp),kp=1,maxpointspec_act),
424!    +                (xm_depd(ksp,kp),kp=1,maxpointspec_act)
425!  endif
426!!! CHANGE
427
428
429
430
431! Check whether concentrations are to be calculated
432!**************************************************
433
434    if ((ldirect*itime.ge.ldirect*loutstart).and. &
435         (ldirect*itime.le.ldirect*loutend)) then ! add to grid
436      if (mod(itime-loutstart,loutsample).eq.0) then
437
438! If we are exactly at the start or end of the concentration averaging interval,
439! give only half the weight to this sample
440!*****************************************************************************
441
442        if ((itime.eq.loutstart).or.(itime.eq.loutend)) then
443          weight=0.5
444        else
445          weight=1.0
446        endif
447        outnum=outnum+weight
448
449        call conccalc(itime,weight)
450
451      endif
452
453! :TODO: MPI output of particle positions;  each process sequentially
454!   access the same file
455      if ((mquasilag.eq.1).and.(itime.eq.(loutstart+loutend)/2)) &
456           call partoutput_short(itime)    ! dump particle positions in extremely compressed format
457
458
459! Output and reinitialization of grid
460! If necessary, first sample of new grid is also taken
461!*****************************************************
462
463      if ((itime.eq.loutend).and.(outnum.gt.0.)) then
464        if ((iout.le.3.).or.(iout.eq.5)) then
465
466! MPI: Root process collects/sums grids
467!**************************************
468          call mpif_tm_reduce_grid
469
470          if (mp_measure_time) call mpif_mtime('iotime',0)
471          if (surf_only.ne.1) then
472            if (lroot) then
473              if (lnetcdfout.eq.1) then
474                call concoutput_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,&
475                     &drygridtotalunc)
476              else
477                call concoutput(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
478              endif
479            else
480! zero arrays on non-root processes
481              gridunc(:,:,:,:,:,:,:)=0.
482              creceptor(:,:)=0.
483            end if
484          else
485            if (lroot) then
486              if (lnetcdfout.eq.1) then
487                call concoutput_surf_netcdf(itime,outnum,gridtotalunc,wetgridtotalunc,&
488                     &drygridtotalunc)
489              else
490                call concoutput_surf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
491              end if
492            else
493! zero arrays on non-root processes
494              gridunc(:,:,:,:,:,:,:)=0.
495              creceptor(:,:)=0.
496            endif
497          endif
498          if (mp_measure_time) call mpif_mtime('iotime',1)
499
500! :TODO: Correct calling of conc_surf above?
501
502!   call concoutput_surf(itime,outnum,gridtotalunc,wetgridtotalunc,drygridtotalunc)
503! endif
504
505          if (nested_output.eq.1) then
506
507! MPI: Root process collects/sums nested grids
508!*********************************************
509            call mpif_tm_reduce_grid_nest
510 
511           if (mp_measure_time) call mpif_mtime('iotime',0)
512
513            if (lnetcdfout.eq.0) then
514              if (surf_only.ne.1) then
515
516                if (lroot) then
517                  call concoutput_nest(itime,outnum)
518                else
519                  griduncn(:,:,:,:,:,:,:)=0.
520                end if
521
522              else  ! :TODO: check for zeroing in the netcdf module
523                call concoutput_surf_nest(itime,outnum)
524
525              end if
526
527            else
528
529              if (surf_only.ne.1) then
530                if (lroot) then             
531                  call concoutput_nest_netcdf(itime,outnum)
532                else
533                  griduncn(:,:,:,:,:,:,:)=0.
534                end if
535              else
536                if (lroot) then
537                  call concoutput_surf_nest_netcdf(itime,outnum)
538                else
539                  griduncn(:,:,:,:,:,:,:)=0.
540                end if
541              endif
542
543
544            end if
545          end if
546         
547
548          outnum=0.
549        endif
550        if ((iout.eq.4).or.(iout.eq.5)) call plumetraj(itime)
551        if (iflux.eq.1) call fluxoutput(itime)
552        if (mp_measure_time) call mpif_mtime('iotime',1)
553
554
555! Decide whether to write an estimate of the number of particles released,
556! or exact number (require MPI reduce operation)
557        numpart_tot_mpi = numpart*mp_partgroup_np
558
559        if (mp_exact_numpart.and..not.(lmpreader.and.lmp_use_reader)) then
560          call MPI_Reduce(numpart, numpart_tot_mpi, 1, MPI_INTEGER, MPI_SUM, id_root, &
561               & mp_comm_used, mp_ierr)
562        endif
563       
564        !CGZ-lifetime: output species lifetime
565        if (lroot) then
566        !   write(*,*) 'Overview species lifetime in days', &
567        !        real((species_lifetime(:,1)/species_lifetime(:,2))/real(3600.0*24.0))
568        !   write(*,*) 'all info:',species_lifetime
569          write(*,45) itime,numpart_tot_mpi,gridtotalunc,&
570               &wetgridtotalunc,drygridtotalunc
571        !   if (verbosity.gt.0) then
572        !     write (*,*) 'timemanager> starting simulation'
573        !   end if
574        end if
575
57645      format(i13,' SECONDS SIMULATED: ',i13, ' PARTICLES:    Uncertainty: ',3f7.3)
57746      format(' Simulated ',f7.1,' hours (',i13,' s), ',i13, ' particles')
578        if (ipout.ge.1) then
579          do ip=0, mp_partgroup_np-1
580            if (ip.eq.mp_partid) call partoutput(itime) ! dump particle positions
581            call mpif_mpi_barrier
582          end do
583        end if
584
585        loutnext=loutnext+loutstep
586        loutstart=loutnext-loutaver/2
587        loutend=loutnext+loutaver/2
588        if (itime.eq.loutstart) then
589          weight=0.5
590          outnum=outnum+weight
591          call conccalc(itime,weight)
592        endif
593
594
595! Check, whether particles are to be split:
596! If so, create new particles and attribute all information from the old
597! particles also to the new ones; old and new particles both get half the
598! mass of the old ones
599!************************************************************************
600       
601        if (ldirect*itime.ge.ldirect*itsplit) then
602          n=numpart
603          do j=1,numpart
604            if (ldirect*itime.ge.ldirect*itrasplit(j)) then
605!                if (n.lt.maxpart) then
606              if (n.lt.maxpart_mpi) then
607                n=n+1
608                itrasplit(j)=2*(itrasplit(j)-itramem(j))+itramem(j)
609                itrasplit(n)=itrasplit(j)
610                itramem(n)=itramem(j)
611                itra1(n)=itra1(j)
612                idt(n)=idt(j)
613                npoint(n)=npoint(j)
614                nclass(n)=nclass(j)
615                xtra1(n)=xtra1(j)
616                ytra1(n)=ytra1(j)
617                ztra1(n)=ztra1(j)
618                uap(n)=uap(j)
619                ucp(n)=ucp(j)
620                uzp(n)=uzp(j)
621                us(n)=us(j)
622                vs(n)=vs(j)
623                ws(n)=ws(j)
624                cbt(n)=cbt(j)
625                do ks=1,nspec
626                  xmass1(j,ks)=xmass1(j,ks)/2.
627                  xmass1(n,ks)=xmass1(j,ks)
628                end do
629              endif
630            endif
631          end do
632          numpart=n
633        endif
634      endif
635    endif
636   
637
638    if (itime.eq.ideltas) exit         ! almost finished
639
640! Compute interval since radioactive decay of deposited mass was computed
641!************************************************************************
642
643    if (itime.lt.loutnext) then
644      ldeltat=itime-(loutnext-loutstep)
645    else                                  ! first half of next interval
646      ldeltat=itime-loutnext
647    endif
648
649
650! Loop over all particles
651!************************
652    if (mp_measure_time) call mpif_mtime('partloop1',0)
653
654!--------------------------------------------------------------------------------
655! various variables for testing reason of CBL scheme, by mc
656    well_mixed_vector=0. !erase vector to test well mixed condition: modified by mc
657    well_mixed_norm=0.   !erase normalization to test well mixed condition: modified by mc
658    avg_ol=0.
659    avg_wst=0.
660    avg_h=0.
661    avg_air_dens=0.  !erase vector to obtain air density at particle positions: modified by mc
662!--------------------------------------------------------------------------------
663
664    do j=1,numpart
665
666! If integration step is due, do it
667!**********************************
668
669      if (itra1(j).eq.itime) then
670
671        if (ioutputforeachrelease.eq.1) then
672          kp=npoint(j)
673        else
674          kp=1
675        endif
676! Determine age class of the particle
677        itage=abs(itra1(j)-itramem(j))
678        do nage=1,nageclass
679          if (itage.lt.lage(nage)) exit
680        end do
681
682! Initialize newly released particle
683!***********************************
684
685        if ((itramem(j).eq.itime).or.(itime.eq.0)) &
686             call initialize(itime,idt(j),uap(j),ucp(j),uzp(j), &
687             us(j),vs(j),ws(j),xtra1(j),ytra1(j),ztra1(j),cbt(j))
688
689! Memorize particle positions
690!****************************
691
692        xold=xtra1(j)
693        yold=ytra1(j)
694        zold=ztra1(j)
695
696! Integrate Lagevin equation for lsynctime seconds
697!*************************************************
698
699        mp_advance_wtime_beg = mpi_wtime()
700
701        call advance(itime,npoint(j),idt(j),uap(j),ucp(j),uzp(j), &
702             us(j),vs(j),ws(j),nstop,xtra1(j),ytra1(j),ztra1(j),prob, &
703             cbt(j))
704
705        mp_advance_wtime_end = mpi_wtime()
706        mp_advance_wtime_total = mp_advance_wtime_total + (mp_advance_wtime_end - &
707             & mp_advance_wtime_beg)
708
709
710! Calculate the gross fluxes across layer interfaces
711!***************************************************
712
713        if (iflux.eq.1) call calcfluxes(nage,j,xold,yold,zold)
714
715
716! Determine, when next time step is due
717! If trajectory is terminated, mark it
718!**************************************
719
720        if (nstop.gt.1) then
721          if (linit_cond.ge.1) call initial_cond_calc(itime,j)
722          itra1(j)=-999999999
723        else
724          itra1(j)=itime+lsynctime
725
726
727! Dry deposition and radioactive decay for each species
728! Also check maximum (of all species) of initial mass remaining on the particle;
729! if it is below a threshold value, terminate particle
730!*****************************************************************************
731
732          xmassfract=0.
733          do ks=1,nspec
734            if (decay(ks).gt.0.) then             ! radioactive decay
735              decfact=exp(-real(abs(lsynctime))*decay(ks))
736            else
737              decfact=1.
738            endif
739
740            if (DRYDEPSPEC(ks)) then        ! dry deposition
741              drydeposit(ks)=xmass1(j,ks)*prob(ks)*decfact
742              xmass1(j,ks)=xmass1(j,ks)*(1.-prob(ks))*decfact
743              if (decay(ks).gt.0.) then   ! correct for decay (see wetdepo)
744                drydeposit(ks)=drydeposit(ks)* &
745                     exp(real(abs(ldeltat))*decay(ks))
746              endif
747            else                           ! no dry deposition
748              xmass1(j,ks)=xmass1(j,ks)*decfact
749            endif
750
751
752            if (mdomainfill.eq.0) then
753              if (xmass(npoint(j),ks).gt.0.)then
754                   xmassfract=max(xmassfract,real(npart(npoint(j)))* &
755                   xmass1(j,ks)/xmass(npoint(j),ks))
756                   
757                   !CGZ-lifetime: Check mass fraction left/save lifetime
758                   ! if(lroot.and.real(npart(npoint(j)))*xmass1(j,ks)/xmass(npoint(j),ks).lt.inv_e.and.checklifetime(j,ks).eq.0.)then
759                       !Mass below 1% of initial >register lifetime
760                   !     checklifetime(j,ks)=abs(itra1(j)-itramem(j))
761
762                   !     species_lifetime(ks,1)=species_lifetime(ks,1)+abs(itra1(j)-itramem(j))
763                   !     species_lifetime(ks,2)= species_lifetime(ks,2)+1
764                   ! endif
765                   !CGZ-lifetime: Check mass fraction left/save lifetime
766                   
767              endif
768            else
769              xmassfract=1.
770            endif
771          end do
772
773          if (xmassfract.lt.minmass) then ! .and. sum(real(npart(npoint(j)))*xmass1(j,:)).lt.1.0) then   ! terminate all particles carrying less mass
774          !            print*,'terminated particle ',j,' for small mass (', sum(real(npart(npoint(j)))* &
775          !         xmass1(j,:)), ' of ', sum(xmass(npoint(j),:)),')'
776            itra1(j)=-999999999
777          endif
778
779!        Sabine Eckhardt, June 2008
780!        don't create depofield for backward runs
781          if (DRYDEP.AND.(ldirect.eq.1)) then
782            call drydepokernel(nclass(j),drydeposit,real(xtra1(j)), &
783                 real(ytra1(j)),nage,kp)
784            if (nested_output.eq.1) call drydepokernel_nest( &
785                 nclass(j),drydeposit,real(xtra1(j)),real(ytra1(j)), &
786                 nage,kp)
787          endif
788
789! Terminate trajectories that are older than maximum allowed age
790!***************************************************************
791
792          if (abs(itra1(j)-itramem(j)).ge.lage(nageclass)) then
793            if (linit_cond.ge.1) &
794                 call initial_cond_calc(itime+lsynctime,j)
795            itra1(j)=-999999999
796            !print*, 'terminated particle ',j,'for age'
797          endif
798        endif
799
800      endif
801
802    end do ! j=1, numpart
803
804    if(mp_measure_time) call mpif_mtime('partloop1',1)
805
806
807! Added by mc: counter of "unstable" particle velocity during a time scale
808!   of maximumtl=20 minutes (defined in com_mod)
809
810    total_nan_intl=0
811    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)
812    sum_nan_count(i_nan)=nan_count
813    if (i_nan > maxtl/lsynctime) i_nan=1 !lsynctime must be <= maxtl
814    do ii_nan=1, (maxtl/lsynctime)
815      total_nan_intl=total_nan_intl+sum_nan_count(ii_nan)
816    end do
817
818! Output to keep track of the numerical instabilities in CBL simulation
819! and if they are compromising the final result (or not):
820    if (cblflag.eq.1) print *,j,itime,'nan_synctime',nan_count,'nan_tl',total_nan_intl 
821
822
823  end do ! itime=0,ideltas,lsynctime
824
825
826! Complete the calculation of initial conditions for particles not yet terminated
827!*****************************************************************************
828
829! eso :TODO: this not implemented yet (transfer particles to PID 0 or rewrite)
830! the tools to do this are already in mpi_mod.f90
831! :CUR: check this
832  if (lroot) then
833    do j=1,numpart
834      if (linit_cond.ge.1) call initial_cond_calc(itime,j)
835    end do
836  end if
837
838
839  if (ipout.eq.2) then
840! MPI process 0 creates the file, the other processes append to it
841    do ip=0, mp_partgroup_np-1
842      if (ip.eq.mp_partid) then
843        !if (mp_dbg_mode) write(*,*) 'call partoutput(itime), proc, mp_partid',ip,mp_partid
844        call partoutput(itime)    ! dump particle positions
845      end if
846      call mpif_mpi_barrier
847    end do
848  end if
849
850! eso :TODO: MPI
851  if (linit_cond.ge.1.and.lroot) call initial_cond_output(itime)   ! dump initial cond. field
852
853
854! De-allocate memory and end
855!***************************
856
857999 if (iflux.eq.1) then
858    deallocate(flux)
859  endif
860  if (OHREA) then
861    deallocate(OH_field,OH_hourly,lonOH,latOH,altOH)
862  endif
863  if (ldirect.gt.0) then
864    deallocate(drygridunc,wetgridunc)
865  endif
866  deallocate(gridunc)
867  deallocate(xpoint1,xpoint2,ypoint1,ypoint2,zpoint1,zpoint2,xmass, checklifetime)
868  deallocate(ireleasestart,ireleaseend,npart,kindz)
869  deallocate(xmasssave)
870  if (nested_output.eq.1) then
871    deallocate(orooutn, arean, volumen)
872    if (ldirect.gt.0) then
873      deallocate(griduncn,drygriduncn,wetgriduncn)
874    endif
875  endif
876  deallocate(outheight,outheighthalf)
877  deallocate(oroout, area, volume)
878
879  if (mp_measure_time) call mpif_mtime('timemanager',1)
880
881end subroutine timemanager
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG