source: flexpart.git/src/mpi_mod.f90 @ 341f4b7

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

SPECIES files converted to namelist format. Fixed-format SPCIES files renamed with .oldformat extension. Added 2 wet depo parameters to old-format SPECIES. Renamed internal variables and parameters used for wet deposition. incloud_ratio increased to 2.2

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