source: flexpart.git/src/timemanager_mpi.f90 @ 18adf60

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

Bug found in the MQUASILAG=1 option, xmass array accessed with out-of-bounds index. Quick fix implemented (see changes in timemanager and advance).

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