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

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

Code for cloud water should be correct if the total cw is stored in field clwc (old scheme) or in field qc (new scheme). Minor edits in some files.

  • Property mode set to 100644
File size: 79.3 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_sp,MPI_IN_PLACE,&
480           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
481      if (mp_ierr /= 0) goto 600
482      call MPI_SCATTER(ucp,numpart_mpi,mp_sp,MPI_IN_PLACE,&
483           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
484      if (mp_ierr /= 0) goto 600
485      call MPI_SCATTER(uzp,numpart_mpi,mp_sp,MPI_IN_PLACE,&
486           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
487      if (mp_ierr /= 0) goto 600
488
489      call MPI_SCATTER(us,numpart_mpi,mp_sp,MPI_IN_PLACE,&
490           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
491      if (mp_ierr /= 0) goto 600
492      call MPI_SCATTER(vs,numpart_mpi,mp_sp,MPI_IN_PLACE,&
493           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
494      if (mp_ierr /= 0) goto 600
495      call MPI_SCATTER(ws,numpart_mpi,mp_sp,MPI_IN_PLACE,&
496           &numpart_mpi,mp_sp,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_sp,MPI_IN_PLACE,&
506           &numpart_mpi,mp_sp,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_sp,MPI_IN_PLACE,&
511             &numpart_mpi,mp_sp,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_sp,uap,&
540           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
541      if (mp_ierr /= 0) goto 600
542      call MPI_SCATTER(ucp,numpart_mpi,mp_sp,ucp,&
543           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
544      if (mp_ierr /= 0) goto 600
545      call MPI_SCATTER(uzp,numpart_mpi,mp_sp,uzp,&
546           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
547      if (mp_ierr /= 0) goto 600
548
549      call MPI_SCATTER(us,numpart_mpi,mp_sp,us,&
550           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
551      if (mp_ierr /= 0) goto 600
552      call MPI_SCATTER(vs,numpart_mpi,mp_sp,vs,&
553           &numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
554      if (mp_ierr /= 0) goto 600
555      call MPI_SCATTER(ws,numpart_mpi,mp_sp,ws,&
556           &numpart_mpi,mp_sp,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_sp,ztra1,&
566           &numpart_mpi,mp_sp,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_sp,xmass1(:,i),&
571             &numpart_mpi,mp_sp,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_sp, uap, &
660           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
661      if (mp_ierr /= 0) goto 600
662      call MPI_GATHER(ucp, numpart_mpi, mp_sp, ucp, &
663           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
664      if (mp_ierr /= 0) goto 600
665      call MPI_GATHER(uzp, numpart_mpi, mp_sp, uzp, &
666           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
667      if (mp_ierr /= 0) goto 600
668
669      call MPI_GATHER(us, numpart_mpi, mp_sp, us, &
670           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
671      if (mp_ierr /= 0) goto 600
672      call MPI_GATHER(vs, numpart_mpi, mp_sp, vs, &
673           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
674      if (mp_ierr /= 0) goto 600
675      call MPI_GATHER(ws, numpart_mpi, mp_sp, ws, &
676           &numpart_mpi, mp_sp, 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_sp, ztra1, &
687           &numpart_mpi, mp_sp, 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_sp, xmass1(:,i), &
692             &numpart_mpi, mp_sp, 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_sp, uap, &
725             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
726        if (mp_ierr /= 0) goto 600
727        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ucp, &
728             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
729        if (mp_ierr /= 0) goto 600
730        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, uzp, &
731             &numpart_mpi, mp_sp, 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_sp, us, &
735             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
736        if (mp_ierr /= 0) goto 600
737        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, vs, &
738             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
739        if (mp_ierr /= 0) goto 600
740        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ws, &
741             &numpart_mpi, mp_sp, 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_sp, ztra1, &
752             &numpart_mpi, mp_sp, 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_sp, xmass1(:,i), &
757               &numpart_mpi, mp_sp, 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_sp, uap, &
787             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
788        if (mp_ierr /= 0) goto 600
789        call MPI_GATHER(ucp, numpart_mpi, mp_sp, ucp, &
790             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
791        if (mp_ierr /= 0) goto 600
792        call MPI_GATHER(uzp, numpart_mpi, mp_sp, uzp, &
793             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
794        if (mp_ierr /= 0) goto 600
795
796        call MPI_GATHER(us, numpart_mpi, mp_sp, us, &
797             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
798        if (mp_ierr /= 0) goto 600
799        call MPI_GATHER(vs, numpart_mpi, mp_sp, vs, &
800             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
801        if (mp_ierr /= 0) goto 600
802        call MPI_GATHER(ws, numpart_mpi, mp_sp, ws, &
803             &numpart_mpi, mp_sp, 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_sp, ztra1, &
814             &numpart_mpi, mp_sp, 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_sp, xmass1(:,i), &
819               &numpart_mpi, mp_sp, 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 cloud water 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_sp,id_read,MPI_COMM_WORLD,mp_ierr)
917      if (mp_ierr /= 0) goto 600
918      call MPI_Bcast(excessoro(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
919      if (mp_ierr /= 0) goto 600
920      call MPI_Bcast(lsm(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
921      if (mp_ierr /= 0) goto 600
922      call MPI_Bcast(xlanduse(:,:,:),d2_size3*numclass,mp_sp,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_sp,id_read,MPI_COMM_WORLD,mp_ierr)
937    if (mp_ierr /= 0) goto 600
938    call MPI_Bcast(vv(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
939    if (mp_ierr /= 0) goto 600
940    call MPI_Bcast(uupol(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
941    if (mp_ierr /= 0) goto 600
942    call MPI_Bcast(vvpol(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
943    if (mp_ierr /= 0) goto 600
944    call MPI_Bcast(ww(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
945    if (mp_ierr /= 0) goto 600
946    call MPI_Bcast(tt(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
947    if (mp_ierr /= 0) goto 600
948    call MPI_Bcast(rho(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
949    if (mp_ierr /= 0) goto 600
950    call MPI_Bcast(drhodz(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
951    if (mp_ierr /= 0) goto 600
952    call MPI_Bcast(tth(:,:,:,li:ui),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
953    if (mp_ierr /= 0) goto 600
954    call MPI_Bcast(qvh(:,:,:,li:ui),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
955    if (mp_ierr /= 0) goto 600
956    call MPI_Bcast(qv(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
957    if (mp_ierr /= 0) goto 600
958    call MPI_Bcast(pv(:,:,:,li:ui),d3s1,mp_sp,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),d2s1*5,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
966      ! if (mp_ierr /= 0) goto 600
967      call MPI_Bcast(clw4(:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
968      if (mp_ierr /= 0) goto 600
969      ! call MPI_Bcast(clwc(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
970      ! if (mp_ierr /= 0) goto 600
971      ! call MPI_Bcast(ciwc(:,:,:,li:ui),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
972      ! if (mp_ierr /= 0) goto 600
973    end if
974
975! 2D fields
976    call MPI_Bcast(cloudsh(:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
977    if (mp_ierr /= 0) goto 600
978    call MPI_Bcast(vdep(:,:,:,li:ui),d2s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
979    if (mp_ierr /= 0) goto 600
980    call MPI_Bcast(ps(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
981    if (mp_ierr /= 0) goto 600
982    call MPI_Bcast(sd(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
983    if (mp_ierr /= 0) goto 600
984    call MPI_Bcast(tcc(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
985    if (mp_ierr /= 0) goto 600
986    call MPI_Bcast(tt2(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
987    if (mp_ierr /= 0) goto 600
988    call MPI_Bcast(td2(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
989    if (mp_ierr /= 0) goto 600
990    call MPI_Bcast(lsprec(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
991    if (mp_ierr /= 0) goto 600
992    call MPI_Bcast(convprec(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
993    if (mp_ierr /= 0) goto 600
994    call MPI_Bcast(ustar(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
995    if (mp_ierr /= 0) goto 600
996    call MPI_Bcast(wstar(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
997    if (mp_ierr /= 0) goto 600
998    call MPI_Bcast(hmix(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
999    if (mp_ierr /= 0) goto 600
1000    call MPI_Bcast(tropopause(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1001    if (mp_ierr /= 0) goto 600
1002    call MPI_Bcast(oli(:,:,:,li:ui),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1003    if (mp_ierr /= 0) goto 600
1004
1005    call MPI_Bcast(z0,numclass,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1006    if (mp_ierr /= 0) goto 600
1007
1008    if (mp_measure_time) call mpif_mtime('commtime',1)
1009
1010    goto 601
1011
1012600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
1013    stop
1014
1015601 end subroutine mpif_gf_send_vars
1016
1017
1018  subroutine mpif_gf_send_vars_nest(memstat)
1019!***********************************************************************
1020! DESCRIPTION
1021!   Distribute 'getfield' variables from reader process to all processes.
1022!   For nested fields
1023!
1024!   Called from timemanager
1025!
1026! NOTE
1027!   This subroutine distributes nested windfields from the reader process to
1028!   all other processes. Usually one set of fields is transfered, but at
1029!   first timestep when there are no fields in memory, both are transfered.
1030!   MPI_Bcast is used, so implicitly all processes are synchronized at this
1031!   step
1032!
1033! TODO
1034!   Transfer cloud water/ice if and when available for nested
1035!
1036!***********************************************************************
1037    use com_mod
1038    use par_mod, only: numwfmem
1039
1040    implicit none
1041
1042    integer, intent(in) :: memstat
1043    integer :: i
1044! li,ui    lower,upper indices of fields to be transfered (1-3, ui>=li)
1045    integer :: li,ui
1046
1047! Common array sizes used for communications
1048    integer :: d3_size1 = nxmaxn*nymaxn*nzmax*maxnests
1049    integer :: d3_size2 = nxmaxn*nymaxn*nuvzmax*maxnests
1050    integer :: d2_size1 = nxmaxn*nymaxn*maxnests
1051    integer :: d2_size2 = nxmaxn*nymaxn*maxspec*maxnests
1052    integer :: d2_size3 = nxmaxn*nymaxn*maxnests
1053
1054    integer :: d3s1,d3s2,d2s1,d2s2
1055
1056! First time routine is called the unchangeable fields will be transfered   
1057    logical :: first_call=.true.
1058
1059!**********************************************************************
1060
1061! Sizes of arrays transfered
1062    d3s1=d3_size1
1063    d3s2=d3_size2
1064    d2s1=d2_size1
1065    d2s2=d2_size2
1066
1067! Decide which field(s) need to be transfered
1068    if (memstat.ge.32) then ! distribute 2 fields
1069      li=1; ui=2
1070      d3s1=2*d3_size1
1071      d3s2=2*d3_size2
1072      d2s1=2*d2_size1
1073      d2s2=2*d2_size2
1074    else if (memstat.eq.17) then ! distribute 1 field, on 1st index
1075      li=1; ui=1
1076    else if (memstat.eq.18) then ! distribute 1 field, on 2nd index
1077      li=2; ui=2
1078    else if (memstat.eq.19) then ! distribute 1 field, on 3nd index
1079      li=3; ui=3
1080    else
1081      write(*,*) '#### mpi_mod::mpif_gf_send_vars_nest> ERROR: ', &
1082           & 'wrong value of memstat, exiting ####', memstat
1083      stop
1084    end if
1085
1086
1087! Time for MPI communications
1088    if (mp_measure_time) call mpif_mtime('commtime',0)
1089
1090! Distribute variables, send from getfield process (id_read) to other
1091! processes
1092!**********************************************************************
1093
1094! Static fields/variables sent only at startup
1095    if (first_call) then
1096      call MPI_Bcast(oron(:,:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1097      if (mp_ierr /= 0) goto 600
1098      call MPI_Bcast(excessoron(:,:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1099      if (mp_ierr /= 0) goto 600
1100      call MPI_Bcast(lsmn(:,:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1101      if (mp_ierr /= 0) goto 600
1102      call MPI_Bcast(xlandusen(:,:,:,:),d2_size3*numclass,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1103      if (mp_ierr /= 0) goto 600
1104      first_call=.false.
1105    end if
1106
1107! MPI prefers contiguous arrays for sending (else a buffer is created),
1108! hence the loop
1109
1110    do i=1, numbnests
1111! 3D fields
1112      call MPI_Bcast(uun(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1113      if (mp_ierr /= 0) goto 600
1114      call MPI_Bcast(vvn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1115      if (mp_ierr /= 0) goto 600
1116      call MPI_Bcast(wwn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1117      if (mp_ierr /= 0) goto 600
1118      call MPI_Bcast(ttn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1119      if (mp_ierr /= 0) goto 600
1120      call MPI_Bcast(rhon(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1121      if (mp_ierr /= 0) goto 600
1122      call MPI_Bcast(drhodzn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1123      if (mp_ierr /= 0) goto 600
1124      call MPI_Bcast(tthn(:,:,:,li:ui,i),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1125      if (mp_ierr /= 0) goto 600
1126      call MPI_Bcast(qvhn(:,:,:,li:ui,i),d3s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1127      if (mp_ierr /= 0) goto 600
1128      call MPI_Bcast(qvn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1129      if (mp_ierr /= 0) goto 600
1130      call MPI_Bcast(pvn(:,:,:,li:ui,i),d3s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1131      if (mp_ierr /= 0) goto 600
1132      call MPI_Bcast(cloudsn(:,:,:,li:ui,i),d3s1,MPI_INTEGER1,id_read,MPI_COMM_WORLD,mp_ierr)
1133      if (mp_ierr /= 0) goto 600
1134
1135! 2D fields
1136      call MPI_Bcast(cloudsnh(:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1137      if (mp_ierr /= 0) goto 600
1138      call MPI_Bcast(vdepn(:,:,:,li:ui,i),d2s2,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1139      if (mp_ierr /= 0) goto 600
1140      call MPI_Bcast(psn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1141      if (mp_ierr /= 0) goto 600
1142      call MPI_Bcast(sdn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1143      if (mp_ierr /= 0) goto 600
1144      call MPI_Bcast(tccn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1145      if (mp_ierr /= 0) goto 600
1146      call MPI_Bcast(tt2n(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1147      if (mp_ierr /= 0) goto 600
1148      call MPI_Bcast(td2n(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1149      if (mp_ierr /= 0) goto 600
1150      call MPI_Bcast(lsprecn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1151      if (mp_ierr /= 0) goto 600
1152      call MPI_Bcast(convprecn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1153      if (mp_ierr /= 0) goto 600
1154      call MPI_Bcast(ustarn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1155      if (mp_ierr /= 0) goto 600
1156      call MPI_Bcast(wstarn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1157      if (mp_ierr /= 0) goto 600
1158      call MPI_Bcast(olin(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1159      if (mp_ierr /= 0) goto 600
1160      call MPI_Bcast(hmixn(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1161      if (mp_ierr /= 0) goto 600
1162      call MPI_Bcast(tropopausen(:,:,:,li:ui,i),d2s1,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1163      if (mp_ierr /= 0) goto 600
1164    end do
1165
1166    if (mp_measure_time) call mpif_mtime('commtime',1)
1167
1168    goto 601
1169
1170600 write(*,*) "mpi_mod> mp_ierr \= 0",mp_ierr
1171    stop
1172
1173601 end subroutine mpif_gf_send_vars_nest
1174
1175
1176  subroutine mpif_gf_send_vars_async(memstat)
1177!*******************************************************************************
1178! DESCRIPTION
1179!   Distribute 'getfield' variables from reader process to all processes.
1180!   Called from timemanager by reader process only
1181!
1182! NOTE
1183!   This version uses asynchronious sends. The newest fields are sent in the
1184!   background, so calculations can continue while
1185!   MPI communications are performed.
1186!
1187!   The way the MPI tags/requests are sequenced means that this routine must
1188!   carefully match routine 'mpif_gf_recv_vars_async'
1189!
1190! VARIABLES
1191!   memstat -- input, for resolving pointer to windfield index being read
1192!   mind    -- index where to place new fields
1193!
1194! TODO
1195!
1196!*******************************************************************************
1197    use com_mod
1198
1199    implicit none
1200
1201    integer, intent(in) :: memstat
1202    integer :: mind
1203    integer :: dest,i
1204
1205! Common array sizes used for communications
1206    integer :: d3s1 = nxmax*nymax*nzmax
1207    integer :: d3s2 = nxmax*nymax*nuvzmax
1208    integer :: d2s1 = nxmax*nymax
1209    integer :: d2s2 = nxmax*nymax*maxspec
1210    integer :: d1s1 = maxwf
1211
1212!*******************************************************************************
1213
1214! At the time the send is posted, the reader process is one step further
1215! in the permutation of memstat compared with the receiving processes
1216
1217    if (memstat.ge.32) then
1218! last read was synchronous, to indices 1 and 2, 3 is free
1219      write(*,*) "#### mpi_mod::mpif_gf_send_vars_async> ERROR: &
1220           & memstat>=32 should never happen here."
1221      stop
1222    else if (memstat.eq.17) then
1223! old fields on 1,2, send 3
1224      mind=3
1225    else if (memstat.eq.18) then
1226! old fields on 2,3, send 1
1227      mind=1
1228    else if (memstat.eq.19) then
1229! old fields on 3,1, send 2
1230      mind=2
1231    else
1232      write(*,*) "#### mpi_mod::mpif_gf_send_vars_async> ERROR: &
1233           & invalid memstat"
1234    end if
1235
1236    if (mp_dev_mode) then
1237      if (mp_pid.ne.id_read) then
1238        write(*,*) 'MPI_DEV_MODE: error calling mpif_gf_send_vars_async'
1239      end if
1240    end if
1241
1242    if (mp_dev_mode) write(*,*) '## in mpif_gf_send_vars_async, memstat:', memstat
1243
1244! Time for MPI communications
1245    if (mp_measure_time) call mpif_mtime('commtime',0)
1246
1247! Loop over receiving processes, initiate data sending
1248!*****************************************************
1249
1250    do dest=0,mp_np-1 ! mp_np-2 will also work if last proc reserved for reading
1251! TODO: use mp_partgroup_np here
1252      if (dest.eq.id_read) cycle
1253      i=dest*nvar_async
1254      call MPI_Isend(uu(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1255      if (mp_ierr /= 0) goto 600
1256      i=i+1
1257      call MPI_Isend(vv(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1258      if (mp_ierr /= 0) goto 600
1259      i=i+1
1260      call MPI_Isend(uupol(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1261      if (mp_ierr /= 0) goto 600
1262      i=i+1
1263      call MPI_Isend(vvpol(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1264      if (mp_ierr /= 0) goto 600
1265      i=i+1
1266      call MPI_Isend(ww(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1267      if (mp_ierr /= 0) goto 600
1268      i=i+1
1269      call MPI_Isend(tt(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1270      if (mp_ierr /= 0) goto 600
1271      i=i+1
1272      call MPI_Isend(rho(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1273      if (mp_ierr /= 0) goto 600
1274      i=i+1
1275      call MPI_Isend(drhodz(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1276      if (mp_ierr /= 0) goto 600
1277      i=i+1
1278      call MPI_Isend(tth(:,:,:,mind),d3s2,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1279      if (mp_ierr /= 0) goto 600
1280      i=i+1
1281      call MPI_Isend(qvh(:,:,:,mind),d3s2,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1282      if (mp_ierr /= 0) goto 600
1283      i=i+1
1284      call MPI_Isend(qv(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1285      if (mp_ierr /= 0) goto 600
1286      i=i+1
1287      call MPI_Isend(pv(:,:,:,mind),d3s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1288      if (mp_ierr /= 0) goto 600
1289      i=i+1
1290      call MPI_Isend(clouds(:,:,:,mind),d3s1,MPI_INTEGER1,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1291      i=i+1
1292      if (mp_ierr /= 0) goto 600
1293
1294      call MPI_Isend(cloudsh(:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1295      if (mp_ierr /= 0) goto 600
1296      i=i+1
1297      call MPI_Isend(vdep(:,:,:,mind),d2s2,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1298      if (mp_ierr /= 0) goto 600
1299      i=i+1
1300      call MPI_Isend(ps(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1301      if (mp_ierr /= 0) goto 600
1302      i=i+1
1303      call MPI_Isend(sd(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1304      if (mp_ierr /= 0) goto 600
1305      i=i+1
1306      call MPI_Isend(tcc(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1307      if (mp_ierr /= 0) goto 600
1308      i=i+1
1309      call MPI_Isend(tt2(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1310      if (mp_ierr /= 0) goto 600
1311      i=i+1
1312      call MPI_Isend(td2(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1313      if (mp_ierr /= 0) goto 600
1314      i=i+1
1315      call MPI_Isend(lsprec(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1316      if (mp_ierr /= 0) goto 600
1317      i=i+1
1318      call MPI_Isend(convprec(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1319      if (mp_ierr /= 0) goto 600
1320      i=i+1
1321      call MPI_Isend(ustar(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1322      if (mp_ierr /= 0) goto 600
1323      i=i+1
1324      call MPI_Isend(wstar(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1325      if (mp_ierr /= 0) goto 600
1326      i=i+1
1327      call MPI_Isend(hmix(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1328      if (mp_ierr /= 0) goto 600
1329      i=i+1
1330      call MPI_Isend(tropopause(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1331      if (mp_ierr /= 0) goto 600
1332      i=i+1
1333      call MPI_Isend(oli(:,:,:,mind),d2s1,mp_sp,dest,tm1,MPI_COMM_WORLD,reqs(i),mp_ierr)
1334      if (mp_ierr /= 0) goto 600
1335
1336! Send cloud water if it exists. Increment counter always (as on receiving end)
1337      if (readclouds) then
1338        i=i+1
1339        ! call MPI_Isend(icloud_stats(:,:,:,mind),d2s1*5,mp_sp,dest,tm1,&
1340        !      &MPI_COMM_WORLD,reqs(i),mp_ierr)
1341        call MPI_Isend(clw4(:,:,mind),d2s1,mp_sp,dest,tm1,&
1342             &MPI_COMM_WORLD,reqs(i),mp_ierr)
1343
1344        if (mp_ierr /= 0) goto 600
1345
1346        ! call MPI_Isend(clwc(:,:,:,mind),d3s1,mp_sp,dest,tm1,&
1347        !      &MPI_COMM_WORLD,reqs(i),mp_ierr)
1348        ! if (mp_ierr /= 0) goto 600
1349        ! i=i+1
1350
1351        ! call MPI_Isend(ciwc(:,:,:,mind),d3s1,mp_sp,dest,tm1,&
1352        !      &MPI_COMM_WORLD,reqs(i),mp_ierr)
1353        ! if (mp_ierr /= 0) goto 600
1354
1355      end if
1356    end do
1357
1358    if (mp_measure_time) call mpif_mtime('commtime',1)
1359
1360    goto 601
1361
1362600 write(*,*) "#### mpi_mod::mpif_gf_send_vars_async> mp_ierr \= 0", mp_ierr
1363    stop
1364
1365601 end subroutine mpif_gf_send_vars_async
1366
1367
1368  subroutine mpif_gf_recv_vars_async(memstat)
1369!*******************************************************************************
1370! DESCRIPTION
1371!   Receive 'getfield' variables from reader process.
1372!   Called from timemanager by all processes except reader
1373!
1374! NOTE
1375!   This version uses asynchronious communications.
1376!
1377! VARIABLES
1378!   memstat -- input, used to resolve windfield index being received
1379!
1380! TODO
1381!
1382!*******************************************************************************
1383    use com_mod
1384
1385    implicit none
1386
1387    integer, intent(in) :: memstat
1388    integer :: mind,j
1389
1390! Common array sizes used for communications
1391    integer :: d3s1 = nxmax*nymax*nzmax
1392    integer :: d3s2 = nxmax*nymax*nuvzmax
1393    integer :: d2s1 = nxmax*nymax
1394    integer :: d2s2 = nxmax*nymax*maxspec
1395    integer :: d1_size1 = maxwf
1396
1397!    integer :: d3s1,d3s2,d2s1,d2s2
1398!*******************************************************************************
1399
1400! At the time this immediate receive is posted, memstat is the state of
1401! windfield indices at the previous/current time. From this, the future
1402! state is deduced.
1403    if (memstat.eq.32) then
1404! last read was synchronous to indices 1 and 2, 3 is free
1405      mind=3
1406    else if (memstat.eq.17) then
1407! last read was asynchronous to index 3, 1 is free
1408      mind=1
1409    else if (memstat.eq.18) then
1410! last read was asynchronous to index 1, 2 is free
1411      mind=2
1412    else if (memstat.eq.19) then
1413! last read was asynchronous to index 2, 3 is free
1414      mind=3
1415    end if
1416
1417    if (mp_dev_mode) then
1418      if (mp_pid.eq.id_read) then
1419        write(*,*) 'MPI_DEV_MODE: error calling mpif_gf_recv_vars_async'
1420      end if
1421    end if
1422
1423! Time for MPI communications
1424    if (mp_measure_time) call mpif_mtime('commtime',0)
1425
1426    if (mp_dev_mode) write(*,*) '## in mpif_gf_send_vars_async, memstat:', memstat
1427
1428! Initiate receiving of data
1429!***************************
1430
1431! Get MPI tags/requests for communications
1432    j=mp_pid*nvar_async
1433    call MPI_Irecv(uu(:,:,:,mind),d3s1,mp_sp,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(vv(:,:,:,mind),d3s1,mp_sp,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(uupol(:,:,:,mind),d3s1,mp_sp,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(vvpol(:,:,:,mind),d3s1,mp_sp,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(ww(:,:,:,mind),d3s1,mp_sp,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(tt(:,:,:,mind),d3s1,mp_sp,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(rho(:,:,:,mind),d3s1,mp_sp,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(drhodz(:,:,:,mind),d3s1,mp_sp,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(tth(:,:,:,mind),d3s2,mp_sp,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    call MPI_Irecv(qvh(:,:,:,mind),d3s2,mp_sp,id_read,MPI_ANY_TAG,&
1470         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1471    if (mp_ierr /= 0) goto 600
1472    j=j+1
1473
1474    call MPI_Irecv(qv(:,:,:,mind),d3s1,mp_sp,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(pv(:,:,:,mind),d3s1,mp_sp,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    call MPI_Irecv(clouds(:,:,:,mind),d3s1,MPI_INTEGER1,id_read,MPI_ANY_TAG,&
1483         &MPI_COMM_WORLD,reqs(j),mp_ierr)   
1484    if (mp_ierr /= 0) goto 600
1485    j=j+1
1486
1487    call MPI_Irecv(cloudsh(:,:,mind),d2s1,mp_sp,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(vdep(:,:,:,mind),d2s2,mp_sp,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(ps(:,:,:,mind),d2s1,mp_sp,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(sd(:,:,:,mind),d2s1,mp_sp,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(tcc(:,:,:,mind),d2s1,mp_sp,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(tt2(:,:,:,mind),d2s1,mp_sp,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(td2(:,:,:,mind),d2s1,mp_sp,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(lsprec(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1516         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1517    if (mp_ierr /= 0) goto 600
1518    j=j+1
1519    call MPI_Irecv(convprec(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1520         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1521    if (mp_ierr /= 0) goto 600
1522
1523    call MPI_Irecv(ustar(:,:,:,mind),d2s1,mp_sp,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(wstar(:,:,:,mind),d2s1,mp_sp,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(hmix(:,:,:,mind),d2s1,mp_sp,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(tropopause(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1536         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1537    if (mp_ierr /= 0) goto 600
1538    j=j+1
1539    call MPI_Irecv(oli(:,:,:,mind),d2s1,mp_sp,id_read,MPI_ANY_TAG,&
1540         &MPI_COMM_WORLD,reqs(j),mp_ierr)
1541    if (mp_ierr /= 0) goto 600
1542
1543
1544! Post request for clwc. These data are possibly not sent, request must then be cancelled
1545! For now assume that data at all steps either have or do not have water
1546    if (readclouds) then
1547      j=j+1
1548
1549      ! call MPI_Irecv(icloud_stats(:,:,:,mind),d2s1*5,mp_sp,id_read,MPI_ANY_TAG,&
1550      !      &MPI_COMM_WORLD,reqs(j),mp_ierr)
1551      call MPI_Irecv(clw4(:,:,mind),d2s1*5,mp_sp,id_read,MPI_ANY_TAG,&
1552           &MPI_COMM_WORLD,reqs(j),mp_ierr)
1553      if (mp_ierr /= 0) goto 600
1554
1555      ! call MPI_Irecv(clwc(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1556      !      &MPI_COMM_WORLD,reqs(j),mp_ierr)   
1557      ! if (mp_ierr /= 0) goto 600
1558      ! j=j+1
1559      ! call MPI_Irecv(ciwc(:,:,:,mind),d3s1,mp_sp,id_read,MPI_ANY_TAG,&
1560      !      &MPI_COMM_WORLD,reqs(j),mp_ierr)   
1561      ! if (mp_ierr /= 0) goto 600
1562
1563    end if
1564
1565
1566    if (mp_measure_time) call mpif_mtime('commtime',1)
1567
1568    goto 601
1569
1570600 write(*,*) "#### mpi_mod::mpif_gf_recv_vars_async> MPI ERROR ", mp_ierr
1571    stop
1572
1573601 end subroutine mpif_gf_recv_vars_async
1574
1575
1576  subroutine mpif_gf_request
1577!*******************************************************************************
1578! DESCRIPTION
1579!   Check for completion of MPI_Isend/MPI_Irecv data transfer.
1580!   
1581!   
1582! NOTE
1583!   implicit synchronisation between all processes takes place here
1584!
1585!
1586!*******************************************************************************
1587    use com_mod, only: readclouds
1588
1589    implicit none
1590
1591
1592    integer :: n_req
1593    integer :: i
1594
1595!***********************************************************************
1596
1597    n_req=nvar_async*mp_np
1598
1599    if (mp_measure_time) call mpif_mtime('commtime',0)
1600
1601!    call MPI_Wait(rm1,MPI_STATUS_IGNORE,mp_ierr)
1602
1603! TODO: cancel recv request if readclouds=.false.
1604!    if (readclouds) then
1605    call MPI_Waitall(n_req,reqs,MPI_STATUSES_IGNORE,mp_ierr)
1606!    endif
1607! else
1608!   do i = 0, nvar_async*mp_np-1
1609!     if (mod(i,27).eq.0 .or. mod(i,28).eq.0) then
1610!       call MPI_Cancel(reqs(i),mp_ierr)
1611!       cycle
1612!     end if
1613!     call MPI_Wait(reqs(i),MPI_STATUS_IGNORE,mp_ierr)
1614!   end do
1615! end if
1616
1617    if (mp_ierr /= 0) goto 600
1618
1619    if (mp_measure_time) call mpif_mtime('commtime',1)
1620
1621    goto 601
1622
1623600 write(*,*) "#### mpi_mod::mpif_gf_request> MPI ERROR ", mp_ierr
1624    stop
1625
1626601 end subroutine mpif_gf_request
1627
1628
1629  subroutine mpif_tm_reduce_grid
1630!***********************************************************************
1631! Collect grid variable to PID 0, adding from all processes.
1632!
1633! NOTE
1634!   - Older MPI implementations may lack the MPI_IN_PLACE operation.
1635!     If so use 1) below and change unc_mod
1636!
1637!***********************************************************************
1638    use com_mod
1639    use unc_mod
1640    use par_mod
1641
1642    implicit none
1643
1644    integer :: grid_size2d,grid_size3d
1645    integer, parameter :: rcpt_size=maxreceptor*maxspec
1646
1647!**********************************************************************
1648    grid_size2d=numxgrid*numygrid*maxspec* &
1649         & maxpointspec_act*nclassunc*maxageclass
1650    grid_size3d=numxgrid*numygrid*numzgrid*maxspec* &
1651         & maxpointspec_act*nclassunc*maxageclass
1652
1653
1654! Time for MPI communications
1655    if (mp_measure_time) call mpif_mtime('commtime',0)
1656
1657! 1) Using a separate grid (gridunc0) for received values
1658! call MPI_Reduce(gridunc, gridunc0, grid_size3d, mp_sp, MPI_SUM, id_root, &
1659!      & mp_comm_used, mp_ierr)
1660! if (mp_ierr /= 0) goto 600
1661
1662! 2) Using in-place reduction
1663    if (lroot) then
1664      call MPI_Reduce(MPI_IN_PLACE, gridunc, grid_size3d, mp_sp, MPI_SUM, id_root, &
1665           & mp_comm_used, mp_ierr)
1666      if (mp_ierr /= 0) goto 600
1667    else
1668      call MPI_Reduce(gridunc, gridunc, grid_size3d, mp_sp, MPI_SUM, id_root, &
1669           & mp_comm_used, mp_ierr)
1670    end if
1671
1672    if ((WETDEP).and.(ldirect.gt.0)) then
1673      call MPI_Reduce(wetgridunc, wetgridunc0, grid_size2d, mp_sp, MPI_SUM, id_root, &
1674           & mp_comm_used, mp_ierr)
1675      if (mp_ierr /= 0) goto 600
1676    end if
1677
1678    if ((DRYDEP).and.(ldirect.gt.0)) then
1679      call MPI_Reduce(drygridunc, drygridunc0, grid_size2d, mp_sp, MPI_SUM, id_root, &
1680           & mp_comm_used, mp_ierr)
1681      if (mp_ierr /= 0) goto 600
1682    end if
1683
1684! Receptor concentrations   
1685    if (lroot) then
1686      call MPI_Reduce(MPI_IN_PLACE,creceptor,rcpt_size,mp_sp,MPI_SUM,id_root, &
1687           & mp_comm_used,mp_ierr)
1688      if (mp_ierr /= 0) goto 600
1689    else
1690      call MPI_Reduce(creceptor,creceptor,rcpt_size,mp_sp,MPI_SUM,id_root, &
1691           & mp_comm_used,mp_ierr)
1692    end if
1693
1694    if (mp_measure_time) call mpif_mtime('commtime',1)
1695
1696    goto 601
1697
1698600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
1699    stop
1700
1701601 end subroutine mpif_tm_reduce_grid
1702
1703
1704  subroutine mpif_tm_reduce_grid_nest
1705!***********************************************************************
1706! Collect nested grids to PID 0, adding from all processes.
1707!
1708! NOTE
1709!   - Compiling with 'fcheck=all' (gfortran) will cause a run-time error
1710!     as wetgriduncn0 is only allocated for root process. Possibly valid
1711!     MPI code but invalid gfortran.
1712!
1713!***********************************************************************
1714    use com_mod
1715    use unc_mod
1716    use par_mod
1717
1718    implicit none
1719
1720    integer :: grid_size2d,grid_size3d
1721
1722!**********************************************************************
1723
1724    grid_size3d=numxgridn*numygridn*numzgrid*maxspec* &
1725         & maxpointspec_act*nclassunc*maxageclass
1726    grid_size2d=numxgridn*numygridn*maxspec* &
1727         & maxpointspec_act*nclassunc*maxageclass
1728
1729
1730! Time for MPI communications
1731    if (mp_measure_time) call mpif_mtime('commtime',0)
1732
1733! Using a separate grid (gridunc0) for received values, for debugging
1734! call MPI_Reduce(griduncn, griduncn0, grid_size3d, mp_sp, MPI_SUM, id_root, &
1735!      & mp_comm_used, mp_ierr)
1736! if (mp_ierr /= 0) goto 600
1737
1738! Using in-place reduction
1739    if (lroot) then
1740      call MPI_Reduce(MPI_IN_PLACE, griduncn, grid_size3d, mp_sp, MPI_SUM, id_root, &
1741           & mp_comm_used, mp_ierr)
1742      if (mp_ierr /= 0) goto 600
1743    else
1744      call MPI_Reduce(griduncn, griduncn, grid_size3d, mp_sp, MPI_SUM, id_root, &
1745           & mp_comm_used, mp_ierr)
1746    end if
1747
1748    if ((WETDEP).and.(ldirect.gt.0)) then
1749      call MPI_Reduce(wetgriduncn, wetgriduncn0, grid_size2d, mp_sp, MPI_SUM, id_root, &
1750           & mp_comm_used, mp_ierr)
1751      if (mp_ierr /= 0) goto 600
1752    end if
1753
1754    if ((DRYDEP).and.(ldirect.gt.0)) then
1755      call MPI_Reduce(drygriduncn, drygriduncn0, grid_size2d, mp_sp, MPI_SUM, id_root, &
1756           & mp_comm_used, mp_ierr)
1757      if (mp_ierr /= 0) goto 600
1758    end if
1759
1760    if (mp_measure_time) call mpif_mtime('commtime',1)
1761
1762    goto 601
1763
1764600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
1765    stop
1766
1767601 end subroutine mpif_tm_reduce_grid_nest
1768
1769
1770  subroutine mpif_mtime(ident,imode)
1771!***********************************************************************
1772! Measure CPU/Wall time in various parts of the code
1773!
1774! VARIABLES
1775!   ident        character, identifies routine to measure
1776!   imode        integer, 0:start clock(s)  1: stop clock(s)
1777!
1778!***********************************************************************
1779    implicit none
1780
1781    character(LEN=*), intent(in) :: ident
1782    integer, intent(in) :: imode
1783
1784!***********************************************************************
1785
1786    select case(ident)
1787
1788    case ('timemanager')
1789      if (imode.eq.0) then
1790        call cpu_time(tm_tot_beg)
1791      else
1792        call cpu_time(tm_tot_end)
1793        tm_tot_total = tm_tot_total + (tm_tot_end-tm_tot_beg)
1794      end if
1795
1796    case ('wetdepo')
1797      if (imode.eq.0) then
1798        mp_wetdepo_wtime_beg = mpi_wtime()
1799        call cpu_time(mp_wetdepo_time_beg)
1800      else
1801        mp_wetdepo_wtime_end = mpi_wtime()
1802        call cpu_time(mp_wetdepo_time_end)
1803
1804        mp_wetdepo_wtime_total = mp_wetdepo_wtime_total + &
1805             &(mp_wetdepo_wtime_end - mp_wetdepo_wtime_beg)
1806        mp_wetdepo_time_total = mp_wetdepo_time_total + &
1807             &(mp_wetdepo_time_end - mp_wetdepo_time_beg)
1808      end if
1809
1810    case ('getfields')
1811      if (imode.eq.0) then
1812        mp_getfields_wtime_beg = mpi_wtime()
1813        call cpu_time(mp_getfields_time_beg)
1814      else
1815        mp_getfields_wtime_end = mpi_wtime()
1816        call cpu_time(mp_getfields_time_end)
1817
1818        mp_getfields_wtime_total = mp_getfields_wtime_total + &
1819             &(mp_getfields_wtime_end - mp_getfields_wtime_beg)
1820        mp_getfields_time_total = mp_getfields_time_total + &
1821             &(mp_getfields_time_end - mp_getfields_time_beg)
1822      end if
1823
1824    case ('partloop1')
1825      if (imode.eq.0) then
1826        call cpu_time(tm_nploop_beg)
1827      else
1828        call cpu_time(tm_nploop_end)
1829        tm_nploop_total = tm_nploop_total + (tm_nploop_end - tm_nploop_beg)
1830      end if
1831
1832    case ('conccalc')
1833      if (imode.eq.0) then
1834        call cpu_time(mp_conccalc_time_beg)
1835      else
1836        call cpu_time(mp_conccalc_time_end)
1837        mp_conccalc_time_total = mp_conccalc_time_total + mp_conccalc_time_end - &
1838             &mp_conccalc_time_beg
1839      end if
1840    case ('rootonly')
1841      if (imode.eq.0) then
1842        call cpu_time(mp_root_time_beg)
1843        mp_root_wtime_beg = mpi_wtime()
1844      else
1845        call cpu_time(mp_root_time_end)
1846        mp_root_wtime_end = mpi_wtime()
1847        mp_root_time_total = mp_root_time_total + &
1848             &(mp_root_time_end - mp_root_time_beg)
1849        mp_root_wtime_total = mp_root_wtime_total + &
1850             &(mp_root_wtime_end - mp_root_wtime_beg)
1851      end if
1852
1853    case ('iotime')
1854      if (imode.eq.0) then
1855        mp_io_wtime_beg = mpi_wtime()
1856        call cpu_time(mp_io_time_beg)
1857      else
1858        mp_io_wtime_end = mpi_wtime()
1859        call cpu_time(mp_io_time_end)
1860
1861        mp_io_wtime_total = mp_io_wtime_total + (mp_io_wtime_end - &
1862             & mp_io_wtime_beg)
1863        mp_io_time_total = mp_io_time_total + (mp_io_time_end - &
1864             & mp_io_time_beg)
1865      end if
1866
1867    case ('verttransform')
1868      if (imode.eq.0) then
1869        mp_vt_wtime_beg = mpi_wtime()
1870        call cpu_time(mp_vt_time_beg)
1871      else
1872        mp_vt_wtime_end = mpi_wtime()
1873        call cpu_time(mp_vt_time_end)
1874
1875        mp_vt_wtime_total = mp_vt_wtime_total + (mp_vt_wtime_end - &
1876             & mp_vt_wtime_beg)
1877        mp_vt_time_total = mp_vt_time_total + (mp_vt_time_end - &
1878             & mp_vt_time_beg)
1879      end if
1880
1881    case ('readwind')
1882      if (imode.eq.0) then
1883        call cpu_time(mp_readwind_time_beg)
1884        mp_readwind_wtime_beg = mpi_wtime()
1885      else
1886        call cpu_time(mp_readwind_time_end)
1887        mp_readwind_wtime_end = mpi_wtime()
1888
1889        mp_readwind_time_total = mp_readwind_time_total + &
1890             &(mp_readwind_time_end - mp_readwind_time_beg)
1891        mp_readwind_wtime_total = mp_readwind_wtime_total + &
1892             &(mp_readwind_wtime_end - mp_readwind_wtime_beg)
1893      end if
1894
1895    case ('commtime')
1896      if (imode.eq.0) then
1897        call cpu_time(mp_comm_time_beg)
1898        mp_comm_wtime_beg = mpi_wtime()
1899      else
1900        call cpu_time(mp_comm_time_end)
1901        mp_comm_wtime_end = mpi_wtime()
1902        mp_comm_time_total = mp_comm_time_total + &
1903             &(mp_comm_time_end - mp_comm_time_beg)
1904        mp_comm_wtime_total = mp_comm_wtime_total + &
1905             &(mp_comm_wtime_end - mp_comm_wtime_beg)
1906      end if
1907
1908    case ('flexpart')
1909      if (imode.eq.0) then
1910        mp_total_wtime_beg=mpi_wtime()
1911      else
1912        mp_total_wtime_end=mpi_wtime()
1913        mp_total_wtime_total = mp_total_wtime_end-mp_total_wtime_beg
1914      end if
1915
1916    case default
1917      write(*,*) 'mpi_mod::mpif_mtime> WARNING: unknown case identifier'
1918
1919    end select
1920
1921
1922  end subroutine mpif_mtime
1923
1924
1925  subroutine mpif_finalize
1926!***********************************************************************
1927! Finalize MPI
1928! Optionally print detailed time diagnostics
1929!***********************************************************************
1930    implicit none
1931
1932    integer :: ip,j,r
1933
1934!***********************************************************************
1935
1936    IF (mp_measure_time) THEN
1937      do ip=0, mp_np-1
1938        call MPI_BARRIER(MPI_COMM_WORLD, mp_ierr)
1939
1940        if (ip == mp_pid) then
1941          write(*,FMT='(72("#"))')
1942          write(*,FMT='(A60,I3)') 'STATISTICS FOR MPI PROCESS:', mp_pid
1943          if (ip == 0) then
1944            write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR ROOT-ONLY (PID 0) &
1945                 &SECTIONS: ', mp_root_time_total
1946            write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR ROOT-ONLY (PID 0) &
1947                 &SECTIONS:', mp_root_wtime_total
1948          end if
1949          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR FLEXPART: ' &
1950               &, mp_total_wtime_total
1951          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR MPI COMMUNICATIONS:',&
1952               & mp_comm_time_total
1953          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR MPI COMMUNICATIONS:',&
1954               & mp_comm_wtime_total
1955          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR MPI BARRIERS:',&
1956               & mp_barrier_time_total
1957          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR MPI BARRIERS:',&
1958               & mp_barrier_wtime_total
1959          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR TIMEMANAGER:',&
1960               & tm_tot_total
1961          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR TIMEMANAGER/NUMPART LOOP:',&
1962               & tm_nploop_total
1963          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR ADVANCE:',&
1964               & mp_advance_wtime_total
1965          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR GETFIELDS:',&
1966               & mp_getfields_wtime_total
1967          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR GETFIELDS:',&
1968               & mp_getfields_time_total
1969          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR READWIND:',&
1970               & mp_readwind_wtime_total
1971          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR READWIND:',&
1972               & mp_readwind_time_total
1973          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR FILE IO:',&
1974               & mp_io_wtime_total
1975          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR FILE IO:',&
1976               & mp_io_time_total
1977          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR WETDEPO:',&
1978               & mp_wetdepo_wtime_total
1979          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR WETDEPO:',&
1980               & mp_wetdepo_time_total
1981          write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR CONCCALC:',&
1982               & mp_conccalc_time_total
1983          ! write(*,FMT='(A60,TR1,F9.2)') 'TOTAL WALL TIME FOR VERTTRANSFORM:',&
1984          !      & mp_vt_wtime_total
1985          ! write(*,FMT='(A60,TR1,F9.2)') 'TOTAL CPU TIME FOR VERTTRANSFORM:',&
1986          !      & mp_vt_time_total
1987! NB: the 'flush' function is possibly a gfortran-specific extension
1988          call flush()
1989        end if
1990      end do
1991    end if
1992
1993! This call to barrier is for correctly formatting output
1994    call MPI_BARRIER(MPI_COMM_WORLD, mp_ierr)
1995
1996    if (lroot.and.mp_measure_time) then
1997      write(*,FMT='(72("#"))')
1998      WRITE(*,*) "To turn off output of time measurements, set "
1999      WRITE(*,*) "    mp_measure_time=.false."
2000      WRITE(*,*) "in file mpi_mod.f90"
2001      write(*,FMT='(72("#"))')
2002    end if
2003
2004! j=mp_pid*nvar_async
2005! In the implementation with 3 fields, the processes may have posted
2006! MPI_Irecv requests that should be cancelled here
2007!! TODO:
2008! if (.not.lmp_sync) then
2009!   r=mp_pid*nvar_async
2010!   do j=r,r+nvar_async-1
2011!     call MPI_Cancel(j,mp_ierr)
2012!     if (mp_ierr /= 0) write(*,*) '#### mpif_finalize::MPI_Cancel> ERROR ####'
2013!   end do
2014! end if
2015
2016    call MPI_FINALIZE(mp_ierr)
2017    if (mp_ierr /= 0) then
2018      write(*,*) '#### mpif_finalize::MPI_FINALIZE> MPI ERROR ', mp_ierr, ' ####'
2019      stop
2020    end if
2021
2022
2023    end subroutine mpif_finalize
2024
2025
2026    subroutine get_lun(my_lun)
2027!***********************************************************************
2028! get_lun:
2029!   Starting from 100, get next free logical unit number
2030!***********************************************************************
2031
2032      implicit none
2033
2034      integer, intent(inout) :: my_lun
2035      integer, save :: free_lun=100
2036      logical :: exists, iopen
2037
2038!***********************************************************************
2039
2040      loop1: do
2041        inquire(UNIT=free_lun, EXIST=exists, OPENED=iopen)
2042        if (exists .and. .not.iopen) exit loop1
2043        free_lun = free_lun+1
2044      end do loop1
2045      my_lun = free_lun
2046
2047    end subroutine get_lun
2048
2049
2050    subroutine write_data_dbg(array_in, array_name, tstep, ident)
2051!***********************************************************************
2052! Write one-dimensional arrays to file (for debugging purposes)
2053!***********************************************************************
2054      implicit none
2055
2056      real, intent(in), dimension(:) :: array_in
2057      integer, intent(in) :: tstep
2058      integer :: lios
2059      character(LEN=*), intent(in) :: ident, array_name
2060
2061      character(LEN=8) :: c_ts
2062      character(LEN=40) :: fn_1, fn_2
2063
2064!***********************************************************************
2065
2066      write(c_ts, FMT='(I8.8,BZ)') tstep
2067      fn_1='-'//trim(adjustl(c_ts))//'-'//trim(ident)
2068      write(c_ts, FMT='(I2.2,BZ)') mp_np
2069      fn_2= trim(adjustl(array_name))//trim(adjustl(fn_1))//'-np'//trim(adjustl(c_ts))//'.dat'
2070
2071      call get_lun(dat_lun)
2072      open(UNIT=dat_lun, FILE=fn_2, IOSTAT=lios, ACTION='WRITE', &
2073           FORM='UNFORMATTED', STATUS='REPLACE')
2074      write(UNIT=dat_lun, IOSTAT=lios) array_in
2075      close(UNIT=dat_lun)
2076
2077    end subroutine write_data_dbg
2078
2079
2080  subroutine set_fields_synthetic()
2081!*******************************************************************************
2082! DESCRIPTION
2083!   Set all meteorological fields to synthetic (usually constant/homogeneous)
2084!   values.
2085!   Used for validation and error-checking
2086!
2087! NOTE
2088!   This version uses asynchronious communications.
2089!
2090! VARIABLES
2091!   
2092!
2093! TODO
2094!
2095!*******************************************************************************
2096    use com_mod
2097
2098    implicit none
2099
2100    integer :: li=1, ui=2 ! wfmem indices (i.e, operate on entire field)
2101
2102    if (.not.lmp_sync) ui=3
2103
2104
2105! Variables transferred at initialization only
2106!*********************************************
2107!    readclouds=readclouds_
2108    oro(:,:)=0.0
2109    excessoro(:,:)=0.0
2110    lsm(:,:)=0.0
2111    xlanduse(:,:,:)=0.0
2112!    wftime
2113!    numbwf
2114!    nmixz
2115!    height
2116
2117! Time-varying fields:
2118    uu(:,:,:,li:ui) = 10.0
2119    vv(:,:,:,li:ui) = 0.0
2120    uupol(:,:,:,li:ui) = 10.0 ! TODO check if ok
2121    vvpol(:,:,:,li:ui)=0.0
2122    ww(:,:,:,li:ui)=0.
2123    tt(:,:,:,li:ui)=300.
2124    rho(:,:,:,li:ui)=1.3
2125    drhodz(:,:,:,li:ui)=.0
2126    tth(:,:,:,li:ui)=0.0
2127    qvh(:,:,:,li:ui)=1.0
2128    qv(:,:,:,li:ui)=1.0
2129
2130    pv(:,:,:,li:ui)=1.0
2131    clouds(:,:,:,li:ui)=0.0
2132
2133    clwc(:,:,:,li:ui)=0.0
2134    ciwc(:,:,:,li:ui)=0.0
2135 
2136! 2D fields
2137
2138    cloudsh(:,:,li:ui)=0.0
2139    vdep(:,:,:,li:ui)=0.0
2140    ps(:,:,:,li:ui)=1.0e5
2141    sd(:,:,:,li:ui)=0.0
2142    tcc(:,:,:,li:ui)=0.0
2143    tt2(:,:,:,li:ui)=300.
2144    td2(:,:,:,li:ui)=250.
2145    lsprec(:,:,:,li:ui)=0.0
2146    convprec(:,:,:,li:ui)=0.0
2147    ustar(:,:,:,li:ui)=1.0
2148    wstar(:,:,:,li:ui)=1.0
2149    hmix(:,:,:,li:ui)=10000.
2150    tropopause(:,:,:,li:ui)=10000.
2151    oli(:,:,:,li:ui)=0.01
2152    z0=1.0
2153 
2154  end subroutine set_fields_synthetic
2155
2156end module mpi_mod
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG