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

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

Completed handling of nested wind fields with cloud water (for wet deposition).

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