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

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

Fixed a bug where receptor arrays were not properly zeroed

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