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

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

Made enabling netCDF output during compilation optional

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