mpi_mod.f90 110 KB
Newer Older
1
2
! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
! SPDX-License-Identifier: GPL-3.0-or-later
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27

module mpi_mod

!*****************************************************************************
!                                                                            *
!  DESCRIPTION                                                               *
!    This module contains subroutines and common variables used for the      *
!    MPI parallelization of FLEXPART.                                        *
!                                                                            *
!  NOTE                                                                      *
!    Depending on the MPI library installed on your system (e.g. mpich2,     *
!    OpenMPI) you may need to choose below in this file between              *
!      use mpi                                                               *
!    (if the MPI library comes with the file 'mpi.mod'); or                  *
!      include 'mpif.h'                                                      *
!                                                                            *
!                                                                            *
!*****************************************************************************
!                                                                            *
! Variables:                                                                 *
!                                                                            *
! mp_ierr                 MPI error code                                     *
! mp_np                   Number of MPI processes                            *
! mp_pid                  Process ID of each MPI process                     *
! mp_seed                 Parameter for random number seed                   *
Espen Sollum's avatar
Espen Sollum committed
28
29
! read_grp_min            Minimum number of processes at which one will be   *
!                         used as reader                                     *
30
31
32
33
34
35
! numpart_mpi,            Number of particles per node                       *
! maxpart_mpi                                                                *
! mp_partid               MPI process ID for particle calculation            *
! mp_partgroup_           Refers to the subset of processors performing      *
!                         loops over particles. Will be all processes        *
!                         unless a dedicated process runs getfields/readwind *
Espen Sollum's avatar
Espen Sollum committed
36
! lmp_sync                If .false., use asynchronous MPI                   *
37
! mp_cp                   Real precision to use for deposition fields        *
38
39
40
41
!                                                                            *
!                                                                            *
!                                                                            *
!                                                                            *
42
43
44
45
46
47
48
!*****************************************************************************
!                                                                            *
! Modification by DJM, 2017-05-09 - added #ifdef USE_MPIINPLACE cpp          *
! directive to mpif_tm_reduce_grid() to insure that MPI_IN_PLACE is          *
! used in the MPI_Reduce() only if specifically compiled with that           *
! directive.                                                                 *
!                                                                            *
49
50
51
!*****************************************************************************

  use mpi 
52
  use par_mod, only: dp,sp
53
54
55
56
57
58
59
60
61
62
  use com_mod, only: lroot

  implicit none

!  include 'mpif.h'

  public

! Set aside a process for reading windfields if using at least these many processes
!==================================================
Espen Sollum's avatar
Espen Sollum committed
63
  integer, parameter, private :: read_grp_min=4
64
65
66
67
68
69
70
71
!==================================================

! Variables for each MPI process in the world group
  integer :: mp_ierr, mp_np, mp_pid, mp_partid
  integer, private :: world_group_id

! Variables for MPI processes in the 'particle' group
  integer, allocatable, dimension(:) :: mp_partgroup_rank
72
  integer, allocatable, dimension(:) :: npart_per_process
73
74
75
76
  integer :: mp_partgroup_comm, mp_partgroup_pid, mp_partgroup_np

  integer :: mp_seed=0
  integer, parameter :: mp_sp=MPI_REAL4, mp_dp=MPI_REAL8
77
  integer :: mp_cp
78
79
80
  integer, parameter :: id_root=0 ! master process

! MPI tags/requests for send/receive operation
Espen Sollum's avatar
Espen Sollum committed
81
  integer :: tm1
82
  integer, parameter :: nvar_async=26
83
!integer, dimension(:), allocatable :: tags
84
85
  integer, dimension(:), allocatable :: reqs

86
87
88
! Status array used for certain MPI operations (MPI_RECV)
  integer, dimension(MPI_STATUS_SIZE) :: mp_status

89
90
91
92
93
94
95
96
97

  integer :: id_read   ! readwind/getfield process
  integer :: numpart_mpi,maxpart_mpi ! number of particles per node
  integer :: tot_numpart=0
  integer :: mp_comm_used ! global or subgroup communicator

  logical :: lmpreader=.false. ! is set to true for reading process(es) only.
  logical :: lmp_use_reader=.false. ! true if separate readwind process is used

98
! .true. if only using synchronous MPI send/recv (default)
99
100
! If setting this to .false., numwfmem must be set to 3
!===============================================================================
101
  logical :: lmp_sync=.true. 
102
103
!===============================================================================

Espen Sollum's avatar
Espen Sollum committed
104
! mp_dbg_mode       Used for debugging MPI.
105
106
107
108
! mp_dev_mode       various checks related to debugging the parallel code
! mp_dbg_out        write some arrays to data file for debugging
! mp_measure_time   Measure cpu/wall time, write out at end of run
! mp_time_barrier   Measure MPI barrier time
109
! mp_exact_numpart  Use an extra MPI communication to give the exact number of particles
110
!                   to standard output (this does not otherwise affect the simulation)
111
112
113
  logical, parameter :: mp_dbg_mode = .false.
  logical, parameter :: mp_dev_mode = .false.
  logical, parameter :: mp_dbg_out = .false.
114
  logical, parameter :: mp_time_barrier=.true.
Espen Sollum's avatar
Espen Sollum committed
115
  logical, parameter :: mp_measure_time=.false.
116
  logical, parameter :: mp_exact_numpart=.true.
117
118

! for measuring CPU/Wall time
119
120
121
122
123
124
125
126
127
128
  real(sp),private :: mp_comm_time_beg, mp_comm_time_end, mp_comm_time_total=0.
  real(dp),private :: mp_comm_wtime_beg, mp_comm_wtime_end, mp_comm_wtime_total=0.
  real(sp),private :: mp_root_time_beg, mp_root_time_end, mp_root_time_total=0.
  real(dp),private :: mp_root_wtime_beg, mp_root_wtime_end, mp_root_wtime_total=0.
  real(sp),private :: mp_barrier_time_beg, mp_barrier_time_end, mp_barrier_time_total=0.
  real(dp),private :: mp_barrier_wtime_beg, mp_barrier_wtime_end, mp_barrier_wtime_total=0.
  real(sp),private :: tm_nploop_beg, tm_nploop_end, tm_nploop_total=0.
  real(sp),private :: tm_tot_beg, tm_tot_end, tm_tot_total=0.
  real(dp),private :: mp_getfields_wtime_beg, mp_getfields_wtime_end, mp_getfields_wtime_total=0.
  real(sp),private :: mp_getfields_time_beg, mp_getfields_time_end, mp_getfields_time_total=0.
129
130
  real(dp),private :: mp_readwind_wtime_beg, mp_readwind_wtime_end, mp_readwind_wtime_total=0.
  real(sp),private :: mp_readwind_time_beg, mp_readwind_time_end, mp_readwind_time_total=0.
131
132
133
134
135
136
137
138
139
  real(dp),private :: mp_io_wtime_beg, mp_io_wtime_end, mp_io_wtime_total=0.
  real(sp),private :: mp_io_time_beg, mp_io_time_end, mp_io_time_total=0.
  real(dp),private :: mp_wetdepo_wtime_beg, mp_wetdepo_wtime_end, mp_wetdepo_wtime_total=0.
  real(sp),private :: mp_wetdepo_time_beg, mp_wetdepo_time_end, mp_wetdepo_time_total=0.
  real(dp),private :: mp_advance_wtime_beg, mp_advance_wtime_end, mp_advance_wtime_total=0.
  real(dp),private :: mp_conccalc_time_beg, mp_conccalc_time_end, mp_conccalc_time_total=0.
  real(dp),private :: mp_total_wtime_beg, mp_total_wtime_end, mp_total_wtime_total=0.
  real(dp),private :: mp_vt_wtime_beg, mp_vt_wtime_end, mp_vt_wtime_total
  real(sp),private :: mp_vt_time_beg, mp_vt_time_end, mp_vt_time_total
140
141
142
143

! dat_lun           logical unit number for i/o
  integer, private :: dat_lun 

144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
! Temporary arrays for particles (allocated and deallocated as needed)
  integer, allocatable, dimension(:) :: nclass_tmp, npoint_tmp, itra1_tmp, idt_tmp, &
       & itramem_tmp, itrasplit_tmp
  real(kind=dp), allocatable, dimension(:) :: xtra1_tmp, ytra1_tmp
  real, allocatable, dimension(:) :: ztra1_tmp 
  real, allocatable, dimension(:,:) :: xmass1_tmp

! mp_redist_fract        Exchange particles between processes if relative difference
!                        is larger. A good value is between 0.0 and 0.5
! mp_maxpart_factor      Allocate more memory per process than strictly needed
!                        (safety factor). Recommended value between 1.5 and 2.5
! mp_min_redist          Do not redistribute particles if below this limit
  real, parameter :: mp_redist_fract=0.2, mp_maxpart_factor=1.5
  integer,parameter :: mp_min_redist=100000


160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
contains

  subroutine mpif_init
!***********************************************************************
! DESCRIPTION
!   Initialize MPI.
!   
!   Create the global communicator MPI_COMM_WORLD
!   If using a separate MPI process for getfields/readwind, a subgroup
!   is created for the other processes.
!
! VARIABLES
!   mpi_mode    default 0, set to 2/3 if running MPI version
!   mp_np       number of running processes, decided at run-time
!***********************************************************************
175
176
    use par_mod, only: maxpart, numwfmem, dep_prec, maxreceptor, maxspec
    use com_mod, only: mpi_mode, verbosity, creceptor0
177
178
179
180
181

    implicit none

    integer :: i,j,s,addmaxpart=0

182
! Each process gets an ID (mp_pid) in the range 0,..,mp_np-1
183
184
185
186
187
188
189
190
    call MPI_INIT(mp_ierr)
    if (mp_ierr /= 0) goto 100
    call MPI_COMM_RANK(MPI_COMM_WORLD, mp_pid, mp_ierr)
    if (mp_ierr /= 0) goto 100
    call MPI_COMM_SIZE(MPI_COMM_WORLD, mp_np, mp_ierr)
    if (mp_ierr /= 0) goto 100


191
! Variable mpi_mode is used to handle subroutines common to parallel/serial version
192
193
194
195
196
197
198
199
200
201
    if (lmp_sync) then
      mpi_mode=2 ! hold 2 windfields in memory
    else
      mpi_mode=3 ! hold 3 windfields in memory
    end if

    if (mp_pid.ne.0) then
      lroot = .false.
    end if

202
203
204
205
! Set MPI precision to use for transferring deposition fields
!************************************************************
    if (dep_prec==dp) then
      mp_cp = MPI_REAL8
206
! TODO: write info message for serial version as well 
Espen Sollum's avatar
Espen Sollum committed
207
      if (lroot.and.verbosity>0) write(*,*) 'Using double precision for deposition fields'
208
209
    else if (dep_prec==sp) then
      mp_cp = MPI_REAL4
Espen Sollum's avatar
Espen Sollum committed
210
      if (lroot.and.verbosity>0) write(*,*) 'Using single precision for deposition fields'
211
212
213
214
215
    else
      write(*,*) 'ERROR: something went wrong setting MPI real precision'
      stop
    end if

216
217
218
! Check for sensible combination of parameters
!*********************************************

219
220
221
222
223
224
225
226
    if (.not.lmp_sync.and.numwfmem.ne.3) then
      if (lroot) then
        write(*,FMT='(80("#"))')
        write(*,*) '#### mpi_mod::mpif_init> ERROR: ', &
             & 'numwfmem must be set to 3 for asyncronous reading ####'
        write(*,FMT='(80("#"))')
      end if
      call MPI_FINALIZE(mp_ierr)
227
228
229
230
      stop
    else if (lmp_sync.and.numwfmem.ne.2.and.lroot) then
      write(*,FMT='(80("#"))')
      write(*,*) '#### mpi_mod::mpif_init> WARNING: ', &
231
           & 'numwfmem should be set to 2 for syncronous'
232
233
      write(*,*) ' reading. Results will still be valid, but unneccesary memory &
           &is allocated.'
234
      write(*,FMT='(80("#"))')
Espen Sollum's avatar
Espen Sollum committed
235
! Force "syncronized" version if all processes will call getfields
236
    else if ((.not.lmp_sync.and.mp_np.lt.read_grp_min).or.(mp_np.eq.1)) then
Espen Sollum's avatar
Espen Sollum committed
237
238
239
240
241
242
      if (lroot) then
        write(*,FMT='(80("#"))')
        write(*,*) '#### mpi_mod::mpif_init> WARNING: ', &
             & 'all procs call getfields. Setting lmp_sync=.true.'
        write(*,FMT='(80("#"))')
      end if
243
      lmp_sync=.true.
244
    end if
245
246

! TODO: Add more warnings for unimplemented flexpart features
247
248
249
250
251
252

! Set ID of process that calls getfield/readwind. 
! Using the last process in the group ensures statistical identical results
! as running with one process less but not using separate read process
!**********************************************************************

253
!    id_read = min(mp_np-1, 1)
254
255
256
257
258
259
    id_read = mp_np-1

    if (mp_pid.eq.id_read) lmpreader=.true.

    call MPI_Comm_group (MPI_COMM_WORLD, world_group_id, mp_ierr)

Espen Sollum's avatar
Espen Sollum committed
260
! Create a MPI group/communicator that will calculate trajectories. 
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
! Skip this step if program is run with only a few processes
!************************************************************************
    allocate(mp_partgroup_rank(0:mp_np-2))

! This allows checking for allocation of arrays when no subgroup is used
    mp_partgroup_pid=0

    if (read_grp_min.lt.2) then
      write(*,*) '#### mpi_mod::mpif_init> ERROR ####', &
           & 'read_grp_min should be at least 2. Exiting'
      stop
    end if

    if (mp_np.ge.read_grp_min) then
      lmp_use_reader = .true.

Espen Sollum's avatar
Espen Sollum committed
277
! Map the subgroup IDs= 0,1,2,...,mp_np-2, skipping reader process
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
      j=-1 
      do i=0, mp_np-2 !loop over all (subgroup) IDs
        j=j+1
        if (i.eq.id_read) then
          j=j+1
        end if
        mp_partgroup_rank(i) = j
      end do

      call MPI_Group_incl (world_group_id, mp_np-1, mp_partgroup_rank, &
           &mp_partgroup_pid, mp_ierr)
      if (mp_ierr /= 0) goto 100
      call MPI_Comm_create (MPI_COMM_WORLD, mp_partgroup_pid, mp_partgroup_comm, mp_ierr)
      if (mp_ierr /= 0) goto 100

      if (mp_pid.ne.id_read) then
        call MPI_Comm_rank (mp_partgroup_comm, mp_partgroup_pid, mp_ierr)
        if (mp_ierr /= 0) goto 100

! Get size of new group (=mp_np-1)
        call MPI_COMM_SIZE(mp_partgroup_comm, mp_partgroup_np, mp_ierr)
        if (mp_ierr /= 0) goto 100
        if (mp_partgroup_np.ne.mp_np-1) then
          write(*,*)  '#### mpi_mod:: mpif_init> ERROR ####&
               & mp_partgroup_np.ne.mp_np-1'
          stop
        endif

      else
        mp_partgroup_pid = -1
      end if
    end if

    if (lmp_use_reader) then
      mp_comm_used = mp_partgroup_comm
      mp_partid = mp_partgroup_pid
      mp_partgroup_np=mp_np-1
    else
      mp_comm_used = MPI_COMM_WORLD
      mp_partgroup_np = mp_np
      mp_partid = mp_pid
    end if

! Set maxpart per process
322
! ESO 08/2016: Increase maxpart per process, in case of unbalanced distribution
323
    maxpart_mpi=int(mp_maxpart_factor*real(maxpart)/real(mp_partgroup_np))
324
    if (mp_np == 1) maxpart_mpi = maxpart
325
326
327
328
329
330
331
332
333
334
335

! Do not allocate particle data arrays for readwind process
    if (lmpreader.and.lmp_use_reader) then
      maxpart_mpi=0
    end if

! Set random seed for each non-root process
    if (mp_pid.gt.0) then
      s = 244
      mp_seed = -abs(mod((s*181)*((mp_pid-83)*359), 104729))
    end if
Espen Sollum's avatar
Espen Sollum committed
336
337
338
339
    if (mp_dbg_mode) then
      mp_seed=0
      if (lroot) write(*,*) 'MPI: setting seed=0'
    end if
340
341
342
343

! Allocate request array for asynchronous MPI
    if (.not.lmp_sync) then
      allocate(reqs(0:nvar_async*mp_np-1))
Espen Sollum's avatar
Espen Sollum committed
344
      reqs(:)=MPI_REQUEST_NULL
345
346
    else
      allocate(reqs(0:1))
Espen Sollum's avatar
Espen Sollum committed
347
      reqs(:)=MPI_REQUEST_NULL
348
349
    end if

350
351
352
! Allocate array for number of particles per process    
    allocate(npart_per_process(0:mp_partgroup_np-1))

353
354
355
356
357
358
359
! Write whether MPI_IN_PLACE is used or not
#ifdef USE_MPIINPLACE
    if (lroot) write(*,*) 'Using MPI_IN_PLACE operations'
#else
    if (lroot) allocate(creceptor0(maxreceptor,maxspec))
    if (lroot) write(*,*) 'Not using MPI_IN_PLACE operations'
#endif
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
    goto 101

100 write(*,*) '#### mpi_mod::mpif_init> ERROR ####', mp_ierr
    stop

101 end subroutine mpif_init


  subroutine mpif_mpi_barrier
!***********************************************************************
! Set MPI barrier (all processes are synchronized here), with the
! possible exception of a separate 'readwind' process.
! Optionally measure cpu/wall time.
!
!***********************************************************************
    implicit none

    if (mp_time_barrier) then
      call cpu_time(mp_barrier_time_beg)
      mp_barrier_wtime_beg = mpi_wtime()
    endif

    call MPI_BARRIER(mp_comm_used, mp_ierr)
    if (mp_ierr /= 0) goto 100

    if (mp_time_barrier) then
      call cpu_time(mp_barrier_time_end)
      mp_barrier_wtime_end = mpi_wtime()

      mp_barrier_time_total=mp_barrier_time_total+&
           &(mp_barrier_time_end-mp_barrier_time_beg)
      mp_barrier_wtime_total=mp_barrier_wtime_total+&
           &(mp_barrier_wtime_end-mp_barrier_wtime_beg)
    end if

    goto 101

100 write(*,*) '#### mpi_mod::mpif_mpi_barrier> ERROR ####', mp_ierr
    stop

101 end subroutine mpif_mpi_barrier


  subroutine mpif_com_mod_allocate
!*******************************************************************************    
! Dynamic allocation of arrays (in serial code these have size maxpart)
!
407
408
409
! Not in use anymore, moved to com_mod for interoperability with serial version
! 
! TODO: remove 
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
!*******************************************************************************
    use com_mod

    implicit none 

    integer :: array_size

    array_size = maxpart_mpi

    allocate(itra1(array_size),npoint(array_size),&
         & nclass(array_size),idt(array_size),&
         & itramem(array_size),itrasplit(array_size),&
         & xtra1(array_size),ytra1(array_size),&
         & ztra1(array_size),xmass1(array_size, maxspec))

    allocate(uap(array_size),ucp(array_size),&
         & uzp(array_size),us(array_size),&
         & vs(array_size),&
         & ws(array_size),cbt(array_size))


  end subroutine mpif_com_mod_allocate


  subroutine mpif_tm_send_dims
!***********************************************************************
! Distribute array dimensions from pid0 to all processes.
! numpart_mpi/numpart is sent to allow dynamic allocation
!
! Currently not in use.
!
! Variables of similar type (integer, real) are placed in an array
! and sent collectively, to avoid the overhead associated with individual
! MPI send/receives
!
!
!***********************************************************************
    use com_mod

    implicit none

    integer :: add_part=0

