source: flexpart.git/src/mpi_mod.f90 @ 62e65c7

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

'incloud_ratio' set to 6.2. 'maxnests' reverted to 0. Updated info messages (readspecies)

  • Property mode set to 100644
File size: 80.8 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
22module mpi_mod
23
24!*****************************************************************************
25!                                                                            *
26!  DESCRIPTION                                                               *
27!    This module contains subroutines and common variables used for the      *
28!    MPI parallelization of FLEXPART.                                        *
29!                                                                            *
30!  NOTE                                                                      *
31!    Depending on the MPI library installed on your system (e.g. mpich2,     *
32!    OpenMPI) you may need to choose below in this file between              *
33!      use mpi                                                               *
34!    (if the MPI library comes with the file 'mpi.mod'); or                  *
35!      include 'mpif.h'                                                      *
36!                                                                            *
37!                                                                            *
38!*****************************************************************************
39!                                                                            *
40! Variables:                                                                 *
41!                                                                            *
42! mp_ierr                 MPI error code                                     *
43! mp_np                   Number of MPI processes                            *
44! mp_pid                  Process ID of each MPI process                     *
45! mp_seed                 Parameter for random number seed                   *
46! read_grp_min            Minimum number of processes at which one will be   *
47!                         used as reader                                     *
48! numpart_mpi,            Number of particles per node                       *
49! maxpart_mpi                                                                *
50! mp_partid               MPI process ID for particle calculation            *
51! mp_partgroup_           Refers to the subset of processors performing      *
52!                         loops over particles. Will be all processes        *
53!                         unless a dedicated process runs getfields/readwind *
54! lmp_sync                If .false., use asynchronous MPI                   *
55! mp_cp                   Real precision to use for deposition fields        *
56!                                                                            *
57!                                                                            *
58!                                                                            *
59!                                                                            *
60!*****************************************************************************
61
62  use mpi
63  use par_mod, only: dp,sp
64  use com_mod, only: lroot
65
66  implicit none
67
68!  include 'mpif.h'
69
70  public
71
72! Set aside a process for reading windfields if using at least these many processes
73!==================================================
74  integer, parameter, private :: read_grp_min=4
75!==================================================
76
77! Variables for each MPI process in the world group
78  integer :: mp_ierr, mp_np, mp_pid, mp_partid
79  integer, private :: world_group_id
80
81! Variables for MPI processes in the 'particle' group
82  integer, allocatable, dimension(:) :: mp_partgroup_rank
83  integer :: mp_partgroup_comm, mp_partgroup_pid, mp_partgroup_np
84
85  integer :: mp_seed=0
86  integer, parameter :: mp_sp=MPI_REAL4, mp_dp=MPI_REAL8
87  integer :: mp_cp
88  integer, parameter :: id_root=0 ! master process
89
90! MPI tags/requests for send/receive operation
91  integer :: tm1
92  integer, parameter :: nvar_async=26 !27 !29 :DBG:
93!integer, dimension(:), allocatable :: tags
94  integer, dimension(:), allocatable :: reqs
95
96
97  integer :: id_read   ! readwind/getfield process
98  integer :: numpart_mpi,maxpart_mpi ! number of particles per node
99  integer :: tot_numpart=0
100  integer :: mp_comm_used ! global or subgroup communicator
101
102  logical :: lmpreader=.false. ! is set to true for reading process(es) only.
103  logical :: lmp_use_reader=.false. ! true if separate readwind process is used
104
105! .true. if only using synchronous MPI send/recv (default)
106! If setting this to .false., numwfmem must be set to 3
107!===============================================================================
108  logical :: lmp_sync=.true.
109!===============================================================================
110
111! mp_dbg_mode       Used for debugging MPI.
112! mp_dev_mode       various checks related to debugging the parallel code
113! mp_dbg_out        write some arrays to data file for debugging
114! mp_measure_time   Measure cpu/wall time, write out at end of run
115! mp_time_barrier   Measure MPI barrier time
116! mp_exact_numpart  Use an extra MPI communication to give the exact number of particles
117!                   to standard output (this does *not* otherwise affect the simulation)
118  logical, parameter :: mp_dbg_mode = .false.
119  logical, parameter :: mp_dev_mode = .false.
120  logical, parameter :: mp_dbg_out = .false.
121  logical, parameter :: mp_time_barrier=.true.
122  logical, parameter :: mp_measure_time=.false.
123  logical, parameter :: mp_exact_numpart=.true.
124
125! for measuring CPU/Wall time
126  real(sp),private :: mp_comm_time_beg, mp_comm_time_end, mp_comm_time_total=0.
127  real(dp),private :: mp_comm_wtime_beg, mp_comm_wtime_end, mp_comm_wtime_total=0.
128  real(sp),private :: mp_root_time_beg, mp_root_time_end, mp_root_time_total=0.
129  real(dp),private :: mp_root_wtime_beg, mp_root_wtime_end, mp_root_wtime_total=0.
130  real(sp),private :: mp_barrier_time_beg, mp_barrier_time_end, mp_barrier_time_total=0.
131  real(dp),private :: mp_barrier_wtime_beg, mp_barrier_wtime_end, mp_barrier_wtime_total=0.
132  real(sp),private :: tm_nploop_beg, tm_nploop_end, tm_nploop_total=0.
133  real(sp),private :: tm_tot_beg, tm_tot_end, tm_tot_total=0.
134  real(dp),private :: mp_getfields_wtime_beg, mp_getfields_wtime_end, mp_getfields_wtime_total=0.
135  real(sp),private :: mp_getfields_time_beg, mp_getfields_time_end, mp_getfields_time_total=0.
136  real(dp),private :: mp_readwind_wtime_beg, mp_readwind_wtime_end, mp_readwind_wtime_total=0.
137  real(sp),private :: mp_readwind_time_beg, mp_readwind_time_end, mp_readwind_time_total=0.
138  real(dp),private :: mp_io_wtime_beg, mp_io_wtime_end, mp_io_wtime_total=0.
139  real(sp),private :: mp_io_time_beg, mp_io_time_end, mp_io_time_total=0.
140  real(dp),private :: mp_wetdepo_wtime_beg, mp_wetdepo_wtime_end, mp_wetdepo_wtime_total=0.
141  real(sp),private :: mp_wetdepo_time_beg, mp_wetdepo_time_end, mp_wetdepo_time_total=0.
142  real(dp),private :: mp_advance_wtime_beg, mp_advance_wtime_end, mp_advance_wtime_total=0.
143  real(dp),private :: mp_conccalc_time_beg, mp_conccalc_time_end, mp_conccalc_time_total=0.
144  real(dp),private :: mp_total_wtime_beg, mp_total_wtime_end, mp_total_wtime_total=0.
145  real(dp),private :: mp_vt_wtime_beg, mp_vt_wtime_end, mp_vt_wtime_total
146  real(sp),private :: mp_vt_time_beg, mp_vt_time_end, mp_vt_time_total
147
148! dat_lun           logical unit number for i/o
149  integer, private :: dat_lun
150
151contains
152
153  subroutine mpif_init
154!***********************************************************************
155! DESCRIPTION
156!   Initialize MPI.
157!   
158!   Create the global communicator MPI_COMM_WORLD
159!   If using a separate MPI process for getfields/readwind, a subgroup
160!   is created for the other processes.
161!
162! VARIABLES
163!   mpi_mode    default 0, set to 2/3 if running MPI version
164!   mp_np       number of running processes, decided at run-time
165!***********************************************************************
166    use par_mod, only: maxpart, numwfmem, dep_prec
167    use com_mod, only: mpi_mode, verbosity
168
169    implicit none
170
171    integer :: i,j,s,addmaxpart=0
172
173! Each process gets an ID (mp_pid) in the range 0,..,mp_np-1
174    call MPI_INIT(mp_ierr)
175    if (mp_ierr /= 0) goto 100
176    call MPI_COMM_RANK(MPI_COMM_WORLD, mp_pid, mp_ierr)
177    if (mp_ierr /= 0) goto 100
178    call MPI_COMM_SIZE(MPI_COMM_WORLD, mp_np, mp_ierr)
179    if (mp_ierr /= 0) goto 100
180
181
182! Variable mpi_mode is used to handle subroutines common to parallel/serial version
183    if (lmp_sync) then
184      mpi_mode=2 ! hold 2 windfields in memory
185    else
186      mpi_mode=3 ! hold 3 windfields in memory
187    end if
188
189    if (mp_pid.ne.0) then
190      lroot = .false.
191    end if
192
193! Set MPI precision to use for transferring deposition fields
194!************************************************************
195    if (dep_prec==dp) then
196      mp_cp = MPI_REAL8
197      ! TODO: write info message for serial version as well
198      if (lroot.and.verbosity>0) write(*,*) 'Using double precision for deposition fields'
199    else if (dep_prec==sp) then
200      mp_cp = MPI_REAL4
201      if (lroot.and.verbosity>0) write(*,*) 'Using single precision for deposition fields'
202    else
203      write(*,*) 'ERROR: something went wrong setting MPI real precision'
204      stop
205    end if
206
207! Check for sensible combination of parameters
208!*********************************************
209
210    if (.not.lmp_sync.and.numwfmem.ne.3) then
211      if (lroot) then
212        write(*,FMT='(80("#"))')
213        write(*,*) '#### mpi_mod::mpif_init> ERROR: ', &
214             & 'numwfmem must be set to 3 for asyncronous reading ####'
215        write(*,FMT='(80("#"))')
216      end if
217      call MPI_FINALIZE(mp_ierr)
218      stop
219    else if (lmp_sync.and.numwfmem.ne.2.and.lroot) then
220      write(*,FMT='(80("#"))')
221      write(*,*) '#### mpi_mod::mpif_init> WARNING: ', &
222           & 'numwfmem should be set to 2 for syncronous'
223      write(*,*) ' reading. Results will still be valid, but unneccesary memory &
224           &is allocated.'
225      write(*,FMT='(80("#"))')
226! Force "syncronized" version if all processes will call getfields
227    else if ((.not.lmp_sync.and.mp_np.lt.read_grp_min).or.(mp_np.eq.1)) then
228      if (lroot) then
229        write(*,FMT='(80("#"))')
230        write(*,*) '#### mpi_mod::mpif_init> WARNING: ', &
231             & 'all procs call getfields. Setting lmp_sync=.true.'
232        write(*,FMT='(80("#"))')
233      end if
234      lmp_sync=.true. ! :DBG: eso fix this...
235    end if
236
237! TODO: Add more warnings for unimplemented flexpart features
238
239! Set ID of process that calls getfield/readwind.
240! Using the last process in the group ensures statistical identical results
241! as running with one process less but not using separate read process
242!**********************************************************************
243
244!    id_read = min(mp_np-1, 1)
245    id_read = mp_np-1
246
247    if (mp_pid.eq.id_read) lmpreader=.true.
248
249    call MPI_Comm_group (MPI_COMM_WORLD, world_group_id, mp_ierr)
250
251! Create a MPI group/communicator that will calculate trajectories.
252! Skip this step if program is run with only a few processes
253!************************************************************************
254    allocate(mp_partgroup_rank(0:mp_np-2))
255
256! This allows checking for allocation of arrays when no subgroup is used
257    mp_partgroup_pid=0
258
259    if (read_grp_min.lt.2) then
260      write(*,*) '#### mpi_mod::mpif_init> ERROR ####', &
261           & 'read_grp_min should be at least 2. Exiting'
262      stop
263    end if
264
265    if (mp_np.ge.read_grp_min) then
266      lmp_use_reader = .true.
267
268! Map the subgroup IDs= 0,1,2,...,mp_np-2, skipping reader process
269      j=-1
270      do i=0, mp_np-2 !loop over all (subgroup) IDs
271        j=j+1
272        if (i.eq.id_read) then
273          j=j+1
274        end if
275        mp_partgroup_rank(i) = j
276      end do
277
278      call MPI_Group_incl (world_group_id, mp_np-1, mp_partgroup_rank, &
279           &mp_partgroup_pid, mp_ierr)
280      if (mp_ierr /= 0) goto 100
281      call MPI_Comm_create (MPI_COMM_WORLD, mp_partgroup_pid, mp_partgroup_comm, mp_ierr)
282      if (mp_ierr /= 0) goto 100
283
284      if (mp_pid.ne.id_read) then
285        call MPI_Comm_rank (mp_partgroup_comm, mp_partgroup_pid, mp_ierr)
286        if (mp_ierr /= 0) goto 100
287
288! Get size of new group (=mp_np-1)
289        call MPI_COMM_SIZE(mp_partgroup_comm, mp_partgroup_np, mp_ierr)
290        if (mp_ierr /= 0) goto 100
291        if (mp_partgroup_np.ne.mp_np-1) then
292          write(*,*)  '#### mpi_mod:: mpif_init> ERROR ####&
293               & mp_partgroup_np.ne.mp_np-1'
294          stop
295        endif
296
297      else
298        mp_partgroup_pid = -1
299      end if
300    end if
301
302    if (lmp_use_reader) then
303      mp_comm_used = mp_partgroup_comm
304      mp_partid = mp_partgroup_pid
305      mp_partgroup_np=mp_np-1
306    else
307      mp_comm_used = MPI_COMM_WORLD
308      mp_partgroup_np = mp_np
309      mp_partid = mp_pid
310    end if
311
312! Set maxpart per process
313    if (mp_partid.lt.mod(maxpart,mp_partgroup_np)) addmaxpart=1
314    maxpart_mpi=int(maxpart/mp_partgroup_np)+addmaxpart
315
316! Do not allocate particle data arrays for readwind process
317    if (lmpreader.and.lmp_use_reader) then
318      maxpart_mpi=0
319    end if
320
321! Set random seed for each non-root process
322    if (mp_pid.gt.0) then
323!    if (mp_pid.ge.0) then
324!      call system_clock(s)
325      s = 244
326      mp_seed = -abs(mod((s*181)*((mp_pid-83)*359), 104729))
327    end if
328    if (mp_dev_mode) write(*,*) 'PID, mp_seed: ',mp_pid, mp_seed
329    if (mp_dbg_mode) then
330! :DBG: For debugging, set all seed to 0 and maxrand to e.g. 20
331      mp_seed=0
332      if (lroot) write(*,*) 'MPI: setting seed=0'
333    end if
334
335! Allocate request array for asynchronous MPI
336    if (.not.lmp_sync) then
337      allocate(reqs(0:nvar_async*mp_np-1))
338      reqs(:)=MPI_REQUEST_NULL
339    else
340      allocate(reqs(0:1))
341      reqs(:)=MPI_REQUEST_NULL
342    end if
343
344    goto 101
345
346100 write(*,*) '#### mpi_mod::mpif_init> ERROR ####', mp_ierr
347    stop
348
349101 end subroutine mpif_init
350
351
352  subroutine mpif_mpi_barrier
353!***********************************************************************
354! Set MPI barrier (all processes are synchronized here), with the
355! possible exception of a separate 'readwind' process.
356! Optionally measure cpu/wall time.
357!
358!***********************************************************************
359    implicit none
360
361    if (mp_time_barrier) then
362      call cpu_time(mp_barrier_time_beg)
363      mp_barrier_wtime_beg = mpi_wtime()
364    endif
365
366    call MPI_BARRIER(mp_comm_used, mp_ierr)
367    if (mp_ierr /= 0) goto 100
368
369    if (mp_time_barrier) then
370      call cpu_time(mp_barrier_time_end)
371      mp_barrier_wtime_end = mpi_wtime()
372
373      mp_barrier_time_total=mp_barrier_time_total+&
374           &(mp_barrier_time_end-mp_barrier_time_beg)
375      mp_barrier_wtime_total=mp_barrier_wtime_total+&
376           &(mp_barrier_wtime_end-mp_barrier_wtime_beg)
377    end if
378
379    goto 101
380
381100 write(*,*) '#### mpi_mod::mpif_mpi_barrier> ERROR ####', mp_ierr
382    stop
383
384101 end subroutine mpif_mpi_barrier
385
386
387  subroutine mpif_com_mod_allocate
388!*******************************************************************************   
389! Dynamic allocation of arrays (in serial code these have size maxpart)
390!
391! Not in use anymore, moved to com_mod for interoperability with serial version
392!
393! TODO: remove
394!*******************************************************************************
395    use com_mod
396
397    implicit none
398
399    integer :: array_size
400
401    array_size = maxpart_mpi
402
403    allocate(itra1(array_size),npoint(array_size),&
404         & nclass(array_size),idt(array_size),&
405         & itramem(array_size),itrasplit(array_size),&
406         & xtra1(array_size),ytra1(array_size),&
407         & ztra1(array_size),xmass1(array_size, maxspec))
408
409    allocate(uap(array_size),ucp(array_size),&
410         & uzp(array_size),us(array_size),&
411         & vs(array_size),&
412         & ws(array_size),cbt(array_size))
413
414
415  end subroutine mpif_com_mod_allocate
416
417
418  subroutine mpif_tm_send_dims
419!***********************************************************************
420! Distribute array dimensions from pid0 to all processes.
421! numpart_mpi/numpart is sent to allow dynamic allocation
422!
423! Currently not in use.
424!
425! Variables of similar type (integer, real) are placed in an array
426! and sent collectively, to avoid the overhead associated with individual
427! MPI send/receives
428!
429!
430!***********************************************************************
431    use com_mod
432
433    implicit none
434
435    integer :: add_part=0
436
437    call MPI_Bcast(numpart, 1, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
438
439! MPI subgroup does the particle-loop
440    if (lmp_use_reader) then
441      if (mod(numpart,mp_partgroup_np).ne.0) add_part=1
442      numpart_mpi=int(numpart/mp_partgroup_np)+add_part
443    else
444      if (mod(numpart,mp_np).ne.0) add_part=1
445      numpart_mpi=int(numpart/mp_np)+add_part
446    end if
447
448
449! redefine numpart as 'numpart per process' throughout the code
450!**************************************************************
451    numpart = numpart_mpi
452
453  end subroutine mpif_tm_send_dims
454
455
456  subroutine mpif_tm_send_vars
457!***********************************************************************
458! Distribute particle variables from pid0 to all processes.
459! Called from timemanager
460! *NOT IN USE* at the moment, but can be useful for debugging
461!
462!***********************************************************************
463    use com_mod
464
465    implicit none
466
467    integer :: i
468
469! Time for MPI communications
470!****************************
471    if (mp_measure_time) call mpif_mtime('commtime',0)
472
473! Distribute variables, send from pid 0 to other processes
474!**********************************************************************
475
476! integers
477    if (lroot) then
478      call MPI_SCATTER(npoint,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
479           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
480      if (mp_ierr /= 0) goto 600
481      call MPI_SCATTER(idt,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
482           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
483      if (mp_ierr /= 0) goto 600
484      call MPI_SCATTER(itra1,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
485           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
486      if (mp_ierr /= 0) goto 600
487      call MPI_SCATTER(nclass,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
488           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
489      if (mp_ierr /= 0) goto 600
490      call MPI_SCATTER(itramem,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
491           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
492      if (mp_ierr /= 0) goto 600
493
494! int2
495      call MPI_SCATTER(cbt,numpart_mpi,MPI_INTEGER2,MPI_IN_PLACE,&
496           &numpart_mpi,MPI_INTEGER2,id_root,mp_comm_used,mp_ierr)
497      if (mp_ierr /= 0) goto 600
498
499! real
500      call MPI_SCATTER(uap,numpart_mpi,mp_sp,MPI_IN_PLACE,&
501           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
502      if (mp_ierr /= 0) goto 600
503      call MPI_SCATTER(ucp,numpart_mpi,mp_sp,MPI_IN_PLACE,&
504           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
505      if (mp_ierr /= 0) goto 600
506      call MPI_SCATTER(uzp,numpart_mpi,mp_sp,MPI_IN_PLACE,&
507           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
508      if (mp_ierr /= 0) goto 600
509
510      call MPI_SCATTER(us,numpart_mpi,mp_sp,MPI_IN_PLACE,&
511           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
512      if (mp_ierr /= 0) goto 600
513      call MPI_SCATTER(vs,numpart_mpi,mp_sp,MPI_IN_PLACE,&
514           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
515      if (mp_ierr /= 0) goto 600
516      call MPI_SCATTER(ws,numpart_mpi,mp_sp,MPI_IN_PLACE,&
517           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
518      if (mp_ierr /= 0) goto 600
519
520      call MPI_SCATTER(xtra1,numpart_mpi,mp_dp,MPI_IN_PLACE,&
521           &numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
522      if (mp_ierr /= 0) goto 600
523      call MPI_SCATTER(ytra1,numpart_mpi,mp_dp,MPI_IN_PLACE,&
524           &numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
525      if (mp_ierr /= 0) goto 600
526      call MPI_SCATTER(ztra1,numpart_mpi,mp_sp,MPI_IN_PLACE,&
527           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
528      if (mp_ierr /= 0) goto 600
529
530      do i=1,nspec
531        call MPI_SCATTER(xmass1(:,i),numpart_mpi,mp_sp,MPI_IN_PLACE,&
532             &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
533        if (mp_ierr /= 0) goto 600
534      end do
535
536    else ! (mp_pid >= 1)
537! integers
538      call MPI_SCATTER(npoint,numpart_mpi,MPI_INTEGER,npoint,&
539           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
540      if (mp_ierr /= 0) goto 600
541      call MPI_SCATTER(idt,numpart_mpi,MPI_INTEGER,idt,&
542           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
543      if (mp_ierr /= 0) goto 600
544      call MPI_SCATTER(itra1,numpart_mpi,MPI_INTEGER,itra1,&
545           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
546      if (mp_ierr /= 0) goto 600
547      call MPI_SCATTER(nclass,numpart_mpi,MPI_INTEGER,nclass,&
548           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
549      if (mp_ierr /= 0) goto 600
550      call MPI_SCATTER(itramem,numpart_mpi,MPI_INTEGER,itramem,&
551           &numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
552      if (mp_ierr /= 0) goto 600
553
554! int2
555      call MPI_SCATTER(cbt,numpart_mpi,MPI_INTEGER2,cbt,&
556           &numpart_mpi,MPI_INTEGER2,id_root,mp_comm_used,mp_ierr)
557      if (mp_ierr /= 0) goto 600
558
559! reals
560      call MPI_SCATTER(uap,numpart_mpi,mp_sp,uap,&
561           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
562      if (mp_ierr /= 0) goto 600
563      call MPI_SCATTER(ucp,numpart_mpi,mp_sp,ucp,&
564           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
565      if (mp_ierr /= 0) goto 600
566      call MPI_SCATTER(uzp,numpart_mpi,mp_sp,uzp,&
567           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
568      if (mp_ierr /= 0) goto 600
569
570      call MPI_SCATTER(us,numpart_mpi,mp_sp,us,&
571           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
572      if (mp_ierr /= 0) goto 600
573      call MPI_SCATTER(vs,numpart_mpi,mp_sp,vs,&
574           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
575      if (mp_ierr /= 0) goto 600
576      call MPI_SCATTER(ws,numpart_mpi,mp_sp,ws,&
577           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
578      if (mp_ierr /= 0) goto 600
579
580      call MPI_SCATTER(xtra1,numpart_mpi,mp_dp,xtra1,&
581           &numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
582      if (mp_ierr /= 0) goto 600
583      call MPI_SCATTER(ytra1,numpart_mpi,mp_dp,ytra1,&
584           &numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
585      if (mp_ierr /= 0) goto 600
586      call MPI_SCATTER(ztra1,numpart_mpi,mp_sp,ztra1,&
587           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
588      if (mp_ierr /= 0) goto 600
589
590      do i=1,nspec
591        call MPI_SCATTER(xmass1(:,i),numpart_mpi,mp_sp,xmass1(:,i),&
592             &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
593        if (mp_ierr /= 0) goto 600
594      end do
595
596    end if !(lroot)
597
598    if (mp_measure_time) call mpif_mtime('commtime',1)
599
600    goto 601
601
602600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
603    stop
604601 end subroutine mpif_tm_send_vars
605
606
607  subroutine mpif_tm_collect_vars
608!*******************************************************************************
609! Collect particle data to PID 0 from all processes.                           *
610!                                                                              *
611! This can be called from timemanager_mpi, before caclulating                  *
612! concentrations/writing output grids.                                         *
613!                                                                              *
614! Currently *not in use*, as each process calculates concentrations            *
615! separately, but can be useful for debugging                                  *
616!                                                                              *
617! To use this routine (together with mpif_tm_send_vars) to transfer data       *
618! to the MPI root process (useful for debugging), insert code like this        *
619! (in timemanager_mpi.f90)                                                     *
620!                                                                              *
621!      if (lroot) tot_numpart=numpart ! root stores total numpart              *
622!      call mpif_tm_send_dims                                                  *
623!      if (numpart>1) then                                                     *
624!        call mpif_tm_send_vars                                                *
625!      end if                                                                  *
626!                                                                              *
627!    To collect data from all processes to root process after a parallel       *
628!    region:                                                                   *
629!                                                                              *
630!      if (numpart>0) then                                                     *
631!          if (lroot) numpart = tot_numpart                                    *
632!       call mpif_tm_collect_vars                                              *
633!      end if                                                                  *
634!                                                                              *
635! Then a section that begins with "if (lroot) ..." will behave like            *
636! serial flexpart                                                              *
637!                                                                              *
638!*******************************************************************************
639    use com_mod !, only: numpart, nspec, itra1, npoint, nclass
640
641    implicit none
642
643    integer :: i
644
645    logical :: use_mp_vars = .false.! use mp_... type allocated vars
646    logical :: use_in_place = .true.! use MPI_IN_PLACE to collect vars.
647
648
649! Time for MPI communications
650    if (mp_measure_time) call mpif_mtime('commtime',0)
651
652! Distribute variables, send from pid 0 to other processes
653!**********************************************************************
654
655    if (.not. use_mp_vars.and..not.use_in_place) then
656
657! Integers:
658      call MPI_GATHER(npoint, numpart_mpi, MPI_INTEGER, npoint, &
659           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
660      if (mp_ierr /= 0) goto 600
661      call MPI_GATHER(idt, numpart_mpi, MPI_INTEGER, idt, &
662           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
663      if (mp_ierr /= 0) goto 600
664      call MPI_GATHER(itra1, numpart_mpi, MPI_INTEGER, itra1, &
665           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
666      if (mp_ierr /= 0) goto 600
667      call MPI_GATHER(nclass, numpart_mpi, MPI_INTEGER, nclass, &
668           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
669      if (mp_ierr /= 0) goto 600
670      call MPI_GATHER(itramem, numpart_mpi, MPI_INTEGER, itramem, &
671           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
672      if (mp_ierr /= 0) goto 600
673
674! int2
675      call MPI_GATHER(cbt, numpart_mpi, MPI_INTEGER2, cbt, &
676           &numpart_mpi, MPI_INTEGER2, id_root, mp_comm_used, mp_ierr)
677      if (mp_ierr /= 0) goto 600
678
679! Reals:
680      call MPI_GATHER(uap, numpart_mpi, mp_sp, uap, &
681           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
682      if (mp_ierr /= 0) goto 600
683      call MPI_GATHER(ucp, numpart_mpi, mp_sp, ucp, &
684           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
685      if (mp_ierr /= 0) goto 600
686      call MPI_GATHER(uzp, numpart_mpi, mp_sp, uzp, &
687           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
688      if (mp_ierr /= 0) goto 600
689
690      call MPI_GATHER(us, numpart_mpi, mp_sp, us, &
691           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
692      if (mp_ierr /= 0) goto 600
693      call MPI_GATHER(vs, numpart_mpi, mp_sp, vs, &
694           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
695      if (mp_ierr /= 0) goto 600
696      call MPI_GATHER(ws, numpart_mpi, mp_sp, ws, &
697           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
698      if (mp_ierr /= 0) goto 600
699
700
701      call MPI_GATHER(xtra1, numpart_mpi, mp_dp, xtra1, &
702           &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
703      if (mp_ierr /= 0) goto 600
704      call MPI_GATHER(ytra1, numpart_mpi, mp_dp, ytra1, &
705           &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
706      if (mp_ierr /= 0) goto 600
707      call MPI_GATHER(ztra1, numpart_mpi, mp_sp, ztra1, &
708           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
709      if (mp_ierr /= 0) goto 600
710
711      do i=1, nspec
712        call MPI_GATHER(xmass1(:,i), numpart_mpi, mp_sp, xmass1(:,i), &
713             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
714        if (mp_ierr /= 0) goto 600
715
716      end do
717
718! Use variables npoint etc directly for communications
719!***********************************************************************
720    else if (use_in_place.and..not.use_mp_vars) then
721      if (lroot) then
722! Integers:
723        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, npoint, &
724             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
725        if (mp_ierr /= 0) goto 600
726        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, idt, &
727             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
728        if (mp_ierr /= 0) goto 600
729        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, itra1, &
730             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
731        if (mp_ierr /= 0) goto 600
732        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, nclass, &
733             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
734        if (mp_ierr /= 0) goto 600
735        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, itramem, &
736             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
737        if (mp_ierr /= 0) goto 600
738
739! int2
740        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER2, cbt, &
741             &numpart_mpi, MPI_INTEGER2, id_root, mp_comm_used, mp_ierr)
742        if (mp_ierr /= 0) goto 600
743
744! Reals:
745        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, uap, &
746             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
747        if (mp_ierr /= 0) goto 600
748        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ucp, &
749             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
750        if (mp_ierr /= 0) goto 600
751        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, uzp, &
752             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
753        if (mp_ierr /= 0) goto 600
754
755        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, us, &
756             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
757        if (mp_ierr /= 0) goto 600
758        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, vs, &
759             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
760        if (mp_ierr /= 0) goto 600
761        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ws, &
762             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
763        if (mp_ierr /= 0) goto 600
764
765
766        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_dp, xtra1, &
767             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
768        if (mp_ierr /= 0) goto 600
769        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_dp, ytra1, &
770             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
771        if (mp_ierr /= 0) goto 600
772        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ztra1, &
773             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
774        if (mp_ierr /= 0) goto 600
775
776        do i=1, nspec
777          call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, xmass1(:,i), &
778               &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
779          if (mp_ierr /= 0) goto 600
780        end do
781
782      else ! (for mp_pid >= 1)
783
784! Integers:
785        call MPI_GATHER(npoint, numpart_mpi, MPI_INTEGER, npoint, &
786             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
787        if (mp_ierr /= 0) goto 600
788        call MPI_GATHER(idt, numpart_mpi, MPI_INTEGER, idt, &
789             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
790        if (mp_ierr /= 0) goto 600
791        call MPI_GATHER(itra1, numpart_mpi, MPI_INTEGER, itra1, &
792             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
793        if (mp_ierr /= 0) goto 600
794        call MPI_GATHER(nclass, numpart_mpi, MPI_INTEGER, nclass, &
795             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
796        if (mp_ierr /= 0) goto 600
797        call MPI_GATHER(itramem, numpart_mpi, MPI_INTEGER, itramem, &
798             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
799        if (mp_ierr /= 0) goto 600
800
801! int2
802        call MPI_GATHER(cbt, numpart_mpi, MPI_INTEGER2, cbt, &
803             &numpart_mpi, MPI_INTEGER2, id_root, mp_comm_used, mp_ierr)
804        if (mp_ierr /= 0) goto 600
805
806! Reals:
807        call MPI_GATHER(uap, numpart_mpi, mp_sp, uap, &
808             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
809        if (mp_ierr /= 0) goto 600
810        call MPI_GATHER(ucp, numpart_mpi, mp_sp, ucp, &
811             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
812        if (mp_ierr /= 0) goto 600
813        call MPI_GATHER(uzp, numpart_mpi, mp_sp, uzp, &
814             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
815        if (mp_ierr /= 0) goto 600
816
817        call MPI_GATHER(us, numpart_mpi, mp_sp, us, &
818             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
819        if (mp_ierr /= 0) goto 600
820        call MPI_GATHER(vs, numpart_mpi, mp_sp, vs, &
821             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
822        if (mp_ierr /= 0) goto 600
823        call MPI_GATHER(ws, numpart_mpi, mp_sp, ws, &
824             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
825        if (mp_ierr /= 0) goto 600
826
827
828        call MPI_GATHER(xtra1, numpart_mpi, mp_dp, xtra1, &
829             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
830        if (mp_ierr /= 0) goto 600
831        call MPI_GATHER(ytra1, numpart_mpi, mp_dp, ytra1, &
832             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
833        if (mp_ierr /= 0) goto 600
834        call MPI_GATHER(ztra1, numpart_mpi, mp_sp, ztra1, &
835             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
836        if (mp_ierr /= 0) goto 600
837
838        do i=1, nspec
839          call MPI_GATHER(xmass1(:,i), numpart_mpi, mp_sp, xmass1(:,i), &
840               &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
841          if (mp_ierr /= 0) goto 600
842        end do
843      end if ! (lroot)
844    end if !  (use_in_place)
845
846    if (mp_measure_time) call mpif_mtime('commtime',1)
847
848    goto 601
849
850600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
851    stop
852601 end subroutine mpif_tm_collect_vars
853
854
855  subroutine mpif_gf_send_vars(memstat)
856!*******************************************************************************
857! DESCRIPTION
858!   Distribute 'getfield' variables from reader process
859!
860!   Called from timemanager
861!
862! NOTE
863!   This subroutine distributes windfields read from the reader process to
864!   all other processes. Usually one set of fields is transfered, but at
865!   first timestep when there are no fields in memory, both are transfered.
866!   MPI_Bcast is used, so implicitly all processes are synchronized at this
867!   step
868!
869!
870!*******************************************************************************
871    use com_mod
872    use par_mod,only: numwfmem
873
874    implicit none
875
876    integer, intent(in) :: memstat
877
878! Common array sizes used for communications
879    integer, parameter :: d3_size1 = nxmax*nymax*nzmax
880    integer, parameter :: d3_size2 = nxmax*nymax*nuvzmax
881    integer, parameter :: d2_size1 = nxmax*nymax
882    integer, parameter :: d2_size2 = nxmax*nymax*maxspec
883    integer, parameter :: d2_size3 = nxmax*nymax
884    integer, parameter :: d1_size1 = maxwf
885
886    integer :: d3s1,d3s2,d2s1,d2s2
887
888! li,ui    lower,upper indices of fields to be transfered (1-3, ui>=li)
889    integer :: li,ui
890
891! First time routine is called the unchangeable fields will be transfered   
892    logical :: first_call=.true.
893
894!**********************************************************************
895
896! Sizes of arrays transfered
897    d3s1=d3_size1
898    d3s2=d3_size2
899    d2s1=d2_size1
900    d2s2=d2_size2
901
902! Decide which field(s) need to be transfered
903    if (memstat.ge.32) then ! distribute 2 fields, to 1st/2nd indices
904      li=1; ui=2
905      d3s1=2*d3_size1
906      d3s2=2*d3_size2
907      d2s1=2*d2_size1
908      d2s2=2*d2_size2
909    else if (memstat.eq.17) then ! distribute 1 field, place on 1st index
910      li=1; ui=1
911    else if (memstat.eq.18) then ! distribute 1 field, place on 2nd index
912      li=2; ui=2
913    else if (memstat.eq.19) then ! distribute 1 field, place on 3nd index
914      li=3; ui=3
915    else
916      write(*,*) '#### mpi_mod::mpif_gf_send_vars> ERROR: ', &
917           & 'wrong value of memstat, exiting ####', memstat
918      stop
919    end if
920
921
922! Time for MPI communications
923    if (mp_measure_time) call mpif_mtime('commtime',0)
924
925! Send variables from getfield process (id_read) to other processes
926!**********************************************************************
927
928! The non-reader processes need to know if cloud water were read.
929    call MPI_Bcast(readclouds,1,MPI_LOGICAL,id_read,MPI_COMM_WORLD,mp_ierr)
930    if (mp_ierr /= 0) goto 600
931
932! Static fields/variables sent only at startup
933    if (first_call) then
934      call MPI_Bcast(oro(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
935      if (mp_ierr /= 0) goto 600
936      call MPI_Bcast(excessoro(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
937      if (mp_ierr /= 0) goto 600
938      call MPI_Bcast(lsm(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
939      if (mp_ierr /= 0) goto 600
940      call MPI_Bcast(xlanduse(:,:,:),d2_size3*numclass,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
941      if (mp_ierr /= 0) goto 600
942      call MPI_Bcast(wftime,d1_size1,MPI_INTEGER,id_read,MPI_COMM_WORLD,mp_ierr)
943      if (mp_ierr /= 0) goto 600
944      call MPI_Bcast(numbwf,1,MPI_INTEGER,id_read,MPI_COMM_WORLD,mp_ierr)
945      if (mp_ierr /= 0) goto 600
946      call MPI_Bcast(nmixz,1,MPI_INTEGER,id_read,MPI_COMM_WORLD,mp_ierr)
947      if (mp_ierr /= 0) goto 600
948      call MPI_Bcast(height,nzmax,MPI_INTEGER,id_read,MPI_COMM_WORLD,mp_ierr)
949      if (mp_ierr /= 0) goto 600
950
951      first_call=.false.
952    endif
953
954    call MPI_Bcast(uu(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
955    if (mp_ierr /= 0) goto 600
956    call MPI_Bcast(vv(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
957    if (mp_ierr /= 0) goto 600
958    call MPI_Bcast(uupol(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
959    if (mp_ierr /= 0) goto 600
960    call MPI_Bcast(vvpol(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
961    if (mp_ierr /= 0) goto 600
962    call MPI_Bcast(ww(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
963    if (mp_ierr /= 0) goto 600
964    call MPI_Bcast(tt(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
965    if (mp_ierr /= 0) goto 600
966    call MPI_Bcast(rho(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
967    if (mp_ierr /= 0) goto 600
968    call MPI_Bcast(drhodz(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
969    if (mp_ierr /= 0) goto 600
970    call MPI_Bcast(tth(:,:,:,li:ui),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
971    if (mp_ierr /= 0) goto 600
972    call MPI_Bcast(qvh(:,:,:,li:ui),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
973    if (mp_ierr /= 0) goto 600
974    call MPI_Bcast(qv(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
975    if (mp_ierr /= 0) goto 600
976    call MPI_Bcast(pv(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
977    if (mp_ierr /= 0) goto 600
978    call MPI_Bcast(clouds(:,:,:,li:ui),d3s1,MPI_INTEGER1,id_read,MPI_COMM_WORLD,mp_ierr)
979    if (mp_ierr /= 0) goto 600
980
981! cloud water/ice:
982    if (readclouds) then
983      ! call MPI_Bcast(icloud_stats(:,:,:,li:ui),d2s1*5,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
984      ! if (mp_ierr /= 0) goto 600
985      call MPI_Bcast(ctwc(:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
986      if (mp_ierr /= 0) goto 600
987      ! call MPI_Bcast(clwc(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
988      ! if (mp_ierr /= 0) goto 600
989      ! call MPI_Bcast(ciwc(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
990      ! if (mp_ierr /= 0) goto 600
991    end if
992
993! 2D fields
994    call MPI_Bcast(cloudsh(:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
995    if (mp_ierr /= 0) goto 600
996    call MPI_Bcast(vdep(:,:,:,li:ui),d2s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
997    if (mp_ierr /= 0) goto 600
998    call MPI_Bcast(ps(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
999    if (mp_ierr /= 0) goto 600
1000    call MPI_Bcast(sd(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1001    if (mp_ierr /= 0) goto 600
1002    call MPI_Bcast(tcc(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1003    if (mp_ierr /= 0) goto 600
1004    call MPI_Bcast(tt2(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1005    if (mp_ierr /= 0) goto 600
1006    call MPI_Bcast(td2(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1007    if (mp_ierr /= 0) goto 600
1008    call MPI_Bcast(lsprec(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1009    if (mp_ierr /= 0) goto 600
1010    call MPI_Bcast(convprec(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1011    if (mp_ierr /= 0) goto 600
1012    call MPI_Bcast(ustar(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1013    if (mp_ierr /= 0) goto 600
1014    call MPI_Bcast(wstar(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1015    if (mp_ierr /= 0) goto 600
1016    call MPI_Bcast(hmix(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1017    if (mp_ierr /= 0) goto 600
1018    call MPI_Bcast(tropopause(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1019    if (mp_ierr /= 0) goto 600
1020    call MPI_Bcast(oli(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1021    if (mp_ierr /= 0) goto 600
1022
1023    if (mp_measure_time) call mpif_mtime('commtime',1)
1024
1025    goto 601
1026
1027600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
1028    stop
1029
1030601 end subroutine mpif_gf_send_vars
1031
1032
1033  subroutine mpif_gf_send_vars_nest(memstat)
1034!***********************************************************************
1035! DESCRIPTION
1036!   Distribute 'getfield' variables from reader process to all processes.
1037!   For nested fields
1038!
1039!   Called from timemanager
1040!
1041! NOTE
1042!   This subroutine distributes nested windfields from the reader process to
1043!   all other processes. Usually one set of fields is transfered, but at
1044!   first timestep when there are no fields in memory, both are transfered.
1045!   MPI_Bcast is used, so implicitly all processes are synchronized at this
1046!   step
1047!
1048!
1049!***********************************************************************
1050    use com_mod
1051    use par_mod,only: numwfmem
1052
1053    implicit none
1054
1055    integer, intent(in) :: memstat
1056    integer :: i
1057! li,ui    lower,upper indices of fields to be transfered (1-3, ui>=li)
1058    integer :: li,ui
1059
1060! Common array sizes used for communications
1061    integer :: d3_size1 = nxmaxn*nymaxn*nzmax
1062    integer :: d3_size2 = nxmaxn*nymaxn*nuvzmax
1063    integer :: d2_size1 = nxmaxn*nymaxn
1064    integer :: d2_size2 = nxmaxn*nymaxn*maxspec
1065    integer :: d2_size3 = nxmaxn*nymaxn
1066
1067    integer :: d3s1,d3s2,d2s1,d2s2
1068
1069! First time routine is called the unchangeable fields will be transfered   
1070    logical :: first_call=.true.
1071
1072!**********************************************************************
1073
1074! Sizes of arrays transfered
1075    d3s1=d3_size1
1076    d3s2=d3_size2
1077    d2s1=d2_size1
1078    d2s2=d2_size2
1079
1080! Decide which field(s) need to be transfered
1081    if (memstat.ge.32) then ! distribute 2 fields
1082      li=1; ui=2
1083      d3s1=2*d3_size1
1084      d3s2=2*d3_size2
1085      d2s1=2*d2_size1
1086      d2s2=2*d2_size2
1087    else if (memstat.eq.17) then ! distribute 1 field, on 1st index
1088      li=1; ui=1
1089    else if (memstat.eq.18) then ! distribute 1 field, on 2nd index
1090      li=2; ui=2
1091    else if (memstat.eq.19) then ! distribute 1 field, on 3nd index
1092      li=3; ui=3
1093    else
1094      write(*,*) '#### mpi_mod::mpif_gf_send_vars_nest> ERROR: ', &
1095           & 'wrong value of memstat, exiting ####', memstat
1096      stop
1097    end if
1098
1099
1100! Time for MPI communications
1101    if (mp_measure_time) call mpif_mtime('commtime',0)
1102
1103! Distribute variables, send from getfield process (id_read) to other
1104! processes
1105!**********************************************************************
1106
1107! The non-reader processes need to know if cloud water were read.
1108    call MPI_Bcast(readclouds_nest,maxnests,MPI_LOGICAL,id_read,MPI_COMM_WORLD,mp_ierr)
1109    if (mp_ierr /= 0) goto 600
1110
1111! Static fields/variables sent only at startup
1112    if (first_call) then
1113      call MPI_Bcast(oron(:,:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1114      if (mp_ierr /= 0) goto 600
1115      call MPI_Bcast(excessoron(:,:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1116      if (mp_ierr /= 0) goto 600
1117      call MPI_Bcast(lsmn(:,:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1118      if (mp_ierr /= 0) goto 600
1119      call MPI_Bcast(xlandusen(:,:,:,:),d2_size3*numclass,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1120      if (mp_ierr /= 0) goto 600
1121      first_call=.false.
1122    end if
1123
1124! MPI prefers contiguous arrays for sending (else a buffer is created),
1125! hence the loop over nests
1126!**********************************************************************
1127    do i=1, numbnests
1128! 3D fields
1129      call MPI_Bcast(uun(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1130      if (mp_ierr /= 0) goto 600
1131      call MPI_Bcast(vvn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1132      if (mp_ierr /= 0) goto 600
1133      call MPI_Bcast(wwn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1134      if (mp_ierr /= 0) goto 600
1135      call MPI_Bcast(ttn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1136      if (mp_ierr /= 0) goto 600
1137      call MPI_Bcast(rhon(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1138      if (mp_ierr /= 0) goto 600
1139      call MPI_Bcast(drhodzn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1140      if (mp_ierr /= 0) goto 600
1141      call MPI_Bcast(tthn(:,:,:,li:ui,i),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1142      if (mp_ierr /= 0) goto 600
1143      call MPI_Bcast(qvhn(:,:,:,li:ui,i),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1144      if (mp_ierr /= 0) goto 600
1145      call MPI_Bcast(qvn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1146      if (mp_ierr /= 0) goto 600
1147      call MPI_Bcast(pvn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1148      if (mp_ierr /= 0) goto 600
1149      call MPI_Bcast(cloudsn(:,:,:,li:ui,i),d3s1,MPI_INTEGER1,id_read,MPI_COMM_WORLD,mp_ierr)
1150      if (mp_ierr /= 0) goto 600
1151
1152! cloud water/ice:
1153    if (readclouds_nest(i)) then
1154      ! call MPI_Bcast(icloud_stats(:,:,:,li:ui),d2s1*5,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1155      ! if (mp_ierr /= 0) goto 600
1156      call MPI_Bcast(ctwcn(:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1157      if (mp_ierr /= 0) goto 600
1158    end if
1159
1160! 2D fields
1161      call MPI_Bcast(cloudshn(:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1162      if (mp_ierr /= 0) goto 600
1163      call MPI_Bcast(vdepn(:,:,:,li:ui,i),d2s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1164      if (mp_ierr /= 0) goto 600
1165      call MPI_Bcast(psn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1166      if (mp_ierr /= 0) goto 600
1167      call MPI_Bcast(sdn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1168      if (mp_ierr /= 0) goto 600
1169      call MPI_Bcast(tccn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1170      if (mp_ierr /= 0) goto 600
1171      call MPI_Bcast(tt2n(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1172      if (mp_ierr /= 0) goto 600
1173      call MPI_Bcast(td2n(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1174      if (mp_ierr /= 0) goto 600
1175      call MPI_Bcast(lsprecn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1176      if (mp_ierr /= 0) goto 600
1177      call MPI_Bcast(convprecn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1178      if (mp_ierr /= 0) goto 600
1179      call MPI_Bcast(ustarn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1180      if (mp_ierr /= 0) goto 600
1181      call MPI_Bcast(wstarn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1182      if (mp_ierr /= 0) goto 600
1183      call MPI_Bcast(olin(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1184      if (mp_ierr /= 0) goto 600
1185      call MPI_Bcast(hmixn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1186      if (mp_ierr /= 0) goto 600
1187      call MPI_Bcast(tropopausen(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1188      if (mp_ierr /= 0) goto 600
1189    end do
1190
1191    if (mp_measure_time) call mpif_mtime('commtime',1)
1192
1193    goto 601
1194
1195600 write(*,*) "mpi_mod> mp_ierr \= 0",mp_ierr
1196    stop
1197
1198601 end subroutine mpif_gf_send_vars_nest
1199
1200
1201  subroutine mpif_gf_send_vars_async(memstat)
1202!*******************************************************************************
1203! DESCRIPTION
1204!   Distribute 'getfield' variables from reader process to all processes.
1205!   Called from timemanager by reader process only
1206!
1207! NOTE
1208!   This version uses asynchronious sends. The newest fields are sent in the
1209!   background, so calculations can continue while
1210!   MPI communications are performed.
1211!
1212!   The way the MPI tags/requests are sequenced means that this routine must
1213!   carefully match routine 'mpif_gf_recv_vars_async'
1214!
1215! VARIABLES
1216!   memstat -- input, for resolving pointer to windfield index being read
1217!   mind    -- index where to place new fields
1218!
1219!
1220!*******************************************************************************
1221    use com_mod
1222
1223    implicit none
1224
1225    integer, intent(in) :: memstat
1226    integer :: mind
1227    integer :: dest,i
1228
1229! Common array sizes used for communications
1230    integer :: d3s1 = nxmax*nymax*nzmax
1231    integer :: d3s2 = nxmax*nymax*nuvzmax
1232    integer :: d2s1 = nxmax*nymax
1233    integer :: d2s2 = nxmax*nymax*maxspec
1234!    integer :: d1s1 = maxwf
1235
1236!*******************************************************************************
1237
1238! At the time the send is posted, the reader process is one step further
1239! in the permutation of memstat compared with the receiving processes
1240
1241    if (memstat.ge.32) then
1242! last read was synchronous, to indices 1 and 2, 3 is free
1243      write(*,*) "#### mpi_mod::mpif_gf_send_vars_async> ERROR: &
1244           & memstat>=32 should never happen here."
1245      stop
1246    else if (memstat.eq.17) then
1247! old fields on 1,2, send 3
1248      mind=3
1249    else if (memstat.eq.18) then
1250! old fields on 2,3, send 1
1251      mind=1
1252    else if (memstat.eq.19) then
1253! old fields on 3,1, send 2
1254      mind=2
1255    else
1256      write(*,*) "#### mpi_mod::mpif_gf_send_vars_async> ERROR: &
1257           & invalid memstat"
1258    end if
1259
1260    if (mp_dev_mode) then
1261      if (mp_pid.ne.id_read) then
1262        write(*,*) 'MPI_DEV_MODE: error calling mpif_gf_send_vars_async'
1263      end if
1264    end if
1265
1266    if (mp_dev_mode) write(*,*) '## in mpif_gf_send_vars_async, memstat:', memstat
1267
1268! Time for MPI communications
1269    if (mp_measure_time) call mpif_mtime('commtime',0)
1270
1271! Loop over receiving processes, initiate data sending
1272!*****************************************************
1273
1274    do dest=0,mp_np-1 ! mp_np-2 will also work if last proc reserved for reading
1275! TODO: use mp_partgroup_np here
1276      if (dest.eq.id_read) cycle
1277      i=dest*nvar_async
1278      call MPI_Isend(uu(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1279      if (mp_ierr /= 0) goto 600
1280      i=i+1
1281      call MPI_Isend(vv(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1282      if (mp_ierr /= 0) goto 600
1283      i=i+1
1284      call MPI_Isend(uupol(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1285      if (mp_ierr /= 0) goto 600
1286      i=i+1
1287      call MPI_Isend(vvpol(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1288      if (mp_ierr /= 0) goto 600
1289      i=i+1
1290      call MPI_Isend(ww(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1291      if (mp_ierr /= 0) goto 600
1292      i=i+1
1293      call MPI_Isend(tt(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1294      if (mp_ierr /= 0) goto 600
1295      i=i+1
1296      call MPI_Isend(rho(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1297      if (mp_ierr /= 0) goto 600
1298      i=i+1
1299      call MPI_Isend(drhodz(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1300      if (mp_ierr /= 0) goto 600
1301      i=i+1
1302      call MPI_Isend(tth(:,:,:,mind),d3s2,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1303      if (mp_ierr /= 0) goto 600
1304      i=i+1
1305      call MPI_Isend(qvh(:,:,:,mind),d3s2,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1306      if (mp_ierr /= 0) goto 600
1307      i=i+1
1308      call MPI_Isend(qv(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1309      if (mp_ierr /= 0) goto 600
1310      i=i+1
1311      call MPI_Isend(pv(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1312      if (mp_ierr /= 0) goto 600
1313      i=i+1
1314      call MPI_Isend(clouds(:,:,:,mind),d3s1,MPI_INTEGER1,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1315      i=i+1
1316      if (mp_ierr /= 0) goto 600
1317
1318      call MPI_Isend(cloudsh(:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1319      if (mp_ierr /= 0) goto 600
1320      i=i+1
1321      call MPI_Isend(vdep(:,:,:,mind),d2s2,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1322      if (mp_ierr /= 0) goto 600
1323      i=i+1
1324      call MPI_Isend(ps(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1325      if (mp_ierr /= 0) goto 600
1326      i=i+1
1327      call MPI_Isend(sd(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1328      if (mp_ierr /= 0) goto 600
1329      i=i+1
1330      call MPI_Isend(tcc(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1331      if (mp_ierr /= 0) goto 600
1332      i=i+1
1333      call MPI_Isend(tt2(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1334      if (mp_ierr /= 0) goto 600
1335      i=i+1
1336      call MPI_Isend(td2(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1337      if (mp_ierr /= 0) goto 600
1338      i=i+1
1339      call MPI_Isend(lsprec(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1340      if (mp_ierr /= 0) goto 600
1341      i=i+1
1342      call MPI_Isend(convprec(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1343      if (mp_ierr /= 0) goto 600
1344      i=i+1
1345      call MPI_Isend(ustar(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1346      if (mp_ierr /= 0) goto 600
1347      i=i+1
1348      call MPI_Isend(wstar(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1349      if (mp_ierr /= 0) goto 600
1350      i=i+1
1351      call MPI_Isend(hmix(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1352      if (mp_ierr /= 0) goto 600
1353      i=i+1
1354      call MPI_Isend(tropopause(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1355      if (mp_ierr /= 0) goto 600
1356      i=i+1
1357      call MPI_Isend(oli(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1358      if (mp_ierr /= 0) goto 600
1359
1360! Send cloud water if it exists. Increment counter always (as on receiving end)
1361      if (readclouds) then
1362        i=i+1
1363        ! call MPI_Isend(icloud_stats(:,:,:,mind),d2s1*5,mp_sp,dest,tm1,&
1364        !      &MPI_COMM_WORLD,reqs(i),mp_ierr)
1365        call MPI_Isend(ctwc(:,:,mind),d2s1,mp_sp,dest,tm1,&
1366             &MPI_COMM_WORLD,reqs(i),mp_ierr)
1367
1368        if (mp_ierr /= 0) goto 600
1369
1370        ! call MPI_Isend(clwc(:,:,:,mind),d3s1,mp_sp,dest,tm1,&
1371        !      &MPI_COMM_WORLD,reqs(i),mp_ierr)
1372        ! if (mp_ierr /= 0) goto 600
1373        ! i=i+1
1374
1375        ! call MPI_Isend(ciwc(:,:,:,mind),d3s1,mp_sp,dest,tm1,&
1376        !      &MPI_COMM_WORLD,reqs(i),mp_ierr)
1377        ! if (mp_ierr /= 0) goto 600
1378
1379      end if
1380    end do
1381
1382    if (mp_measure_time) call mpif_mtime('commtime',1)
1383
1384    goto 601
1385
1386600 write(*,*) "#### mpi_mod::mpif_gf_send_vars_async> mp_ierr \= 0", mp_ierr
1387    stop
1388
1389601 end subroutine mpif_gf_send_vars_async
1390
1391
1392  subroutine mpif_gf_recv_vars_async(memstat)
1393!*******************************************************************************
1394! DESCRIPTION
1395!   Receive 'getfield' variables from reader process.
1396!   Called from timemanager by all processes except reader
1397!
1398! NOTE
1399!   This version uses asynchronious communications.
1400!
1401! VARIABLES
1402!   memstat -- input, used to resolve windfield index being received
1403!
1404!
1405!*******************************************************************************
1406    use com_mod
1407
1408    implicit none
1409
1410    integer, intent(in) :: memstat
1411    integer :: mind,j
1412
1413! Common array sizes used for communications
1414    integer :: d3s1 = nxmax*nymax*nzmax
1415    integer :: d3s2 = nxmax*nymax*nuvzmax
1416    integer :: d2s1 = nxmax*nymax
1417    integer :: d2s2 = nxmax*nymax*maxspec
1418    !integer :: d1_size1 = maxwf
1419
1420!    integer :: d3s1,d3s2,d2s1,d2s2
1421!*******************************************************************************
1422
1423! At the time this immediate receive is posted, memstat is the state of
1424! windfield indices at the previous/current time. From this, the future
1425! state is deduced.
1426    if (memstat.eq.32) then
1427! last read was synchronous to indices 1 and 2, 3 is free
1428      mind=3
1429    else if (memstat.eq.17) then
1430! last read was asynchronous to index 3, 1 is free
1431      mind=1
1432    else if (memstat.eq.18) then
1433! last read was asynchronous to index 1, 2 is free
1434      mind=2
1435    else if (memstat.eq.19) then
1436! last read was asynchronous to index 2, 3 is free
1437      mind=3
1438    else
1439! illegal state
1440      write(*,*) 'mpi_mod> FLEXPART ERROR: Illegal memstat value. Exiting.'
1441      stop
1442    end if
1443
1444    if (mp_dev_mode) then
1445      if (mp_pid.eq.id_read) then
1446        write(*,*) 'MPI_DEV_MODE: error calling mpif_gf_recv_vars_async'
1447      end if
1448    end if
1449
1450! Time for MPI communications
1451    if (mp_measure_time) call mpif_mtime('commtime',0)
1452
1453    if (mp_dev_mode) write(*,*) '## in mpif_gf_send_vars_async, memstat:', memstat
1454
1455! Initiate receiving of data
1456!***************************
1457
1458! Get MPI tags/requests for communications
1459    j=mp_pid*nvar_async
1460    call MPI_Irecv(uu(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1461         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1462    if (mp_ierr /= 0) goto 600
1463    j=j+1
1464    call MPI_Irecv(vv(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1465         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1466    if (mp_ierr /= 0) goto 600
1467    j=j+1
1468    call MPI_Irecv(uupol(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1469         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1470    if (mp_ierr /= 0) goto 600
1471    j=j+1
1472    call MPI_Irecv(vvpol(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1473         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1474    if (mp_ierr /= 0) goto 600
1475    j=j+1
1476    call MPI_Irecv(ww(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1477         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1478    if (mp_ierr /= 0) goto 600
1479    j=j+1
1480    call MPI_Irecv(tt(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1481         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1482    if (mp_ierr /= 0) goto 600
1483    j=j+1
1484    call MPI_Irecv(rho(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1485         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1486    if (mp_ierr /= 0) goto 600
1487    j=j+1
1488    call MPI_Irecv(drhodz(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1489         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1490    if (mp_ierr /= 0) goto 600
1491    j=j+1
1492    call MPI_Irecv(tth(:,:,:,mind),d3s2,mp_sp,id_read,MPI_ANY_TAG,&
1493         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1494    if (mp_ierr /= 0) goto 600
1495    j=j+1
1496    call MPI_Irecv(qvh(:,:,:,mind),d3s2,mp_sp,id_read,MPI_ANY_TAG,&
1497         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1498    if (mp_ierr /= 0) goto 600
1499    j=j+1
1500
1501    call MPI_Irecv(qv(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1502         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1503    if (mp_ierr /= 0) goto 600
1504    j=j+1
1505    call MPI_Irecv(pv(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1506         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1507    if (mp_ierr /= 0) goto 600
1508    j=j+1
1509    call MPI_Irecv(clouds(:,:,:,mind),d3s1,MPI_INTEGER1,id_read,MPI_ANY_TAG,&
1510         &MPI_COMM_WORLD,reqs(j),mp_ierr)   
1511    if (mp_ierr /= 0) goto 600
1512    j=j+1
1513
1514    call MPI_Irecv(cloudsh(:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1515         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1516    if (mp_ierr /= 0) goto 600
1517    j=j+1
1518    call MPI_Irecv(vdep(:,:,:,mind),d2s2,mp_sp,id_read,MPI_ANY_TAG,&
1519         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1520    if (mp_ierr /= 0) goto 600
1521    j=j+1
1522    call MPI_Irecv(ps(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1523         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1524    if (mp_ierr /= 0) goto 600
1525    j=j+1
1526    call MPI_Irecv(sd(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1527         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1528    if (mp_ierr /= 0) goto 600
1529    j=j+1
1530    call MPI_Irecv(tcc(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1531         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1532    if (mp_ierr /= 0) goto 600
1533    j=j+1
1534    call MPI_Irecv(tt2(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1535         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1536    if (mp_ierr /= 0) goto 600
1537    j=j+1
1538    call MPI_Irecv(td2(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1539         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1540    if (mp_ierr /= 0) goto 600
1541    j=j+1
1542    call MPI_Irecv(lsprec(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1543         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1544    if (mp_ierr /= 0) goto 600
1545    j=j+1
1546    call MPI_Irecv(convprec(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1547         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1548    if (mp_ierr /= 0) goto 600
1549
1550    call MPI_Irecv(ustar(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1551         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1552    if (mp_ierr /= 0) goto 600
1553    j=j+1
1554    call MPI_Irecv(wstar(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1555         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1556    if (mp_ierr /= 0) goto 600
1557    j=j+1
1558    call MPI_Irecv(hmix(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1559         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1560    if (mp_ierr /= 0) goto 600
1561    j=j+1
1562    call MPI_Irecv(tropopause(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1563         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1564    if (mp_ierr /= 0) goto 600
1565    j=j+1
1566    call MPI_Irecv(oli(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1567         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1568    if (mp_ierr /= 0) goto 600
1569
1570
1571! Post request for clwc. These data are possibly not sent, request must then be cancelled
1572! For now assume that data at all steps either have or do not have water
1573    if (readclouds) then
1574      j=j+1
1575
1576      ! call MPI_Irecv(icloud_stats(:,:,:,mind),d2s1*5,mp_sp,id_read,MPI_ANY_TAG,&
1577      !      &MPI_COMM_WORLD,reqs(j),mp_ierr)
1578      call MPI_Irecv(ctwc(:,:,mind),d2s1*5,mp_sp,id_read,MPI_ANY_TAG,&
1579           &MPI_COMM_WORLD,reqs(j),mp_ierr)
1580      if (mp_ierr /= 0) goto 600
1581
1582      ! call MPI_Irecv(clwc(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1583      !      &MPI_COMM_WORLD,reqs(j),mp_ierr)   
1584      ! if (mp_ierr /= 0) goto 600
1585      ! j=j+1
1586      ! call MPI_Irecv(ciwc(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1587      !      &MPI_COMM_WORLD,reqs(j),mp_ierr)   
1588      ! if (mp_ierr /= 0) goto 600
1589
1590    end if
1591
1592
1593    if (mp_measure_time) call mpif_mtime('commtime',1)
1594
1595    goto 601
1596
1597600 write(*,*) "#### mpi_mod::mpif_gf_recv_vars_async> MPI ERROR ", mp_ierr
1598    stop
1599
1600601 end subroutine mpif_gf_recv_vars_async
1601
1602
1603  subroutine mpif_gf_request
1604!*******************************************************************************
1605! DESCRIPTION
1606!   Check for completion of MPI_Isend/MPI_Irecv data transfer.
1607!   
1608!   
1609! NOTE
1610!   implicit synchronisation between all processes takes place here
1611!
1612!
1613!*******************************************************************************
1614!    use com_mod,only: readclouds
1615
1616    implicit none
1617
1618
1619    integer :: n_req !,i
1620
1621!***********************************************************************
1622
1623    n_req=nvar_async*mp_np
1624
1625    if (mp_measure_time) call mpif_mtime('commtime',0)
1626
1627!    call MPI_Wait(rm1,MPI_STATUS_IGNORE,mp_ierr)
1628
1629! TODO: cancel recv request if readclouds=.false.
1630!    if (readclouds) then
1631    call MPI_Waitall(n_req,reqs,MPI_STATUSES_IGNORE,mp_ierr)
1632!    endif
1633! else
1634!   do i = 0, nvar_async*mp_np-1
1635!     if (mod(i,27).eq.0 .or. mod(i,28).eq.0) then
1636!       call MPI_Cancel(reqs(i),mp_ierr)
1637!       cycle
1638!     end if
1639!     call MPI_Wait(reqs(i),MPI_STATUS_IGNORE,mp_ierr)
1640!   end do
1641! end if
1642
1643    if (mp_ierr /= 0) goto 600
1644
1645    if (mp_measure_time) call mpif_mtime('commtime',1)
1646
1647    goto 601
1648
1649600 write(*,*) "#### mpi_mod::mpif_gf_request> MPI ERROR ", mp_ierr
1650    stop
1651
1652601 end subroutine mpif_gf_request
1653
1654
1655  subroutine mpif_tm_reduce_grid
1656!***********************************************************************
1657! Collect grid variable to PID 0, adding from all processes.
1658!
1659! NOTE
1660!   - Older MPI implementations may lack the MPI_IN_PLACE operation.
1661!     If so use 1) below and change unc_mod
1662!
1663!***********************************************************************
1664    use com_mod
1665    use unc_mod
1666    use par_mod
1667
1668    implicit none
1669
1670    integer :: grid_size2d,grid_size3d
1671    integer, parameter :: rcpt_size=maxreceptor*maxspec
1672
1673!**********************************************************************
1674    grid_size2d=numxgrid*numygrid*maxspec* &
1675         & maxpointspec_act*nclassunc*maxageclass
1676    grid_size3d=numxgrid*numygrid*numzgrid*maxspec* &
1677         & maxpointspec_act*nclassunc*maxageclass
1678
1679
1680! Time for MPI communications
1681    if (mp_measure_time) call mpif_mtime('commtime',0)
1682
1683! 1) Using a separate grid (gridunc0) for received values
1684! call MPI_Reduce(gridunc, gridunc0, grid_size3d, mp_sp, MPI_SUM, id_root, &
1685!      & mp_comm_used, mp_ierr)
1686! if (mp_ierr /= 0) goto 600
1687
1688! 2) Using in-place reduction
1689    if (lroot) then
1690      call MPI_Reduce(MPI_IN_PLACE, gridunc, grid_size3d, mp_sp, MPI_SUM, id_root, &
1691           & mp_comm_used, mp_ierr)
1692      if (mp_ierr /= 0) goto 600
1693    else
1694      call MPI_Reduce(gridunc, gridunc, grid_size3d, mp_sp, MPI_SUM, id_root, &
1695           & mp_comm_used, mp_ierr)
1696    end if
1697
1698    if ((WETDEP).and.(ldirect.gt.0)) then
1699      call MPI_Reduce(wetgridunc, wetgridunc0, grid_size2d, mp_cp, MPI_SUM, id_root, &
1700           & mp_comm_used, mp_ierr)
1701      if (mp_ierr /= 0) goto 600
1702    end if
1703
1704    if ((DRYDEP).and.(ldirect.gt.0)) then
1705      call MPI_Reduce(drygridunc, drygridunc0, grid_size2d, mp_cp, MPI_SUM, id_root, &
1706           & mp_comm_used, mp_ierr)
1707      if (mp_ierr /= 0) goto 600
1708    end if
1709
1710! Receptor concentrations   
1711    if (lroot) then
1712      call MPI_Reduce(MPI_IN_PLACE,creceptor,rcpt_size,mp_sp,MPI_SUM,id_root, &
1713           & mp_comm_used,mp_ierr)
1714      if (mp_ierr /= 0) goto 600
1715    else
1716      call MPI_Reduce(creceptor,creceptor,rcpt_size,mp_sp,MPI_SUM,id_root, &
1717           & mp_comm_used,mp_ierr)
1718    end if
1719
1720    if (mp_measure_time) call mpif_mtime('commtime',1)
1721
1722    goto 601
1723
1724600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
1725    stop
1726
1727601 end subroutine mpif_tm_reduce_grid
1728
1729
1730  subroutine mpif_tm_reduce_grid_nest
1731!***********************************************************************
1732! Collect nested grids to PID 0, adding from all processes.
1733!
1734! NOTE
1735!   - Compiling with 'fcheck=all' (gfortran) will cause a run-time error
1736!     as wetgriduncn0 is only allocated for root process. Possibly valid
1737!     MPI code but invalid gfortran.
1738!
1739!***********************************************************************
1740    use com_mod
1741    use unc_mod
1742    use par_mod
1743
1744    implicit none
1745
1746    integer :: grid_size2d,grid_size3d
1747
1748!**********************************************************************
1749
1750    grid_size3d=numxgridn*numygridn*numzgrid*maxspec* &
1751         & maxpointspec_act*nclassunc*maxageclass
1752    grid_size2d=numxgridn*numygridn*maxspec* &
1753         & maxpointspec_act*nclassunc*maxageclass
1754
1755
1756! Time for MPI communications
1757    if (mp_measure_time) call mpif_mtime('commtime',0)
1758
1759! Using a separate grid (gridunc0) for received values, for debugging
1760! call MPI_Reduce(griduncn, griduncn0, grid_size3d, mp_sp, MPI_SUM, id_root, &
1761!      & mp_comm_used, mp_ierr)
1762! if (mp_ierr /= 0) goto 600
1763
1764! Using in-place reduction
1765    if (lroot) then
1766      call MPI_Reduce(MPI_IN_PLACE, griduncn, grid_size3d, mp_sp, MPI_SUM, id_root, &
1767           & mp_comm_used, mp_ierr)
1768      if (mp_ierr /= 0) goto 600
1769    else
1770      call MPI_Reduce(griduncn, griduncn, grid_size3d, mp_sp, MPI_SUM, id_root, &
1771           & mp_comm_used, mp_ierr)
1772    end if
1773
1774    if ((WETDEP).and.(ldirect.gt.0)) then
1775      call MPI_Reduce(wetgriduncn, wetgriduncn0, grid_size2d, mp_cp, MPI_SUM, id_root, &
1776           & mp_comm_used, mp_ierr)
1777      if (mp_ierr /= 0) goto 600
1778    end if
1779
1780    if ((DRYDEP).and.(ldirect.gt.0)) then
1781      call MPI_Reduce(drygriduncn, drygriduncn0, grid_size2d, mp_cp, MPI_SUM, id_root, &
1782           & mp_comm_used, mp_ierr)
1783      if (mp_ierr /= 0) goto 600
1784    end if
1785
1786    if (mp_measure_time) call mpif_mtime('commtime',1)
1787
1788    goto 601
1789
1790600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
1791    stop
1792
1793601 end subroutine mpif_tm_reduce_grid_nest
1794
1795
1796  subroutine mpif_mtime(ident,imode)
1797!***********************************************************************
1798! Measure CPU/Wall time in various parts of the code
1799!
1800! VARIABLES
1801!   ident        character, identifies routine to measure
1802!   imode        integer, 0:start clock(s)  1: stop clock(s)
1803!
1804!***********************************************************************
1805    implicit none
1806
1807    character(LEN=*), intent(in) :: ident
1808    integer, intent(in) :: imode
1809
1810!***********************************************************************
1811
1812    select case(ident)
1813
1814    case ('timemanager')
1815      if (imode.eq.0) then
1816        call cpu_time(tm_tot_beg)
1817      else
1818        call cpu_time(tm_tot_end)
1819        tm_tot_total = tm_tot_total + (tm_tot_end-tm_tot_beg)
1820      end if
1821
1822    case ('wetdepo')
1823      if (imode.eq.0) then
1824        mp_wetdepo_wtime_beg = mpi_wtime()
1825        call cpu_time(mp_wetdepo_time_beg)
1826      else
1827        mp_wetdepo_wtime_end = mpi_wtime()
1828        call cpu_time(mp_wetdepo_time_end)
1829
1830        mp_wetdepo_wtime_total = mp_wetdepo_wtime_total + &
1831             &(mp_wetdepo_wtime_end - mp_wetdepo_wtime_beg)
1832        mp_wetdepo_time_total = mp_wetdepo_time_total + &
1833             &(mp_wetdepo_time_end - mp_wetdepo_time_beg)
1834      end if
1835
1836    case ('advance')
1837      if (imode.eq.0) then
1838        mp_advance_wtime_beg = mpi_wtime()
1839      else
1840        mp_advance_wtime_end = mpi_wtime()
1841
1842        mp_advance_wtime_total = mp_advance_wtime_total + &
1843             &(mp_advance_wtime_end - mp_advance_wtime_beg)
1844      end if
1845
1846    case ('getfields')
1847      if (imode.eq.0) then
1848        mp_getfields_wtime_beg = mpi_wtime()
1849        call cpu_time(mp_getfields_time_beg)
1850      else
1851        mp_getfields_wtime_end = mpi_wtime()
1852        call cpu_time(mp_getfields_time_end)
1853
1854        mp_getfields_wtime_total = mp_getfields_wtime_total + &
1855             &(mp_getfields_wtime_end - mp_getfields_wtime_beg)
1856        mp_getfields_time_total = mp_getfields_time_total + &
1857             &(mp_getfields_time_end - mp_getfields_time_beg)
1858      end if
1859
1860    case ('partloop1')
1861      if (imode.eq.0) then
1862        call cpu_time(tm_nploop_beg)
1863      else
1864        call cpu_time(tm_nploop_end)
1865        tm_nploop_total = tm_nploop_total + (tm_nploop_end - tm_nploop_beg)
1866      end if
1867
1868    case ('conccalc')
1869      if (imode.eq.0) then
1870        call cpu_time(mp_conccalc_time_beg)
1871      else
1872        call cpu_time(mp_conccalc_time_end)
1873        mp_conccalc_time_total = mp_conccalc_time_total + mp_conccalc_time_end - &
1874             &mp_conccalc_time_beg
1875      end if
1876
1877    case ('rootonly')
1878      if (imode.eq.0) then
1879        call cpu_time(mp_root_time_beg)
1880        mp_root_wtime_beg = mpi_wtime()
1881      else
1882        call cpu_time(mp_root_time_end)
1883        mp_root_wtime_end = mpi_wtime()
1884        mp_root_time_total = mp_root_time_total + &
1885             &(mp_root_time_end - mp_root_time_beg)
1886        mp_root_wtime_total = mp_root_wtime_total + &
1887             &(mp_root_wtime_end - mp_root_wtime_beg)
1888      end if
1889
1890    case ('iotime')
1891      if (imode.eq.0) then
1892        mp_io_wtime_beg = mpi_wtime()
1893        call cpu_time(mp_io_time_beg)
1894      else
1895        mp_io_wtime_end = mpi_wtime()
1896        call cpu_time(mp_io_time_end)
1897
1898        mp_io_wtime_total = mp_io_wtime_total + (mp_io_wtime_end - &
1899             & mp_io_wtime_beg)
1900        mp_io_time_total = mp_io_time_total + (mp_io_time_end - &
1901             & mp_io_time_beg)
1902      end if
1903
1904    case ('verttransform')
1905      if (imode.eq.0) then
1906        mp_vt_wtime_beg = mpi_wtime()
1907        call cpu_time(mp_vt_time_beg)
1908      else
1909        mp_vt_wtime_end = mpi_wtime()
1910        call cpu_time(mp_vt_time_end)
1911
1912        mp_vt_wtime_total = mp_vt_wtime_total + (mp_vt_wtime_end - &
1913             & mp_vt_wtime_beg)
1914        mp_vt_time_total = mp_vt_time_total + (mp_vt_time_end - &
1915             & mp_vt_time_beg)
1916      end if
1917
1918    case ('readwind')
1919      if (imode.eq.0) then
1920        call cpu_time(mp_readwind_time_beg)
1921        mp_readwind_wtime_beg = mpi_wtime()
1922      else
1923        call cpu_time(mp_readwind_time_end)
1924        mp_readwind_wtime_end = mpi_wtime()
1925
1926        mp_readwind_time_total = mp_readwind_time_total + &
1927             &(mp_readwind_time_end - mp_readwind_time_beg)
1928        mp_readwind_wtime_total = mp_readwind_wtime_total + &
1929             &(mp_readwind_wtime_end - mp_readwind_wtime_beg)
1930      end if
1931
1932    case ('commtime')
1933      if (imode.eq.0) then
1934        call cpu_time(mp_comm_time_beg)
1935        mp_comm_wtime_beg = mpi_wtime()
1936      else
1937        call cpu_time(mp_comm_time_end)
1938        mp_comm_wtime_end = mpi_wtime()
1939        mp_comm_time_total = mp_comm_time_total + &
1940             &(mp_comm_time_end - mp_comm_time_beg)
1941        mp_comm_wtime_total = mp_comm_wtime_total + &
1942             &(mp_comm_wtime_end - mp_comm_wtime_beg)
1943      end if
1944
1945    case ('flexpart')
1946      if (imode.eq.0) then
1947        mp_total_wtime_beg=mpi_wtime()
1948      else
1949        mp_total_wtime_end=mpi_wtime()
1950        mp_total_wtime_total = mp_total_wtime_end-mp_total_wtime_beg
1951      end if
1952
1953    case default
1954      write(*,*) 'mpi_mod::mpif_mtime> WARNING: unknown case identifier'
1955
1956    end select
1957
1958
1959  end subroutine mpif_mtime
1960
1961
1962  subroutine mpif_finalize
1963!***********************************************************************
1964! Finalize MPI
1965! Optionally print detailed time diagnostics
1966!***********************************************************************
1967    implicit none
1968
1969    integer :: ip !,j,r
1970
1971!***********************************************************************
1972
1973    IF (mp_measure_time) THEN
1974      do ip=0, mp_np-1
1975        call MPI_BARRIER(MPI_COMM_WORLD, mp_ierr)
1976
1977        if (ip == mp_pid) then
1978          write(*,FMT='(72("#"))')
1979          write(*,FMT='(A60,I3)') 'STATISTICS FOR MPI PROCESS:', mp_pid
1980          if (ip == 0) then
1981            write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR ROOT-ONLY (PID 0) &
1982                 &SECTIONS: ', mp_root_time_total
1983            write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR ROOT-ONLY (PID 0) &
1984                 &SECTIONS:', mp_root_wtime_total
1985          end if
1986          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR FLEXPART: ' &
1987               &, mp_total_wtime_total
1988          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR MPI COMMUNICATIONS:',&
1989               & mp_comm_time_total
1990          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR MPI COMMUNICATIONS:',&
1991               & mp_comm_wtime_total
1992          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR MPI BARRIERS:',&
1993               & mp_barrier_time_total
1994          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR MPI BARRIERS:',&
1995               & mp_barrier_wtime_total
1996          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR TIMEMANAGER:',&
1997               & tm_tot_total
1998          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR TIMEMANAGER/NUMPART LOOP:',&
1999               & tm_nploop_total
2000          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR ADVANCE:',&
2001               & mp_advance_wtime_total
2002          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR GETFIELDS:',&
2003               & mp_getfields_wtime_total
2004          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR GETFIELDS:',&
2005               & mp_getfields_time_total
2006          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR READWIND:',&
2007               & mp_readwind_wtime_total
2008          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR READWIND:',&
2009               & mp_readwind_time_total
2010          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR FILE IO:',&
2011               & mp_io_wtime_total
2012          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR FILE IO:',&
2013               & mp_io_time_total
2014          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR WETDEPO:',&
2015               & mp_wetdepo_wtime_total
2016          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR WETDEPO:',&
2017               & mp_wetdepo_time_total
2018          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR CONCCALC:',&
2019               & mp_conccalc_time_total
2020          ! write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR VERTTRANSFORM:',&
2021          !      & mp_vt_wtime_total
2022          ! write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR VERTTRANSFORM:',&
2023          !      & mp_vt_time_total
2024! NB: the 'flush' function is possibly a gfortran-specific extension
2025          call flush()
2026        end if
2027      end do
2028    end if
2029
2030! This call to barrier is for correctly formatting output
2031    call MPI_BARRIER(MPI_COMM_WORLD, mp_ierr)
2032
2033    if (lroot.and.mp_measure_time) then
2034      write(*,FMT='(72("#"))')
2035      WRITE(*,*) "To turn off output of time measurements, set "
2036      WRITE(*,*) "    mp_measure_time=.false."
2037      WRITE(*,*) "in file mpi_mod.f90"
2038      write(*,FMT='(72("#"))')
2039    end if
2040
2041! j=mp_pid*nvar_async
2042! In the implementation with 3 fields, the processes may have posted
2043! MPI_Irecv requests that should be cancelled here
2044! if (.not.lmp_sync) then
2045!   r=mp_pid*nvar_async
2046!   do j=r,r+nvar_async-1
2047!     call MPI_Cancel(j,mp_ierr)
2048!     if (mp_ierr /= 0) write(*,*) '#### mpif_finalize::MPI_Cancel> ERROR ####'
2049!   end do
2050! end if
2051
2052    call MPI_FINALIZE(mp_ierr)
2053    if (mp_ierr /= 0) then
2054      write(*,*) '#### mpif_finalize::MPI_FINALIZE> MPI ERROR ', mp_ierr, ' ####'
2055      stop
2056    end if
2057
2058
2059    end subroutine mpif_finalize
2060
2061
2062    subroutine get_lun(my_lun)
2063!***********************************************************************
2064! get_lun:
2065!   Starting from 100, get next free logical unit number
2066!***********************************************************************
2067
2068      implicit none
2069
2070      integer, intent(inout) :: my_lun
2071      integer, save :: free_lun=100
2072      logical :: exists, iopen
2073
2074!***********************************************************************
2075
2076      loop1: do
2077        inquire(UNIT=free_lun, EXIST=exists, OPENED=iopen)
2078        if (exists .and. .not.iopen) exit loop1
2079        free_lun = free_lun+1
2080      end do loop1
2081      my_lun = free_lun
2082
2083    end subroutine get_lun
2084
2085
2086    subroutine write_data_dbg(array_in, array_name, tstep, ident)
2087!***********************************************************************
2088! Write one-dimensional arrays to file (for debugging purposes)
2089!***********************************************************************
2090      implicit none
2091
2092      real, intent(in), dimension(:) :: array_in
2093      integer, intent(in) :: tstep
2094      integer :: lios
2095      character(LEN=*), intent(in) :: ident, array_name
2096
2097      character(LEN=8) :: c_ts
2098      character(LEN=40) :: fn_1, fn_2
2099
2100!***********************************************************************
2101
2102      write(c_ts, FMT='(I8.8,BZ)') tstep
2103      fn_1='-'//trim(adjustl(c_ts))//'-'//trim(ident)
2104      write(c_ts, FMT='(I2.2,BZ)') mp_np
2105      fn_2= trim(adjustl(array_name))//trim(adjustl(fn_1))//'-np'//trim(adjustl(c_ts))//'.dat'
2106
2107      call get_lun(dat_lun)
2108      open(UNIT=dat_lun, FILE=fn_2, IOSTAT=lios, ACTION='WRITE', &
2109           FORM='UNFORMATTED', STATUS='REPLACE')
2110      write(UNIT=dat_lun, IOSTAT=lios) array_in
2111      close(UNIT=dat_lun)
2112
2113    end subroutine write_data_dbg
2114
2115
2116  subroutine set_fields_synthetic()
2117!*******************************************************************************
2118! DESCRIPTION
2119!   Set all meteorological fields to synthetic (usually constant/homogeneous)
2120!   values.
2121!   Used for validation and error-checking
2122!
2123! NOTE
2124!   This version uses asynchronious communications.
2125!
2126! VARIABLES
2127!   
2128!
2129!
2130!*******************************************************************************
2131    use com_mod
2132
2133    implicit none
2134
2135    integer :: li=1, ui=2 ! wfmem indices (i.e, operate on entire field)
2136
2137    if (.not.lmp_sync) ui=3
2138
2139
2140! Variables transferred at initialization only
2141!*********************************************
2142!    readclouds=readclouds_
2143    oro(:,:)=0.0
2144    excessoro(:,:)=0.0
2145    lsm(:,:)=0.0
2146    xlanduse(:,:,:)=0.0
2147!    wftime
2148!    numbwf
2149!    nmixz
2150!    height
2151
2152! Time-varying fields:
2153    uu(:,:,:,li:ui) = 10.0
2154    vv(:,:,:,li:ui) = 0.0
2155    uupol(:,:,:,li:ui) = 10.0
2156    vvpol(:,:,:,li:ui)=0.0
2157    ww(:,:,:,li:ui)=0.
2158    tt(:,:,:,li:ui)=300.
2159    rho(:,:,:,li:ui)=1.3
2160    drhodz(:,:,:,li:ui)=.0
2161    tth(:,:,:,li:ui)=0.0
2162    qvh(:,:,:,li:ui)=1.0
2163    qv(:,:,:,li:ui)=1.0
2164
2165    pv(:,:,:,li:ui)=1.0
2166    clouds(:,:,:,li:ui)=0
2167
2168    clwc(:,:,:,li:ui)=0.0
2169    ciwc(:,:,:,li:ui)=0.0
2170 
2171! 2D fields
2172
2173    cloudsh(:,:,li:ui)=0
2174    vdep(:,:,:,li:ui)=0.0
2175    ps(:,:,:,li:ui)=1.0e5
2176    sd(:,:,:,li:ui)=0.0
2177    tcc(:,:,:,li:ui)=0.0
2178    tt2(:,:,:,li:ui)=300.
2179    td2(:,:,:,li:ui)=250.
2180    lsprec(:,:,:,li:ui)=0.0
2181    convprec(:,:,:,li:ui)=0.0
2182    ustar(:,:,:,li:ui)=1.0
2183    wstar(:,:,:,li:ui)=1.0
2184    hmix(:,:,:,li:ui)=10000.
2185    tropopause(:,:,:,li:ui)=10000.
2186    oli(:,:,:,li:ui)=0.01
2187 
2188  end subroutine set_fields_synthetic
2189
2190end module mpi_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG