source: flexpart.git/src/timemanager_mpi.f90 @ 78e62dc

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 78e62dc was 78e62dc, checked in by flexpart <>, 9 years ago

New OH parameter in SPECIES files (now 3 instead of 2). New path to OH binariy files.

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