Espen Sollum's avatar
Espen Sollum committed
453
    call MPI_Bcast(numpart, 1, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471

! MPI subgroup does the particle-loop 
    if (lmp_use_reader) then
      if (mod(numpart,mp_partgroup_np).ne.0) add_part=1
      numpart_mpi=int(numpart/mp_partgroup_np)+add_part
    else
      if (mod(numpart,mp_np).ne.0) add_part=1
      numpart_mpi=int(numpart/mp_np)+add_part
    end if


! redefine numpart as 'numpart per process' throughout the code
!**************************************************************
    numpart = numpart_mpi

  end subroutine mpif_tm_send_dims


472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
  subroutine mpif_send_part_properties(num_part)
!***********************************************************************
! Distribute particle properties from root to all processes.
!  
! Used by init_domainfill_mpi
!
! Variables:
! 
! num_part        input, number of particles per process (rounded up)
!
!***********************************************************************
    use com_mod

    implicit none

    integer,intent(in) :: num_part
    integer :: i,jj, addone

490
491
492
493
494
495
496
497
498
499
500
501
502
! Exit if too many particles
    if (num_part.gt.maxpart_mpi) then
      write(*,*) '#####################################################'
      write(*,*) '#### ERROR - TOTAL NUMBER OF PARTICLES REQUIRED  ####'
      write(*,*) '#### EXCEEDS THE MAXIMUM ALLOWED NUMBER. REDUCE  ####'
      write(*,*) '#### EITHER NUMBER OF PARTICLES PER RELEASE POINT####'
      write(*,*) '#### OR INCREASE MAXPART.                        ####'
      write(*,*) '#####################################################'
!      call MPI_FINALIZE(mp_ierr)
      stop
    end if


503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
! Time for MPI communications
!****************************
    if (mp_measure_time) call mpif_mtime('commtime',0)

! Distribute variables, send from pid 0 to other processes (including itself)
!****************************************************************************

    call MPI_SCATTER(nclass_tmp,num_part,MPI_INTEGER,nclass,&
         &num_part,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(npoint_tmp,num_part,MPI_INTEGER,npoint,&
         &num_part,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(itra1_tmp,num_part,MPI_INTEGER,itra1,&
         &num_part,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(idt_tmp,num_part,MPI_INTEGER,idt,&
         &num_part,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(itramem_tmp,num_part,MPI_INTEGER,itramem,&
         &num_part,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(itrasplit_tmp,num_part,MPI_INTEGER,itrasplit,&
         &num_part,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(xtra1_tmp,num_part,mp_dp,xtra1,&
         &num_part,mp_dp,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(ytra1_tmp,num_part,mp_dp,ytra1,&
         &num_part,mp_dp,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    call MPI_SCATTER(ztra1_tmp,num_part,mp_sp,ztra1,&
         &num_part,mp_sp,id_root,mp_comm_used,mp_ierr)
    if (mp_ierr /= 0) goto 600 
    do i=1,nspec
      call MPI_SCATTER(xmass1_tmp(:,i),num_part,mp_sp,xmass1(:,i),&
           &num_part,mp_sp,id_root,mp_comm_used,mp_ierr)
      if (mp_ierr /= 0) goto 600 
    end do

    if (mp_measure_time) call mpif_mtime('commtime',1)

    goto 601

600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
    stop

! After the transfer of particles it it possible that the value of
551
! "numpart" is larger than the actual used, so we reduce it if there are
552
! invalid particles at the end of the arrays
553

554
601 do i=numpart, 1, -1
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
      if (itra1(i).eq.-999999999) then
        numpart=numpart-1
      else
        exit
      end if
    end do


  end subroutine mpif_send_part_properties


  subroutine mpif_calculate_part_redist(itime)
!***********************************************************************
! Calculate number of particles to redistribute between processes. This routine 
! can be called at regular time intervals to keep a level number of
! particles on each process.
!
! First, all processes report their local "numpart" to each other, which is
! stored in array "numpart_mpi(np)". This array is sorted from low to
! high values, along with a corresponding process ID array "idx_arr(np)".
! If the relative difference between the highest and lowest "numpart_mpi" value
! is above a threshold, particles are transferred from process idx_arr(np-1)
! to process idx_arr(0). Repeat for processes idx_arr(np-i) and idx_arr(i)
! for all valid i. 
! Note: If np is an odd number, the process with the median
! number of particles is left unchanged
!
! VARIABLES
! itime        input, current time
!***********************************************************************
    use com_mod

    implicit none
    
    integer, intent(in) :: itime
    real :: pmin,z
    integer :: i,jj,nn, num_part=1,m,imin, num_trans
    logical :: first_iter
593
    integer,allocatable,dimension(:) :: idx_arr
594
595
596
597
598
599
600
601
    real,allocatable,dimension(:) :: sorted ! TODO: we don't really need this

! Exit if running with only 1 process
!************************************************************************
    if (mp_np.eq.1) return

! All processes exchange information on number of particles
!****************************************************************************
602
    allocate( idx_arr(0:mp_partgroup_np-1), sorted(0:mp_partgroup_np-1))
603

604
    call MPI_Allgather(numpart, 1, MPI_INTEGER, npart_per_process, &
605
606
607
608
609
         & 1, MPI_INTEGER, mp_comm_used, mp_ierr)


! Sort from lowest to highest
! Initial guess: correct order
610
    sorted(:) = npart_per_process(:)
611
612
613
614
    do i=0, mp_partgroup_np-1
      idx_arr(i) = i
    end do

615
616
617
! Do not rebalance particles for ipout=3    
    if (ipout.eq.3) return

618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
! For each successive element in index array, see if a lower value exists
    do i=0, mp_partgroup_np-2
      pmin=sorted(i)
      imin=idx_arr(i)
      m=i+1
      do jj=m, mp_partgroup_np-1
        if (pmin.le.sorted(jj)) cycle
        z=pmin
        pmin=sorted(jj)
        sorted(jj)=z

        nn=imin
        imin=idx_arr(jj)
        idx_arr(jj)=nn
      end do
      sorted(i)=pmin
      idx_arr(i)=imin
    end do

! For each pair of processes, transfer particles if the difference is above a
! limit, and numpart at sending process large enough

    m=mp_partgroup_np-1 ! index for last sorted process (most particles)
    do i=0,mp_partgroup_np/2-1
642
      num_trans = npart_per_process(idx_arr(m)) - npart_per_process(idx_arr(i))
643
      if (mp_partid.eq.idx_arr(m).or.mp_partid.eq.idx_arr(i)) then
644
645
        if ( npart_per_process(idx_arr(m)).gt.mp_min_redist.and.&
             & real(num_trans)/real(npart_per_process(idx_arr(m))).gt.mp_redist_fract) then
646
! DBG
647
648
          ! write(*,*) 'mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, npart_per_process', &
          !      &mp_partid, idx_arr(m), idx_arr(i), mp_min_redist, num_trans, npart_per_process
649
! DBG
650
651
652
653
654
655
          call mpif_redist_part(itime, idx_arr(m), idx_arr(i), num_trans/2)
        end if
      end if
      m=m-1
    end do

656
    deallocate(idx_arr, sorted)
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682

  end subroutine mpif_calculate_part_redist


  subroutine mpif_redist_part(itime, src_proc, dest_proc, num_trans)
!***********************************************************************
! Transfer particle properties between two arbitrary processes.
! 
! VARIABLES
!
! itime           input, current time
! src_proc        input, ID of source process
! dest_proc       input, ID of destination process
! num_trans       input, number of particles to transfer
!
!************************************************************************
    use com_mod !TODO:  ,only: nclass etc

    implicit none 

    integer, intent(in) :: itime, src_proc, dest_proc, num_trans
    integer :: ll, ul ! lower and upper indices in arrays
    integer :: arr_sz ! size of temporary arrays
    integer :: mtag   ! MPI message tag
    integer :: i, j, minpart, ipart, maxnumpart
 
683
684
685
686
687
688
689
690
! Check for invalid input arguments
!**********************************
 if (src_proc.eq.dest_proc) then
   write(*,*) '<mpi_mod::mpif_redist_part>: Error: &
        &src_proc.eq.dest_proc' 
   stop
 end if

691
692
693
694
695
696
697
698
699
700
701
702
703
! Measure time for MPI communications
!************************************
    if (mp_measure_time) call mpif_mtime('commtime',0)

! Send particles
!***************
    if (mp_partid.eq.src_proc) then
      mtag = 2000

! Array indices for the transferred particles
      ll=numpart-num_trans+1
      ul=numpart

704
705
706
707
708
      if (mp_dev_mode) then
        write(*,FMT='(72("#"))')
        write(*,*) "Sending ", num_trans, "particles (from/to)", src_proc, dest_proc
        write(*,*) "numpart before: ", numpart
      end if
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746

      call MPI_SEND(nclass(ll:ul),num_trans,&
           & MPI_INTEGER,dest_proc,mtag+1,mp_comm_used,mp_ierr)

      call MPI_SEND(npoint(ll:ul),num_trans,&
           & MPI_INTEGER,dest_proc,mtag+2,mp_comm_used,mp_ierr)

      call MPI_SEND(itra1(ll:ul),num_trans, &
           & MPI_INTEGER,dest_proc,mtag+3,mp_comm_used,mp_ierr)

      call MPI_SEND(idt(ll:ul),num_trans, &
           & MPI_INTEGER,dest_proc,mtag+4,mp_comm_used,mp_ierr)

      call MPI_SEND(itramem(ll:ul),num_trans, &
           & MPI_INTEGER,dest_proc,mtag+5,mp_comm_used,mp_ierr)

      call MPI_SEND(itrasplit(ll:ul),num_trans,&
           & MPI_INTEGER,dest_proc,mtag+6,mp_comm_used,mp_ierr)

      call MPI_SEND(xtra1(ll:ul),num_trans, &
           & mp_dp,dest_proc,mtag+7,mp_comm_used,mp_ierr)

      call MPI_SEND(ytra1(ll:ul),num_trans,&
           & mp_dp,dest_proc,mtag+8,mp_comm_used,mp_ierr)

      call MPI_SEND(ztra1(ll:ul),num_trans,&
           & mp_sp,dest_proc,mtag+9,mp_comm_used,mp_ierr)

      do j=1,nspec
        call MPI_SEND(xmass1(ll:ul,j),num_trans,&
             & mp_sp,dest_proc,mtag+(9+j),mp_comm_used,mp_ierr)
      end do

! Terminate transferred particles, update numpart
      itra1(ll:ul) = -999999999

      numpart = numpart-num_trans

747
748
749
750
      if (mp_dev_mode) then
        write(*,*) "numpart after: ", numpart
        write(*,FMT='(72("#"))')
      end if
751
752
753
754
755
756
757
758
759
760

    else if (mp_partid.eq.dest_proc) then

! Receive particles
!******************
      mtag = 2000
      ! Array indices for the transferred particles
      ll=numpart+1
      ul=numpart+num_trans

761
762
763
764
765
      if (mp_dev_mode) then
        write(*,FMT='(72("#"))')
        write(*,*) "Receiving ", num_trans, "particles (from/to)", src_proc, dest_proc
        write(*,*) "numpart before: ", numpart
      end if
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814

! We could receive the data directly at nclass(ll:ul) etc., but this leaves
! vacant spaces at lower indices. Using temporary arrays instead.
      arr_sz = ul-ll+1
      allocate(itra1_tmp(arr_sz),npoint_tmp(arr_sz),nclass_tmp(arr_sz),&
           & idt_tmp(arr_sz),itramem_tmp(arr_sz),itrasplit_tmp(arr_sz),&
           & xtra1_tmp(arr_sz),ytra1_tmp(arr_sz),ztra1_tmp(arr_sz),&
           & xmass1_tmp(arr_sz, maxspec))
      
      call MPI_RECV(nclass_tmp,num_trans,&
           & MPI_INTEGER,src_proc,mtag+1,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(npoint_tmp,num_trans,&
           & MPI_INTEGER,src_proc,mtag+2,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(itra1_tmp,num_trans, &
           & MPI_INTEGER,src_proc,mtag+3,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(idt_tmp,num_trans, &
           & MPI_INTEGER,src_proc,mtag+4,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(itramem_tmp,num_trans, &
           & MPI_INTEGER,src_proc,mtag+5,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(itrasplit_tmp,num_trans,&
           & MPI_INTEGER,src_proc,mtag+6,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(xtra1_tmp,num_trans, &
           & mp_dp,src_proc,mtag+7,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(ytra1_tmp,num_trans,&
           & mp_dp,src_proc,mtag+8,mp_comm_used,mp_status,mp_ierr)

      call MPI_RECV(ztra1_tmp,num_trans,&
           & mp_sp,src_proc,mtag+9,mp_comm_used,mp_status,mp_ierr)

      do j=1,nspec
        call MPI_RECV(xmass1_tmp(:,j),num_trans,&
             & mp_sp,src_proc,mtag+(9+j),mp_comm_used,mp_status,mp_ierr)
      end do

! This is the maximum value numpart can possibly have after the transfer
      maxnumpart=numpart+num_trans

! Search for vacant space and copy from temporary storage
!********************************************************
      minpart=1
      do i=1, num_trans
! Take into acount that we may have transferred invalid particles
815
        if (itra1_tmp(i).ne.itime) cycle
816
817
        do ipart=minpart,maxnumpart
          if (itra1(ipart).ne.itime) then
818
819
820
821
822
823
824
825
826
827
            itra1(ipart) = itra1_tmp(i)
            npoint(ipart) = npoint_tmp(i)
            nclass(ipart) = nclass_tmp(i)
            idt(ipart) = idt_tmp(i)
            itramem(ipart) = itramem_tmp(i)
            itrasplit(ipart) = itrasplit_tmp(i)
            xtra1(ipart) = xtra1_tmp(i)
            ytra1(ipart) = ytra1_tmp(i)
            ztra1(ipart) = ztra1_tmp(i)
            xmass1(ipart,:) = xmass1_tmp(i,:)
828
829
            goto 200 ! Storage space has been found, stop searching
          end if
830
! :TODO: add check for if too many particles requiried
831
        end do
832
200     minpart=ipart+1
833
      end do
834
835
! Increase numpart, if necessary
      numpart=max(numpart,ipart)
836

837
838
      deallocate(itra1_tmp,npoint_tmp,nclass_tmp,idt_tmp,itramem_tmp,&
           &itrasplit_tmp,xtra1_tmp,ytra1_tmp,ztra1_tmp,xmass1_tmp)
839

840
841
842
843
      if (mp_dev_mode) then
        write(*,*) "numpart after: ", numpart
        write(*,FMT='(72("#"))')
      end if
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858

    else
! This routine should only be called by the two participating processes
      write(*,*) "ERROR: wrong process has entered mpi_mod::mpif_redist_part"
      stop
      return
    end if

! Measure time for MPI communications
!************************************
    if (mp_measure_time) call mpif_mtime('commtime',1)

  end subroutine mpif_redist_part


859
860
861
862
  subroutine mpif_tm_send_vars
!***********************************************************************
! Distribute particle variables from pid0 to all processes.
! Called from timemanager
863
! *NOT IN USE* at the moment, but can be useful for debugging
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
!
!***********************************************************************
    use com_mod

    implicit none

    integer :: i

! Time for MPI communications
!****************************
    if (mp_measure_time) call mpif_mtime('commtime',0)

! Distribute variables, send from pid 0 to other processes
!**********************************************************************

! integers
    if (lroot) then
      call MPI_SCATTER(npoint,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
882
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
883
884
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(idt,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
885
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
886
887
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(itra1,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
888
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
889
890
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(nclass,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
891
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
892
893
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(itramem,numpart_mpi,MPI_INTEGER,MPI_IN_PLACE,&
894
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
895
896
897
898
      if (mp_ierr /= 0) goto 600 

! int2
      call MPI_SCATTER(cbt,numpart_mpi,MPI_INTEGER2,MPI_IN_PLACE,&
899
           & numpart_mpi,MPI_INTEGER2,id_root,mp_comm_used,mp_ierr)
900
901
902
      if (mp_ierr /= 0) goto 600

! real
903
      call MPI_SCATTER(uap,numpart_mpi,mp_sp,MPI_IN_PLACE,&
904
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
905
      if (mp_ierr /= 0) goto 600
906
      call MPI_SCATTER(ucp,numpart_mpi,mp_sp,MPI_IN_PLACE,&
907
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
908
      if (mp_ierr /= 0) goto 600
909
      call MPI_SCATTER(uzp,numpart_mpi,mp_sp,MPI_IN_PLACE,&
910
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
911
912
      if (mp_ierr /= 0) goto 600

913
      call MPI_SCATTER(us,numpart_mpi,mp_sp,MPI_IN_PLACE,&
914
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
915
      if (mp_ierr /= 0) goto 600
916
      call MPI_SCATTER(vs,numpart_mpi,mp_sp,MPI_IN_PLACE,&
917
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
918
      if (mp_ierr /= 0) goto 600
919
      call MPI_SCATTER(ws,numpart_mpi,mp_sp,MPI_IN_PLACE,&
920
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
921
922
923
      if (mp_ierr /= 0) goto 600

      call MPI_SCATTER(xtra1,numpart_mpi,mp_dp,MPI_IN_PLACE,&
924
           & numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
925
926
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(ytra1,numpart_mpi,mp_dp,MPI_IN_PLACE,&
927
           & numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
928
      if (mp_ierr /= 0) goto 600 
929
      call MPI_SCATTER(ztra1,numpart_mpi,mp_sp,MPI_IN_PLACE,&
930
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
931
932
933
      if (mp_ierr /= 0) goto 600 

      do i=1,nspec
934
        call MPI_SCATTER(xmass1(:,i),numpart_mpi,mp_sp,MPI_IN_PLACE,&
935
             & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
936
937
938
939
940
941
        if (mp_ierr /= 0) goto 600 
      end do

    else ! (mp_pid >= 1)
! integers
      call MPI_SCATTER(npoint,numpart_mpi,MPI_INTEGER,npoint,&
942
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
943
944
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(idt,numpart_mpi,MPI_INTEGER,idt,&
945
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
946
947
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(itra1,numpart_mpi,MPI_INTEGER,itra1,&
948
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
949
950
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(nclass,numpart_mpi,MPI_INTEGER,nclass,&
951
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
952
953
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(itramem,numpart_mpi,MPI_INTEGER,itramem,&
954
           & numpart_mpi,MPI_INTEGER,id_root,mp_comm_used,mp_ierr)
955
956
957
958
      if (mp_ierr /= 0) goto 600 

! int2
      call MPI_SCATTER(cbt,numpart_mpi,MPI_INTEGER2,cbt,&
959
           & numpart_mpi,MPI_INTEGER2,id_root,mp_comm_used,mp_ierr)
960
961
962
      if (mp_ierr /= 0) goto 600

! reals
963
      call MPI_SCATTER(uap,numpart_mpi,mp_sp,uap,&
964
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
965
      if (mp_ierr /= 0) goto 600
966
      call MPI_SCATTER(ucp,numpart_mpi,mp_sp,ucp,&
967
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
968
      if (mp_ierr /= 0) goto 600
969
      call MPI_SCATTER(uzp,numpart_mpi,mp_sp,uzp,&
970
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
971
972
      if (mp_ierr /= 0) goto 600

973
      call MPI_SCATTER(us,numpart_mpi,mp_sp,us,&
974
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
975
      if (mp_ierr /= 0) goto 600
976
      call MPI_SCATTER(vs,numpart_mpi,mp_sp,vs,&
977
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
978
      if (mp_ierr /= 0) goto 600
979
      call MPI_SCATTER(ws,numpart_mpi,mp_sp,ws,&
980
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
981
982
983
      if (mp_ierr /= 0) goto 600

      call MPI_SCATTER(xtra1,numpart_mpi,mp_dp,xtra1,&
984
           & numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
985
986
      if (mp_ierr /= 0) goto 600 
      call MPI_SCATTER(ytra1,numpart_mpi,mp_dp,ytra1,&
987
           & numpart_mpi,mp_dp,id_root,mp_comm_used,mp_ierr)
988
      if (mp_ierr /= 0) goto 600 
989
      call MPI_SCATTER(ztra1,numpart_mpi,mp_sp,ztra1,&
990
           & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
991
992
993
      if (mp_ierr /= 0) goto 600 

      do i=1,nspec
994
        call MPI_SCATTER(xmass1(:,i),numpart_mpi,mp_sp,xmass1(:,i),&
995
             & numpart_mpi,mp_sp,id_root,mp_comm_used,mp_ierr)
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
        if (mp_ierr /= 0) goto 600 
      end do

    end if !(lroot)

    if (mp_measure_time) call mpif_mtime('commtime',1)

    goto 601

600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
    stop
601 end subroutine mpif_tm_send_vars


  subroutine mpif_tm_collect_vars
!*******************************************************************************
! Collect particle data to PID 0 from all processes.                           *
!                                                                              *
! This can be called from timemanager_mpi, before caclulating                  *
! concentrations/writing output grids.                                         *
!                                                                              *
1017
1018
! Currently *not in use*, as each process calculates concentrations            *
! separately, but can be useful for debugging                                  *
1019
1020
1021
1022
1023
!                                                                              *
! To use this routine (together with mpif_tm_send_vars) to transfer data       *
! to the MPI root process (useful for debugging), insert code like this        *
! (in timemanager_mpi.f90)                                                     *
!                                                                              *
1024
!      if (lroot) tot_numpart=numpart ! root stores total numpart              *
1025
1026
1027
1028
1029
1030
1031
1032
1033
!      call mpif_tm_send_dims                                                  *
!      if (numpart>1) then                                                     *
!        call mpif_tm_send_vars                                                *
!      end if                                                                  *
!                                                                              *
!    To collect data from all processes to root process after a parallel       *
!    region:                                                                   *
!                                                                              *
!      if (numpart>0) then                                                     *
1034
!          if (lroot) numpart = tot_numpart                                    *
1035
1036
1037
!       call mpif_tm_collect_vars                                              *
!      end if                                                                  *
!                                                                              *
1038
1039
! Then a section that begins with "if (lroot) ..." will behave like            *
! serial flexpart                                                              *
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
!                                                                              *
!*******************************************************************************
    use com_mod !, only: numpart, nspec, itra1, npoint, nclass

    implicit none

    integer :: i

    logical :: use_mp_vars = .false.! use mp_... type allocated vars
    logical :: use_in_place = .true.! use MPI_IN_PLACE to collect vars.


! Time for MPI communications
    if (mp_measure_time) call mpif_mtime('commtime',0)

! Distribute variables, send from pid 0 to other processes
!**********************************************************************

    if (.not. use_mp_vars.and..not.use_in_place) then

! Integers:
      call MPI_GATHER(npoint, numpart_mpi, MPI_INTEGER, npoint, &
           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 
      call MPI_GATHER(idt, numpart_mpi, MPI_INTEGER, idt, &
           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 
      call MPI_GATHER(itra1, numpart_mpi, MPI_INTEGER, itra1, &
           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 
      call MPI_GATHER(nclass, numpart_mpi, MPI_INTEGER, nclass, &
           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 
      call MPI_GATHER(itramem, numpart_mpi, MPI_INTEGER, itramem, &
           &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 

! int2
      call MPI_GATHER(cbt, numpart_mpi, MPI_INTEGER2, cbt, &
           &numpart_mpi, MPI_INTEGER2, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 

! Reals:
1083
1084
      call MPI_GATHER(uap, numpart_mpi, mp_sp, uap, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1085
      if (mp_ierr /= 0) goto 600 
1086
1087
      call MPI_GATHER(ucp, numpart_mpi, mp_sp, ucp, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1088
      if (mp_ierr /= 0) goto 600 
1089
1090
      call MPI_GATHER(uzp, numpart_mpi, mp_sp, uzp, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1091
1092
      if (mp_ierr /= 0) goto 600 

1093
1094
      call MPI_GATHER(us, numpart_mpi, mp_sp, us, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1095
      if (mp_ierr /= 0) goto 600
1096
1097
      call MPI_GATHER(vs, numpart_mpi, mp_sp, vs, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1098
      if (mp_ierr /= 0) goto 600 
1099
1100
      call MPI_GATHER(ws, numpart_mpi, mp_sp, ws, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1101
1102
1103
1104
1105
1106
1107
1108
1109
      if (mp_ierr /= 0) goto 600


      call MPI_GATHER(xtra1, numpart_mpi, mp_dp, xtra1, &
           &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 
      call MPI_GATHER(ytra1, numpart_mpi, mp_dp, ytra1, &
           &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
      if (mp_ierr /= 0) goto 600 
1110
1111
      call MPI_GATHER(ztra1, numpart_mpi, mp_sp, ztra1, &
           &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1112
1113
1114
      if (mp_ierr /= 0) goto 600 

      do i=1, nspec
1115
1116
        call MPI_GATHER(xmass1(:,i), numpart_mpi, mp_sp, xmass1(:,i), &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
        if (mp_ierr /= 0) goto 600 

      end do

! Use variables npoint etc directly for communications
!***********************************************************************
    else if (use_in_place.and..not.use_mp_vars) then
      if (lroot) then
! Integers:
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, npoint, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, idt, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, itra1, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, nclass, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER, itramem, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 

! int2
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, MPI_INTEGER2, cbt, &
             &numpart_mpi, MPI_INTEGER2, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 

! Reals:
1148
1149
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, uap, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1150
        if (mp_ierr /= 0) goto 600 
1151
1152
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ucp, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1153
        if (mp_ierr /= 0) goto 600 
1154
1155
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, uzp, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1156
1157
        if (mp_ierr /= 0) goto 600 

1158
1159
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, us, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1160
        if (mp_ierr /= 0) goto 600
1161
1162
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, vs, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1163
        if (mp_ierr /= 0) goto 600 
1164
1165
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ws, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1166
1167
1168
1169
1170
1171
1172
1173
1174
        if (mp_ierr /= 0) goto 600


        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_dp, xtra1, &
             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_dp, ytra1, &
             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
1175
1176
        call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, ztra1, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1177
1178
1179
        if (mp_ierr /= 0) goto 600 

        do i=1, nspec
1180
1181
          call MPI_GATHER(MPI_IN_PLACE, numpart_mpi, mp_sp, xmass1(:,i), &
               &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
          if (mp_ierr /= 0) goto 600 
        end do

      else ! (for mp_pid >= 1)

! Integers:
        call MPI_GATHER(npoint, numpart_mpi, MPI_INTEGER, npoint, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(idt, numpart_mpi, MPI_INTEGER, idt, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(itra1, numpart_mpi, MPI_INTEGER, itra1, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(nclass, numpart_mpi, MPI_INTEGER, nclass, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(itramem, numpart_mpi, MPI_INTEGER, itramem, &
             &numpart_mpi, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 

! int2
        call MPI_GATHER(cbt, numpart_mpi, MPI_INTEGER2, cbt, &
             &numpart_mpi, MPI_INTEGER2, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 

! Reals:
1210
1211
        call MPI_GATHER(uap, numpart_mpi, mp_sp, uap, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1212
        if (mp_ierr /= 0) goto 600 
1213
1214
        call MPI_GATHER(ucp, numpart_mpi, mp_sp, ucp, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1215
        if (mp_ierr /= 0) goto 600 
1216
1217
        call MPI_GATHER(uzp, numpart_mpi, mp_sp, uzp, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1218
1219
        if (mp_ierr /= 0) goto 600 

1220
1221
        call MPI_GATHER(us, numpart_mpi, mp_sp, us, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1222
        if (mp_ierr /= 0) goto 600
1223
1224
        call MPI_GATHER(vs, numpart_mpi, mp_sp, vs, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1225
        if (mp_ierr /= 0) goto 600 
1226
1227
        call MPI_GATHER(ws, numpart_mpi, mp_sp, ws, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1228
1229
1230
1231
1232
1233
1234
1235
1236
        if (mp_ierr /= 0) goto 600


        call MPI_GATHER(xtra1, numpart_mpi, mp_dp, xtra1, &
             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
        call MPI_GATHER(ytra1, numpart_mpi, mp_dp, ytra1, &
             &numpart_mpi, mp_dp, id_root, mp_comm_used, mp_ierr)
        if (mp_ierr /= 0) goto 600 
1237
1238
        call MPI_GATHER(ztra1, numpart_mpi, mp_sp, ztra1, &
             &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1239
1240
1241
        if (mp_ierr /= 0) goto 600 

        do i=1, nspec
1242
1243
          call MPI_GATHER(xmass1(:,i), numpart_mpi, mp_sp, xmass1(:,i), &
               &numpart_mpi, mp_sp, id_root, mp_comm_used, mp_ierr)
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
          if (mp_ierr /= 0) goto 600 
        end do
      end if ! (lroot)
    end if !  (use_in_place)

    if (mp_measure_time) call mpif_mtime('commtime',1)

    goto 601

600 write(*,*) "mpi_mod> mp_ierr \= 0", mp_ierr
    stop
601 end subroutine mpif_tm_collect_vars


  subroutine mpif_gf_send_vars(memstat)
!*******************************************************************************
! DESCRIPTION
Espen Sollum's avatar
Espen Sollum committed
1261
!   Distribute 'getfield' variables from reader process
1262
1263
1264
1265
1266
! 
!   Called from timemanager
!
! NOTE
!   This subroutine distributes windfields read from the reader process to
1267
!   all other processes. Usually one set of fields are transfered, but at
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
!   first timestep when there are no fields in memory, both are transfered.
!   MPI_Bcast is used, so implicitly all processes are synchronized at this
!   step
!
!
!*******************************************************************************
    use com_mod
    use par_mod,only: numwfmem

    implicit none

    integer, intent(in) :: memstat

! Common array sizes used for communications
    integer, parameter :: d3_size1 = nxmax*nymax*nzmax
    integer, parameter :: d3_size2 = nxmax*nymax*nuvzmax
    integer, parameter :: d2_size1 = nxmax*nymax
    integer, parameter :: d2_size2 = nxmax*nymax*maxspec
    integer, parameter :: d2_size3 = nxmax*nymax
    integer, parameter :: d1_size1 = maxwf

    integer :: d3s1,d3s2,d2s1,d2s2

! li,ui    lower,upper indices of fields to be transfered (1-3, ui>=li) 
    integer :: li,ui

! First time routine is called the unchangeable fields will be transfered    
    logical :: first_call=.true.

!**********************************************************************

! Sizes of arrays transfered
    d3s1=d3_size1
    d3s2=d3_size2
    d2s1=d2_size1 
    d2s2=d2_size2

! Decide which field(s) need to be transfered
    if (memstat.ge.32) then ! distribute 2 fields, to 1st/2nd indices 
      li=1; ui=2
      d3s1=2*d3_size1
      d3s2=2*d3_size2
      d2s1=2*d2_size1 
      d2s2=2*d2_size2
    else if (memstat.eq.17) then ! distribute 1 field, place on 1st index
      li=1; ui=1
    else if (memstat.eq.18) then ! distribute 1 field, place on 2nd index
      li=2; ui=2
    else if (memstat.eq.19) then ! distribute 1 field, place on 3nd index
      li=3; ui=3
    else
      write(*,*) '#### mpi_mod::mpif_gf_send_vars> ERROR: ', &
           & 'wrong value of memstat, exiting ####', memstat
      stop
    end if


! Time for MPI communications
    if (mp_measure_time) call mpif_mtime('commtime',0)

! Send variables from getfield process (id_read) to other processes
!**********************************************************************

1331
! The non-reader processes need to know if cloud water was read.
Espen Sollum's avatar
Espen Sollum committed
1332
1333
1334
    call MPI_Bcast(readclouds,1,MPI_LOGICAL,id_read,MPI_COMM_WORLD,mp_ierr)
    if (mp_ierr /= 0) goto 600 

1335
1336
! Static fields/variables sent only at startup
    if (first_call) then
1337
      call MPI_Bcast(oro(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1338
      if (mp_ierr /= 0) goto 600 
1339
      call MPI_Bcast(excessoro(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1340
      if (mp_ierr /= 0) goto 600 
1341
      call MPI_Bcast(lsm(:,:),d2_size3,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1342
      if (mp_ierr /= 0) goto 600 
1343
      call MPI_Bcast(xlanduse(:,:,:),d2_size3*numclass,mp_sp,id_read,MPI_COMM_WORLD,mp_ierr)
1344
1345
1346
1347
1348
1349
1350
1351
      if (mp_ierr /= 0) goto 600 
      call MPI_Bcast(wftime,d1_size1,MPI_INTEGER,id_read,MPI_COMM_WORLD,mp_ierr)
      if (mp_ierr /=