source: flexpart.git/src/mpi_mod.f90 @ d6a0977

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

Updates to Henrik's wet depo scheme

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