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