source: flexpart.git/src/mpi_mod.f90 @ 6b22af9

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

Local branch for espen, working with OpenMP version of readwind

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