source: flexpart.git/src/timemanager_mpi.f90 @ 0a94e13

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

Added ipout=3 option for time averaged particle output

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