1module midasMpi_mod
2 !
3 ! MODULE midasMpi_mod (prefix='mmpi' category='8. Low-level utilities and constants')
4 !
5 !:Purpose: Subroutine and public variables related to general aspects of mpi.
6 ! Also, subroutine and public variables related to the mpi decomposition
7 ! specific to the MIDAS code.
8 !
9 use utilities_mod
10 implicit none
11 save
12 private
13
14 ! public variables
15 public :: mmpi_myid,mmpi_nprocs,mmpi_myidx,mmpi_myidy,mmpi_npex,mmpi_npey,mmpi_numthread
16 public :: mmpi_comm_EW, mmpi_comm_NS, mmpi_comm_GRID, mmpi_doBarrier
17 public :: mmpi_datyp_real4, mmpi_datyp_real8, mmpi_datyp_int
18 ! public procedures
19 public :: mmpi_initialize,mmpi_getptopo
20 public :: mmpi_allreduce_sumreal8scalar,mmpi_allgather_string
21 public :: mmpi_allreduce_sumR8_1d, mmpi_allreduce_sumR8_2d
22 public :: mmpi_reduce_sumR8_1d, mmpi_reduce_sumR8_2d, mmpi_reduce_sumR8_3d
23 public :: mmpi_setup_latbands, mmpi_setup_lonbands
24 public :: mmpi_setup_m, mmpi_setup_n
25 public :: mmpi_setup_levels
26 public :: mmpi_setup_varslevels
27 public :: mmpi_myidXfromLon, mmpi_myidYfromLat
28
29 integer :: mmpi_myid = 0
30 integer :: mmpi_nprocs = 0
31 integer :: mmpi_myidx = 0
32 integer :: mmpi_myidy = 0
33 integer :: mmpi_npex = 0
34 integer :: mmpi_npey = 0
35 integer :: mmpi_numthread = 0
36
37 integer :: mmpi_comm_EW, mmpi_comm_NS, mmpi_comm_GRID
38 integer :: mmpi_datyp_real4, mmpi_datyp_real8, mmpi_datyp_int
39
40 logical :: mmpi_doBarrier = .true.
41
42 contains
43
44 subroutine mmpi_initialize()
45 implicit none
46
47 ! Locals:
48 integer :: mythread,numthread,omp_get_thread_num,omp_get_num_threads,rpn_comm_mype
49 integer :: ierr
50 integer :: rpn_comm_comm, rpn_comm_datyp
51
52 ! Namelist variables
53 integer :: npex ! number of MPI tasks in 'x' direction (set automatically by launch script)
54 integer :: npey ! number of MPI tasks in 'y' direction (set automatically by launch script)
55
56 ! Initilize MPI
57 npex=0
58 npey=0
59 call rpn_comm_init(mmpi_getptopo,mmpi_myid,mmpi_nprocs,npex,npey)
60 if(mmpi_nprocs.lt.1) then
61 mmpi_nprocs=1
62 mmpi_npex=1
63 mmpi_npey=1
64 mmpi_myid=0
65 mmpi_myidx=0
66 mmpi_myidy=0
67 else
68 ierr = rpn_comm_mype(mmpi_myid,mmpi_myidx,mmpi_myidy)
69 mmpi_npex=npex
70 mmpi_npey=npey
71 endif
72
73 write(*,*) 'mmpi_initialize: mmpi_myid, mmpi_myidx, mmpi_myidy = ', mmpi_myid, mmpi_myidx, mmpi_myidy
74
75 !$OMP PARALLEL PRIVATE(numthread,mythread)
76 mythread=omp_get_thread_num()
77 numthread=omp_get_num_threads()
78 if(mythread.eq.0) then
79 write(*,*) 'mmpi_initialize: NUMBER OF THREADS=',numthread
80 mmpi_numthread = numthread
81 end if
82 !$OMP END PARALLEL
83
84 ! create standard mpi handles to rpn_comm mpi communicators to facilitate
85 ! use of standard mpi routines
86 mmpi_comm_EW = rpn_comm_comm('EW')
87 mmpi_comm_NS = rpn_comm_comm('NS')
88 mmpi_comm_GRID = rpn_comm_comm('GRID')
89
90 mmpi_datyp_real4 = rpn_comm_datyp('MPI_REAL4')
91 mmpi_datyp_real8 = rpn_comm_datyp('MPI_REAL8')
92 mmpi_datyp_int = rpn_comm_datyp('MPI_INTEGER')
93
94 write(*,*) ' '
95 if(mmpi_doBarrier) then
96 write(*,*) 'mmpi_initialize: MPI_BARRIERs will be done to help with interpretation of timings'
97 else
98 write(*,*) 'mmpi_initialize: no MPI_BARRIERs will be done'
99 endif
100 write(*,*) ' '
101
102 end subroutine mmpi_initialize
103
104
105 subroutine mmpi_getptopo( npex, npey )
106
107 implicit none
108
109 ! Arguments:
110 integer, intent(out) :: npex
111 integer, intent(out) :: npey
112
113 ! Locals:
114 integer :: ierr
115 integer :: nulnam,fnom,fclos
116 namelist /ptopo/npex,npey
117
118 npex=1
119 npey=1
120
121 nulnam=0
122 ierr=fnom(nulnam,'ptopo_nml','FTN+SEQ+R/O',0)
123 if(ierr.ne.0) call utl_abort('mpi_getptopo: Error opening file ptopo_nml')
124 read(nulnam,nml=ptopo,iostat=ierr)
125 if(ierr.ne.0) call utl_abort('mpi_getptopo: Error reading namelist')
126 write(*,nml=ptopo)
127 ierr=fclos(nulnam)
128
129 end subroutine mmpi_getptopo
130
131
132 subroutine mmpi_allreduce_sumreal8scalar( sendRecvValue, comm )
133
134 implicit none
135
136 ! Arguments:
137 real(8), intent(inout) :: sendRecvValue ! value to be summed over all mpi tasks
138 character(len=*), intent(in) :: comm ! rpn_comm communicator
139
140 ! Locals:
141 integer :: nsize, ierr, root, rank
142 real(8), allocatable :: allvalues(:)
143
144 ! do a barrier so that timing on reduce operation is accurate
145 call utl_tmg_start(171,'low-level--mpi_allreduce_barr')
146 if(mmpi_doBarrier) call rpn_comm_barrier(comm,ierr)
147 call utl_tmg_stop(171)
148
149 call utl_tmg_start(170,'low-level--mpi_allreduce_sum8')
150
151 ! determine number of processors in the communicating group
152 call rpn_comm_size(comm,nsize,ierr)
153
154 ! determine where to gather the values: first task in group
155 call rpn_comm_rank(comm,rank,ierr)
156 call rpn_comm_allreduce(rank,root,1,"MPI_INTEGER","MPI_MIN",comm,ierr)
157
158 ! gather values to be added onto 1 processor
159 allocate(allvalues(nsize))
160 call rpn_comm_gather(sendRecvValue, 1, "MPI_DOUBLE_PRECISION", allvalues, 1, "MPI_DOUBLE_PRECISION", root, comm, ierr)
161
162 ! sum the values on the "root" mpi task and broadcast to group
163 if(rank.eq.root) sendRecvValue = sum(allvalues(:))
164 deallocate(allvalues)
165 call rpn_comm_bcast(sendRecvValue, 1, "MPI_DOUBLE_PRECISION", root, comm, ierr)
166
167 call utl_tmg_stop(170)
168
169 end subroutine mmpi_allreduce_sumreal8scalar
170
171
172 subroutine mmpi_allreduce_sumR8_1d( sendRecvVector, comm )
173 !
174 !:Purpose: Perform sum of 1d array over all MPI tasks.
175 !
176 implicit none
177
178 ! Arguments:
179 real(8) , intent(inout) :: sendRecvVector(:) ! 1-D vector to be summed over all mpi tasks
180 character(len=*), intent(in) :: comm ! rpn_comm communicator
181
182 ! Locals:
183 integer :: nprocs_mpi, numElements, ierr, root, rank
184 real(8), allocatable :: all_sendRecvVector(:,:)
185
186 ! do a barrier so that timing on reduce operation is accurate
187 call utl_tmg_start(171,'low-level--mpi_allreduce_barr')
188 if ( mmpi_doBarrier ) call rpn_comm_barrier(comm,ierr)
189 call utl_tmg_stop(171)
190
191 call utl_tmg_start(170,'low-level--mpi_allreduce_sum8')
192
193 numElements = size(sendRecvVector)
194
195 ! determine number of processors in the communicating group
196 call rpn_comm_size(comm,nprocs_mpi,ierr)
197
198 ! determine where to gather the values: first task in group
199 call rpn_comm_rank(comm,rank,ierr)
200 call rpn_comm_allreduce(rank,root,1,"mpi_integer","mpi_min",comm,ierr)
201
202 ! gather vectors to be added onto 1 processor
203 allocate(all_sendRecvVector(numElements,0:nprocs_mpi-1))
204 call rpn_comm_gather(sendRecvVector , numElements, "mpi_double_precision", &
205 all_sendRecvVector, numElements, "mpi_double_precision", &
206 root, comm, ierr)
207
208 ! sum the values on the "root" mpi task and broadcast to group
209 if ( rank == root ) sendRecvVector(:) = sum(all_sendRecvVector(:,:),2)
210 deallocate(all_sendRecvVector)
211 call rpn_comm_bcast(sendRecvVector, numElements, "mpi_double_precision", &
212 root, comm, ierr)
213
214 call utl_tmg_stop(170)
215
216 end subroutine mmpi_allreduce_sumR8_1d
217
218
219 subroutine mmpi_allreduce_sumR8_2d( sendRecvVector, comm )
220 !
221 !:Purpose: Perform sum of 2d array over all MPI tasks.
222 !
223 implicit none
224
225 ! Arguments:
226 real(8) , intent(inout) :: sendRecvVector(:,:) ! 2-D vector to be summed over all mpi tasks
227 character(len=*), intent(in) :: comm ! rpn_comm communicator
228
229 ! Locals:
230 integer :: nprocs_mpi, numElements1, numElements2, ierr, root, rank
231 real(8), allocatable :: all_sendRecvVector(:,:,:)
232
233 ! do a barrier so that timing on reduce operation is accurate
234 call utl_tmg_start(171,'low-level--mpi_allreduce_barr')
235 if ( mmpi_doBarrier ) call rpn_comm_barrier(comm,ierr)
236 call utl_tmg_stop(171)
237
238 call utl_tmg_start(170,'low-level--mpi_allreduce_sum8')
239
240 numElements1 = size(sendRecvVector,1)
241 numElements2 = size(sendRecvVector,2)
242
243 ! determine number of processors in the communicating group
244 call rpn_comm_size(comm,nprocs_mpi,ierr)
245
246 ! determine where to gather the values: first task in group
247 call rpn_comm_rank(comm,rank,ierr)
248 call rpn_comm_allreduce(rank,root,1,"mpi_integer","mpi_min",comm,ierr)
249
250 ! gather vectors to be added onto 1 processor
251 allocate(all_sendRecvVector(numElements1,numElements2,0:nprocs_mpi-1))
252 call rpn_comm_gather(sendRecvVector , numElements1*numElements2, "mpi_double_precision", &
253 all_sendRecvVector, numElements1*numElements2, "mpi_double_precision", &
254 root, comm, ierr)
255
256 ! sum the values on the "root" mpi task and broadcast to group
257 if ( rank == root ) sendRecvVector(:,:) = sum(all_sendRecvVector(:,:,:),3)
258 deallocate(all_sendRecvVector)
259 call rpn_comm_bcast(sendRecvVector, numElements1*numElements2, "mpi_double_precision", &
260 root, comm, ierr)
261
262 call utl_tmg_stop(170)
263
264 end subroutine mmpi_allreduce_sumR8_2d
265
266
267 subroutine mmpi_reduce_sumR8_1d( sendVector, recvVector, root, comm )
268 !
269 !:Purpose: Perform sum of 1d array over all MPI tasks.
270 !
271 implicit none
272
273 ! Arguments:
274 real(8) , intent(in) :: sendVector(:) ! 1-D vector to be summed over all mpi tasks
275 real(8) , intent(out) :: recvVector(:) ! 1-D vector to be summed over all mpi tasks
276 integer , intent(in) :: root ! mpi task id where data is put
277 character(len=*), intent(in) :: comm ! rpn_comm communicator
278
279 ! Locals:
280 integer :: nprocs_mpi, numElements, ierr, rank
281 real(8), allocatable :: all_sendRecvVector(:,:)
282
283 ! do a barrier so that timing on reduce operation is accurate
284 call utl_tmg_start(171,'low-level--mpi_allreduce_barr')
285 if ( mmpi_doBarrier ) call rpn_comm_barrier(comm,ierr)
286 call utl_tmg_stop(171)
287
288 call utl_tmg_start(170,'low-level--mpi_allreduce_sum8')
289
290 numElements = size(sendVector)
291
292 ! determine number of processors in the communicating group
293 call rpn_comm_size(comm,nprocs_mpi,ierr)
294
295 ! determine rank of group
296 call rpn_comm_rank(comm,rank,ierr)
297
298 ! gather vectors to be added onto 1 processor
299 if ( rank == root ) then
300 allocate(all_sendRecvVector(numElements,0:nprocs_mpi-1))
301 else
302 allocate(all_sendRecvVector(1,1))
303 end if
304 call rpn_comm_gather(sendVector , numElements, "mpi_double_precision", &
305 all_sendRecvVector, numElements, "mpi_double_precision", &
306 root, comm, ierr)
307
308 ! sum the values on the "root" mpi task
309 if ( rank == root ) recvVector(:) = sum(all_sendRecvVector(:,:),2)
310 deallocate(all_sendRecvVector)
311
312 call utl_tmg_stop(170)
313
314 end subroutine mmpi_reduce_sumR8_1d
315
316
317 subroutine mmpi_reduce_sumR8_2d( sendVector, recvVector, root, comm )
318 !
319 !:Purpose: Perform sum of 2d array over all MPI tasks.
320 !
321 implicit none
322
323 ! Arguments:
324 real(8) , intent(in) :: sendVector(:,:) ! 2-D vector to be summed over all mpi tasks
325 real(8) , intent(out) :: recvVector(:,:) ! 2-D vector to be summed over all mpi tasks
326 integer , intent(in) :: root ! mpi task id where data will be put
327 character(len=*), intent(in) :: comm ! rpn_comm communicator
328
329 ! Locals:
330 integer :: nprocs_mpi, numElements1, numElements2, ierr, rank
331 real(8), allocatable :: all_sendRecvVector(:,:,:)
332
333 ! do a barrier so that timing on reduce operation is accurate
334 call utl_tmg_start(171,'low-level--mpi_allreduce_barr')
335 if ( mmpi_doBarrier ) call rpn_comm_barrier(comm,ierr)
336 call utl_tmg_stop(171)
337
338 call utl_tmg_start(170,'low-level--mpi_allreduce_sum8')
339
340 numElements1 = size(sendVector,1)
341 numElements2 = size(sendVector,2)
342
343 ! determine number of processors in the communicating group
344 call rpn_comm_size(comm,nprocs_mpi,ierr)
345
346 ! determine rank of group
347 call rpn_comm_rank(comm,rank,ierr)
348
349 ! gather vectors to be added onto 1 processor
350 if ( rank == root ) then
351 allocate(all_sendRecvVector(numElements1,numElements2,0:nprocs_mpi-1))
352 else
353 allocate(all_sendRecvVector(1,1,1))
354 end if
355 call rpn_comm_gather(sendVector , numElements1*numElements2, "mpi_double_precision", &
356 all_sendRecvVector, numElements1*numElements2, "mpi_double_precision", &
357 root, comm, ierr)
358
359 ! sum the values on the "root" mpi task
360 if ( rank == root ) recvVector(:,:) = sum(all_sendRecvVector(:,:,:),3)
361 deallocate(all_sendRecvVector)
362
363 call utl_tmg_stop(170)
364
365 end subroutine mmpi_reduce_sumR8_2d
366
367
368 subroutine mmpi_reduce_sumR8_3d( sendVector, recvVector, root, comm )
369 !
370 !:Purpose: Perform sum of 3d array over all MPI tasks.
371 !
372 implicit none
373
374 ! Arguments:
375 real(8) , intent(in) :: sendVector(:,:,:) ! 3-D vector to be summed over all mpi tasks
376 real(8) , intent(out) :: recvVector(:,:,:) ! 3-D vector to be summed over all mpi tasks
377 integer , intent(in) :: root ! mpi task id where data is put
378 character(len=*), intent(in) :: comm ! rpn_comm communicator
379
380 ! Locals:
381 integer :: nprocs_mpi, numElements1, numElements2, numElements3, ierr, rank
382 real(8), allocatable :: all_sendRecvVector(:,:,:,:)
383
384 ! do a barrier so that timing on reduce operation is accurate
385 call utl_tmg_start(171,'low-level--mpi_allreduce_barr')
386 if ( mmpi_doBarrier ) call rpn_comm_barrier(comm,ierr)
387 call utl_tmg_stop(171)
388
389 call utl_tmg_start(170,'low-level--mpi_allreduce_sum8')
390
391 numElements1 = size(sendVector,1)
392 numElements2 = size(sendVector,2)
393 numElements3 = size(sendVector,3)
394
395 ! determine number of processors in the communicating group
396 call rpn_comm_size(comm,nprocs_mpi,ierr)
397
398 ! determine rank of group
399 call rpn_comm_rank(comm,rank,ierr)
400
401 ! gather vectors to be added onto 1 processor
402 if ( rank == root ) then
403 allocate(all_sendRecvVector(numElements1,numElements2,numElements3,0:nprocs_mpi-1))
404 else
405 allocate(all_sendRecvVector(1,1,1,1))
406 end if
407 call rpn_comm_gather(sendVector , numElements1*numElements2*numElements3, "mpi_double_precision", &
408 all_sendRecvVector, numElements1*numElements2*numElements3, "mpi_double_precision", &
409 root, comm, ierr)
410
411 ! sum the values on the "root" mpi task
412 if ( rank == root ) recvVector(:,:,:) = sum(all_sendRecvVector(:,:,:,:),4)
413 deallocate(all_sendRecvVector)
414
415 call utl_tmg_stop(170)
416
417 end subroutine mmpi_reduce_sumR8_3d
418
419
420 subroutine mmpi_allgather_string( str_list, str_list_all, nlist, nchar, nproc, comm, ierr )
421 !
422 !:Purpose: Performs the MPI 'allgather' routine for an array of strings
423 !
424 implicit none
425
426 ! Arguments:
427 character(len=nchar), intent(in) :: str_list(nlist)
428 character(len=*) , intent(in) :: comm
429 integer , intent(in) :: nlist
430 integer , intent(in) :: nchar
431 integer , intent(in) :: nproc
432 character(len=nchar), intent(out) :: str_list_all(nlist,nproc)
433 integer , intent(out) :: ierr
434
435 ! Locals:
436 integer :: num_list(nlist*nchar),num_list_all(nlist*nchar,nproc)
437 integer :: ilist,ichar,iproc
438
439 ! Convert strings to integer sequences
440
441 do ilist=1,nlist
442 do ichar=1,nchar
443 num_list((ilist-1)*nchar+ichar) = iachar(str_list(ilist)(ichar:ichar))
444 end do
445 end do
446
447 ! Perform allgather with converted integer sequences
448
449 call rpn_comm_allgather(num_list,nlist*nchar,"MPI_INTEGER",num_list_all,nlist*nchar,"MPI_INTEGER",comm,ierr)
450
451 ! Convert integer sequences to stnid character strings
452
453 do iproc=1,nproc
454 do ilist=1,nlist
455 do ichar=1,nchar
456 str_list_all(ilist,iproc)(ichar:ichar) = achar(num_list_all((ilist-1)*nchar+ichar,iproc))
457 end do
458 end do
459 end do
460
461 end subroutine mmpi_allgather_string
462
463
464 subroutine mmpi_setup_latbands(nj, latPerPE, latPerPEmax, myLatBeg, myLatEnd, &
465 myLatHalfBeg_opt, myLatHalfEnd_opt, divisible_opt)
466 !:Purpose: compute parameters that define the mpi distribution of
467 ! latitudes over tasks in Y direction (npey)
468 implicit none
469
470 ! Arguments:
471 integer , intent(in) :: nj
472 integer , intent(out) :: latPerPE
473 integer , intent(out) :: latPerPEmax
474 integer , intent(out) :: myLatBeg
475 integer , intent(out) :: myLatEnd
476 integer, optional, intent(out) :: myLatHalfBeg_opt
477 integer, optional, intent(out) :: myLatHalfEnd_opt
478 logical, optional, intent(out) :: divisible_opt
479
480 ! Locals:
481 integer :: latPerPEmin, njlath, ierr
482 logical, save :: firstCall = .true.
483
484 latPerPEmin = floor(real(nj) / real(mmpi_npey))
485 myLatBeg = 1 + (mmpi_myidy * latPerPEmin)
486 if( mmpi_myidy < (mmpi_npey-1) ) then
487 myLatEnd = (1 + mmpi_myidy) * latPerPEmin
488 else
489 myLatEnd = nj
490 end if
491 latPerPE = myLatEnd - myLatBeg + 1
492 call rpn_comm_allreduce(latPerPE,latPerPEmax,1,'MPI_INTEGER','MPI_MAX','NS',ierr)
493
494 if( firstCall ) then
495 write(*,'(a,4i8)') 'mmpi_setup_latbands: latPerPE, latPerPEmax, myLatBeg, myLatEnd = ', &
496 latPerPE, latPerPEmax, myLatBeg, myLatEnd
497 firstCall = .false.
498 end if
499
500 if (present(myLatHalfBeg_opt).and.present(myLatHalfEnd_opt)) then
501 njlath = (nj + 1) / 2
502 if (myLatBeg <= njlath .and. myLatEnd <= njlath) then
503 myLatHalfBeg_opt = myLatBeg
504 myLatHalfEnd_opt = myLatEnd
505 elseif (myLatBeg >= njlath .and. myLatEnd >= njlath) then
506 myLatHalfBeg_opt = 1 + nj - myLatEnd
507 myLatHalfEnd_opt = 1 + nj - myLatBeg
508 else
509 myLatHalfBeg_opt = min(myLatBeg, 1 + nj - myLatEnd)
510 myLatHalfEnd_opt = njlath
511 end if
512 end if
513
514 if( present(divisible_opt) ) then
515 divisible_opt = (latPerPEmin * mmpi_npey == nj)
516 end if
517
518 end subroutine mmpi_setup_latbands
519
520
521 function mmpi_myidYfromLat(latIndex, nj) result(IP_y)
522 !:Purpose: use same logic as setup_latbands to compute myidy
523 ! corresponding to a latitude grid index
524 implicit none
525
526 ! Arguments:
527 integer, intent(in) :: latIndex
528 integer, intent(in) :: nj
529 ! Result:
530 integer :: IP_y
531
532 IP_y = (latIndex-1) / floor( real(nj) / real(mmpi_npey) )
533 IP_y = min( mmpi_npey-1, IP_y )
534
535 end function mmpi_myidYfromLat
536
537
538 subroutine mmpi_setup_lonbands(ni, lonPerPE, lonPerPEmax, myLonBeg, myLonEnd, divisible_opt)
539 !:Purpose: compute parameters that define the mpi distribution of
540 ! longitudes over tasks in X direction (npex)
541 implicit none
542
543 ! Arguments:
544 integer , intent(in) :: ni
545 integer , intent(out) :: lonPerPE
546 integer , intent(out) :: lonPerPEmax
547 integer , intent(out) :: myLonBeg
548 integer , intent(out) :: myLonEnd
549 logical, optional, intent(out) :: divisible_opt
550
551 ! Locals:
552 integer :: lonPerPEmin, ierr
553 logical, save :: firstCall = .true.
554
555 lonPerPEmin = floor(real(ni) / real(mmpi_npex))
556 myLonBeg = 1 + (mmpi_myidx * lonPerPEmin)
557 if( mmpi_myidx < (mmpi_npex-1) ) then
558 myLonEnd = (1 + mmpi_myidx) * lonPerPEmin
559 else
560 myLonEnd = ni
561 end if
562 lonPerPE = myLonEnd - myLonBeg + 1
563 call rpn_comm_allreduce(lonPerPE,lonPerPEmax,1,'MPI_INTEGER','MPI_MAX','EW',ierr)
564
565 if( firstCall ) then
566 write(*,'(a,4i8)') 'mmpi_setup_lonbands: lonPerPE, lonPerPEmax, myLonBeg, myLonEnd = ', &
567 lonPerPE, lonPerPEmax, myLonBeg, myLonEnd
568 firstCall = .false.
569 end if
570
571 if( present(divisible_opt) ) then
572 divisible_opt = (lonPerPEmin * mmpi_npex == ni)
573 end if
574
575 end subroutine mmpi_setup_lonbands
576
577
578 function mmpi_myidXfromLon(lonIndex, ni) result(IP_x)
579 !:Purpose: use same logic as setup_lonbands to compute myidx
580 ! corresponding to a longitude grid index
581 implicit none
582
583 ! Arguments:
584 integer, intent(in) :: lonIndex
585 integer, intent(in) :: ni
586 ! Result:
587 integer :: IP_x
588
589 IP_x = (lonIndex-1) / floor( real(ni) / real(mmpi_npex) )
590 IP_x = min( mmpi_npex-1, IP_x )
591
592 end function mmpi_myidXfromLon
593
594
595 subroutine mmpi_setup_m(ntrunc, mymBeg, mymEnd, mymSkip, mymCount)
596 !:Purpose: compute parameters that define the mpi distribution of
597 ! wavenumber m over tasks in Y direction (npey)
598 implicit none
599
600 ! Arguments:
601 integer, intent(in) :: ntrunc
602 integer, intent(out) :: mymBeg
603 integer, intent(out) :: mymEnd
604 integer, intent(out) :: mymSkip
605 integer, intent(out) :: mymCount
606
607 ! Locals:
608 integer :: jm
609
610 mymBeg = mmpi_myidy
611 mymEnd = ntrunc
612 mymSkip = mmpi_npey
613 mymCount = 0
614 do jm = mymBeg, mymEnd, mymSkip
615 mymCount = mymCount + 1
616 end do
617
618 write(*,'(a,4i8)') 'mmpi_setup_m: mymBeg, mymEnd, mymSkip, mymCount = ', mymBeg, mymEnd, mymSkip, mymCount
619
620 end subroutine mmpi_setup_m
621
622
623 subroutine mmpi_setup_n(ntrunc, mynBeg, mynEnd, mynSkip, mynCount)
624 !:Purpose: compute parameters that define the mpi distribution of
625 ! wavenumber n over tasks in X direction (npex)
626 implicit none
627
628 ! Arguments:
629 integer, intent(in) :: ntrunc
630 integer, intent(out) :: mynBeg
631 integer, intent(out) :: mynEnd
632 integer, intent(out) :: mynSkip
633 integer, intent(out) :: mynCount
634
635 ! Locals:
636 integer :: jn
637
638 mynBeg = mmpi_myidx
639 mynEnd = ntrunc
640 mynSkip = mmpi_npex
641 mynCount = 0
642 do jn = mynBeg, mynEnd, mynSkip
643 mynCount = mynCount + 1
644 end do
645
646 write(*,'(a,4i8)') 'mmpi_setup_n: mynBeg, mynEnd, mynSkip, mynCount = ', mynBeg, mynEnd, mynSkip, mynCount
647
648 end subroutine mmpi_setup_n
649
650
651 subroutine mmpi_setup_levels(numlevels, myLevBeg, myLevEnd, myLevCount)
652 !:Purpose: compute parameters that define the mpi distribution of
653 ! levels over tasks in X direction (npex)
654 implicit none
655
656 ! Arguments:
657 integer, intent(in) :: numlevels
658 integer, intent(out) :: myLevBeg
659 integer, intent(out) :: myLevEnd
660 integer, intent(out) :: myLevCount
661
662 ! Locals:
663 integer :: jlev
664 integer :: jproc
665 integer :: factor
666 integer :: myLevCounts(mmpi_npex)
667
668 ! when possible, always divide into even number of levels per MPI task
669 if(mod(numlevels, 2) /= 0) then
670 write(*,*) 'mmpi_setup_levels: total number of levels is not even, now=', numlevels
671 write(*,*) ' therefore, if global grid, may not be able to do '
672 write(*,*) ' transforms of vor/div <-> u/v'
673 factor = 1
674 else
675 factor = 2
676 end if
677
678 myLevCounts(:) = 0
679 do jproc = 1, mmpi_npex
680 do jlev = jproc, (numlevels / factor), mmpi_npex
681 myLevCounts(jproc) = myLevCounts(jproc) + 1
682 end do
683 end do
684 do jproc = 1, mmpi_npex
685 myLevCounts(jproc) = myLevCounts(jproc) * factor
686 end do
687
688 myLevCount = myLevCounts(mmpi_myidx + 1)
689
690 if (myLevCount > 0) then
691 myLevBeg = 1
692 do jproc = 1, mmpi_myidx
693 myLevBeg = myLevBeg + myLevCounts(jproc)
694 end do
695 myLevEnd = myLevBeg + myLevCount - 1
696 else
697 myLevBeg = 1
698 myLevEnd = 0
699 end if
700
701 write(*,'(a,3i8)') 'mmpi_setup_levels: myLevBeg, myLevEnd, myLevCount = ', &
702 myLevBeg, myLevEnd, myLevCount
703
704 end subroutine mmpi_setup_levels
705
706
707 subroutine mmpi_setup_varslevels(numk, mykBeg, mykEnd, mykCount)
708 !:Purpose: compute parameters that define the mpi distribution of
709 ! variables/levels (i.e. 1->nk) over all tasks (nprocs)
710 implicit none
711
712 ! Arguments:
713 integer, intent(in) :: numk
714 integer, intent(out) :: mykBeg
715 integer, intent(out) :: mykEnd
716 integer, intent(out) :: mykCount
717
718 ! Locals:
719 integer :: jk
720 integer :: jproc
721 integer :: mykCounts(mmpi_nprocs)
722
723 mykCounts(:) = 0
724 do jproc = 1, mmpi_nprocs
725 do jk = jproc, numk, mmpi_nprocs
726 mykCounts(jproc) = mykCounts(jproc) + 1
727 end do
728 end do
729
730 mykCount = mykCounts(mmpi_myid + 1)
731
732 mykBeg = 1
733 do jproc = 1, mmpi_myid
734 mykBeg = mykBeg + mykCounts(jproc)
735 end do
736 mykEnd = mykBeg + mykCount - 1
737
738 write(*,'(a,3i8)') 'mmpi_setup_varslevels: mykBeg, mykEnd, mykCount = ', &
739 mykBeg, mykEnd, mykCount
740
741 end subroutine mmpi_setup_varslevels
742
743end module midasMpi_mod