source: flexpart.git/src/timemanager_mpi.f90 @ 16b61a5

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

Reworked the domain-filling option (MPI). Fixed a slow loop which had errors in loop counter (MPI)

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