source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/mct/mpeu/m_FcComms.F90 @ 4775

Last change on this file since 4775 was 4775, checked in by aclsce, 4 years ago
  • Imported oasis3-mct from Cerfacs svn server (not suppotred anymore).

The version has been extracted from https://oasis3mct.cerfacs.fr/svn/branches/OASIS3-MCT_2.0_branch/oasis3-mct@1818

File size: 19.5 KB
Line 
1!BOP -------------------------------------------------------------------
2!
3! !MODULE: m_FcComms - MPI collective communication operators
4!                      with explict flow control
5!
6! !DESCRIPTION:
7!
8! This module includes implementations of MPI collective operators that
9! have proven problematic on certain systems when run at scale. By
10! introducing additonal flow control, these problems (exhausting internal
11! system resources) can be avoided. These routines were ported from
12! the Community Atmosphere Model's spmd_utils.F90.
13!
14! !INTERFACE:
15!
16! Workaround for performance issue with rsend on cray systems with
17! gemini interconnect
18!
19#ifdef _NO_MPI_RSEND
20#define MPI_RSEND MPI_SEND
21#define mpi_rsend mpi_send
22#define MPI_IRSEND MPI_ISEND
23#define mpi_irsend mpi_isend
24#endif
25
26 module m_FcComms
27
28      implicit none
29
30      private   ! except
31
32      public :: fc_gather_int  ! flow control version of mpi_gather for integer vectors
33      public :: fc_gather_fp   ! flow control version of mpi_gather for FP vectors
34      public :: fc_gatherv_int ! flow control version of mpi_gatherv for integer vectors
35      public :: fc_gatherv_fp  ! flow control version of mpi_gatherv for integer vectors
36      public :: get_fcblocksize ! get current value of max_gather_block_size
37      public :: set_fcblocksize ! set current value of max_gather_block_size
38
39
40! !REVISION HISTORY:
41! 30Jan09 - P.H. Worley <worleyph@ornl.gov> - imported routines
42!           from CAM's spmd_utils to create this module.
43
44  integer, public :: max_gather_block_size = 64
45  character(len=*),parameter :: myname='MCT(MPEU)::m_FcComms'
46
47 contains
48
49!BOP -------------------------------------------------------------------
50!
51! !IROUTINE: fc_gather_int - Gather an array of type integer
52!
53! !DESCRIPTION:
54! This routine gathers a {\em distributed} array of type {\em integer}
55! to the {\tt root} process. Explicit handshaking messages are used
56! to control the number of processes communicating with the root
57! at any one time.
58!
59! If flow_cntl optional parameter
60!    < 0 : use MPI_Gather
61!    >= 0: use point-to-point with handshaking messages and
62!          preposting receive requests up to
63!          min(max(1,flow_cntl),max_gather_block_size)
64!          ahead if optional flow_cntl parameter is present.
65!          Otherwise, max_gather_block_size is used in its place.
66!    Default value is max_gather_block_size.
67! !INTERFACE:
68!
69   subroutine fc_gather_int (sendbuf, sendcnt, sendtype, &
70                             recvbuf, recvcnt, recvtype, &
71                             root, comm, flow_cntl )
72!
73! !USES:
74!
75      use m_die
76      use m_mpif90
77!
78! !INPUT PARAMETERS:
79!
80      integer, intent(in) :: sendbuf(*)
81      integer, intent(in) :: sendcnt
82      integer, intent(in) :: sendtype
83      integer, intent(in) :: recvcnt
84      integer, intent(in) :: recvtype
85      integer, intent(in) :: root
86      integer, intent(in) :: comm
87      integer, optional, intent(in) :: flow_cntl
88
89! !OUTPUT PARAMETERS:
90!
91      integer, intent(out) :: recvbuf(*)
92
93! !REVISION HISTORY:
94! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
95!EOP ___________________________________________________________________
96
97 character(len=*),parameter :: myname_=myname//'::fc_gather_int'
98
99 integer :: signal
100 logical fc_gather         ! use explicit flow control?
101 integer gather_block_size ! number of preposted receive requests
102
103 integer :: mytid, mysize, mtag, p, i, count, displs
104 integer :: preposts, head, tail
105 integer :: rcvid(max_gather_block_size)
106 integer :: status(MP_STATUS_SIZE)
107 integer :: ier ! MPI error code
108
109 signal = 1
110 if ( present(flow_cntl) ) then
111    if (flow_cntl >= 0) then
112       gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
113       fc_gather = .true.
114    else
115       fc_gather = .false.
116    endif
117 else
118    gather_block_size = max(1,max_gather_block_size)
119    fc_gather = .true.
120 endif
121
122 if (fc_gather) then
123 
124    call mpi_comm_rank (comm, mytid, ier)
125    call mpi_comm_size (comm, mysize, ier)
126    mtag = 0
127    if (root .eq. mytid) then
128
129! prepost gather_block_size irecvs, and start receiving data
130       preposts = min(mysize-1, gather_block_size)
131       head = 0
132       count = 0
133       do p=0, mysize-1
134          if (p .ne. root) then
135             if (recvcnt > 0) then
136                count = count + 1
137                if (count > preposts) then
138                   tail = mod(head,preposts) + 1
139                   call mpi_wait (rcvid(tail), status, ier)
140                end if
141                head = mod(head,preposts) + 1
142                displs = p*recvcnt
143                call mpi_irecv ( recvbuf(displs+1), recvcnt, &
144                                 recvtype, p, mtag, comm, rcvid(head), &
145                                 ier )
146                call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
147             end if
148          end if
149       end do
150
151! copy local data
152       displs = mytid*recvcnt
153       do i=1,sendcnt
154          recvbuf(displs+i) = sendbuf(i)
155       enddo
156
157! wait for final data
158       do i=1,min(count,preposts)
159          call mpi_wait (rcvid(i), status, ier)
160       enddo
161
162    else
163
164       if (sendcnt > 0) then
165          call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
166                          status, ier )
167          call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
168                           comm, ier )
169       end if
170
171    endif
172    if (ier /= 0) then
173       call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
174    end if
175
176 else
177 
178    call mpi_gather (sendbuf, sendcnt, sendtype, &
179                     recvbuf, recvcnt, recvtype, &
180                     root, comm, ier)
181    if (ier /= 0) then
182       call MP_perr_die(myname_,':: MPI_GATHER',ier)
183    end if
184
185 endif
186
187 return
188 end subroutine fc_gather_int
189
190!BOP -------------------------------------------------------------------
191!
192! !IROUTINE: fc_gather_fp - Gather an array of type FP
193!
194! !DESCRIPTION:
195! This routine gathers a {\em distributed} array of type {\em FP} to
196! the {\tt root} process. Explicit handshaking messages are used
197! to control the number of processes communicating with the root
198! at any one time.
199!
200! If flow_cntl optional parameter
201!    < 0 : use MPI_Gather
202!    >= 0: use point-to-point with handshaking messages and
203!          preposting receive requests up to
204!          min(max(1,flow_cntl),max_gather_block_size)
205!          ahead if optional flow_cntl parameter is present.
206!          Otherwise, max_gather_block_size is used in its place.
207!    Default value is max_gather_block_size.
208! !INTERFACE:
209!
210   subroutine fc_gather_fp (sendbuf, sendcnt, sendtype, &
211                            recvbuf, recvcnt, recvtype, &
212                             root, comm, flow_cntl )
213!
214! !USES:
215!
216      use m_realkinds, only : FP
217      use m_die
218      use m_mpif90
219!
220! !INPUT PARAMETERS:
221!
222      real (FP), intent(in)  :: sendbuf(*)
223      integer, intent(in) :: sendcnt
224      integer, intent(in) :: sendtype
225      integer, intent(in) :: recvcnt
226      integer, intent(in) :: recvtype
227      integer, intent(in) :: root
228      integer, intent(in) :: comm
229      integer, optional, intent(in) :: flow_cntl
230
231! !OUTPUT PARAMETERS:
232!
233      real (FP), intent(out) :: recvbuf(*)
234
235! !REVISION HISTORY:
236! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
237!EOP ___________________________________________________________________
238
239 character(len=*),parameter :: myname_=myname//'::fc_gather_fp'
240
241 real (FP) :: signal
242 logical fc_gather         ! use explicit flow control?
243 integer gather_block_size ! number of preposted receive requests
244
245 integer :: mytid, mysize, mtag, p, i, count, displs
246 integer :: preposts, head, tail
247 integer :: rcvid(max_gather_block_size)
248 integer :: status(MP_STATUS_SIZE)
249 integer :: ier ! MPI error code
250
251 signal = 1.0
252 if ( present(flow_cntl) ) then
253    if (flow_cntl >= 0) then
254       gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
255       fc_gather = .true.
256    else
257       fc_gather = .false.
258    endif
259 else
260    gather_block_size = max(1,max_gather_block_size)
261    fc_gather = .true.
262 endif
263
264 if (fc_gather) then
265 
266    call mpi_comm_rank (comm, mytid, ier)
267    call mpi_comm_size (comm, mysize, ier)
268    mtag = 0
269    if (root .eq. mytid) then
270
271! prepost gather_block_size irecvs, and start receiving data
272       preposts = min(mysize-1, gather_block_size)
273       head = 0
274       count = 0
275       do p=0, mysize-1
276          if (p .ne. root) then
277             if (recvcnt > 0) then
278                count = count + 1
279                if (count > preposts) then
280                   tail = mod(head,preposts) + 1
281                   call mpi_wait (rcvid(tail), status, ier)
282                end if
283                head = mod(head,preposts) + 1
284                displs = p*recvcnt
285                call mpi_irecv ( recvbuf(displs+1), recvcnt, &
286                                 recvtype, p, mtag, comm, rcvid(head), &
287                                 ier )
288                call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
289             end if
290          end if
291       end do
292
293! copy local data
294       displs = mytid*recvcnt
295       do i=1,sendcnt
296          recvbuf(displs+i) = sendbuf(i)
297       enddo
298
299! wait for final data
300       do i=1,min(count,preposts)
301          call mpi_wait (rcvid(i), status, ier)
302       enddo
303
304    else
305
306       if (sendcnt > 0) then
307          call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
308                          status, ier )
309          call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
310                           comm, ier )
311       end if
312
313    endif
314    if (ier /= 0) then
315       call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
316    end if
317
318 else
319 
320    call mpi_gather (sendbuf, sendcnt, sendtype, &
321                     recvbuf, recvcnt, recvtype, &
322                      root, comm, ier)
323    if (ier /= 0) then
324       call MP_perr_die(myname_,':: MPI_GATHER',ier)
325    end if
326
327 endif
328
329 return
330 end subroutine fc_gather_fp
331
332!BOP -------------------------------------------------------------------
333!
334! !IROUTINE: fc_gatherv_int - Gather an array of type integer
335!
336! !DESCRIPTION:
337! This routine gathers a {\em distributed} array of type {\em integer}
338! to the {\tt root} process. Explicit handshaking messages are used
339! to control the number of processes communicating with the root
340! at any one time.
341!
342! If flow_cntl optional parameter
343!    < 0 : use MPI_Gatherv
344!    >= 0: use point-to-point with handshaking messages and
345!          preposting receive requests up to
346!          min(max(1,flow_cntl),max_gather_block_size)
347!          ahead if optional flow_cntl parameter is present.
348!          Otherwise, max_gather_block_size is used in its place.
349!    Default value is max_gather_block_size.
350! !INTERFACE:
351!
352   subroutine fc_gatherv_int (sendbuf, sendcnt, sendtype, &
353                              recvbuf, recvcnts, displs, recvtype, &
354                              root, comm, flow_cntl )
355!
356! !USES:
357!
358      use m_die
359      use m_mpif90
360!
361! !INPUT PARAMETERS:
362!
363      integer, intent(in) :: sendbuf(*)
364      integer, intent(in) :: sendcnt
365      integer, intent(in) :: sendtype
366      integer, intent(in) :: recvcnts(*)
367      integer, intent(in) :: displs(*)
368      integer, intent(in) :: recvtype
369      integer, intent(in) :: root
370      integer, intent(in) :: comm
371      integer, optional, intent(in) :: flow_cntl
372
373! !OUTPUT PARAMETERS:
374!
375      integer, intent(out) :: recvbuf(*)
376
377! !REVISION HISTORY:
378! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
379!EOP ___________________________________________________________________
380
381 character(len=*),parameter :: myname_=myname//'::fc_gatherv_int'
382
383 integer :: signal
384 logical fc_gather         ! use explicit flow control?
385 integer gather_block_size ! number of preposted receive requests
386
387 integer :: mytid, mysize, mtag, p, q, i, count
388 integer :: preposts, head, tail
389 integer :: rcvid(max_gather_block_size)
390 integer :: status(MP_STATUS_SIZE)
391 integer :: ier ! MPI error code
392
393 signal = 1
394 if ( present(flow_cntl) ) then
395    if (flow_cntl >= 0) then
396       gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
397       fc_gather = .true.
398    else
399       fc_gather = .false.
400    endif
401 else
402    gather_block_size = max(1,max_gather_block_size)
403    fc_gather = .true.
404 endif
405
406 if (fc_gather) then
407 
408    call mpi_comm_rank (comm, mytid, ier)
409    call mpi_comm_size (comm, mysize, ier)
410    mtag = 0
411    if (root .eq. mytid) then
412
413! prepost gather_block_size irecvs, and start receiving data
414       preposts = min(mysize-1, gather_block_size)
415       head = 0
416       count = 0
417       do p=0, mysize-1
418          if (p .ne. root) then
419             q = p+1
420             if (recvcnts(q) > 0) then
421                count = count + 1
422                if (count > preposts) then
423                   tail = mod(head,preposts) + 1
424                   call mpi_wait (rcvid(tail), status, ier)
425                end if
426                head = mod(head,preposts) + 1
427                call mpi_irecv ( recvbuf(displs(q)+1), recvcnts(q), &
428                                 recvtype, p, mtag, comm, rcvid(head), &
429                                 ier )
430                call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
431             end if
432          end if
433       end do
434
435! copy local data
436       q = mytid+1
437       do i=1,sendcnt
438          recvbuf(displs(q)+i) = sendbuf(i)
439       enddo
440
441! wait for final data
442       do i=1,min(count,preposts)
443          call mpi_wait (rcvid(i), status, ier)
444       enddo
445
446    else
447
448       if (sendcnt > 0) then
449          call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
450                          status, ier )
451          call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
452                           comm, ier )
453       end if
454
455    endif
456    if (ier /= 0) then
457       call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
458    end if
459
460 else
461 
462    call mpi_gatherv (sendbuf, sendcnt, sendtype, &
463                      recvbuf, recvcnts, displs, recvtype, &
464                      root, comm, ier)
465    if (ier /= 0) then
466       call MP_perr_die(myname_,':: MPI_GATHERV',ier)
467    end if
468
469 endif
470
471 return
472 end subroutine fc_gatherv_int
473
474!BOP -------------------------------------------------------------------
475!
476! !IROUTINE: fc_gatherv_fp - Gather an array of type FP
477!
478! !DESCRIPTION:
479! This routine gathers a {\em distributed} array of type {\em FP} to
480! the {\tt root} process. Explicit handshaking messages are used
481! to control the number of processes communicating with the root
482! at any one time.
483!
484! If flow_cntl optional parameter
485!    < 0 : use MPI_Gatherv
486!    >= 0: use point-to-point with handshaking messages and
487!          preposting receive requests up to
488!          min(max(1,flow_cntl),max_gather_block_size)
489!          ahead if optional flow_cntl parameter is present.
490!          Otherwise, max_gather_block_size is used in its place.
491!    Default value is max_gather_block_size.
492! !INTERFACE:
493!
494   subroutine fc_gatherv_fp (sendbuf, sendcnt, sendtype, &
495                             recvbuf, recvcnts, displs, recvtype, &
496                             root, comm, flow_cntl )
497!
498! !USES:
499!
500      use m_realkinds, only : FP
501      use m_die
502      use m_mpif90
503!
504! !INPUT PARAMETERS:
505!
506      real (FP), intent(in)  :: sendbuf(*)
507      integer, intent(in) :: sendcnt
508      integer, intent(in) :: sendtype
509      integer, intent(in) :: recvcnts(*)
510      integer, intent(in) :: displs(*)
511      integer, intent(in) :: recvtype
512      integer, intent(in) :: root
513      integer, intent(in) :: comm
514      integer, optional, intent(in) :: flow_cntl
515
516! !OUTPUT PARAMETERS:
517!
518      real (FP), intent(out) :: recvbuf(*)
519
520! !REVISION HISTORY:
521! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
522!EOP ___________________________________________________________________
523
524 character(len=*),parameter :: myname_=myname//'::fc_gatherv_fp'
525
526 real (FP) :: signal
527 logical fc_gather         ! use explicit flow control?
528 integer gather_block_size ! number of preposted receive requests
529
530 integer :: mytid, mysize, mtag, p, q, i, count
531 integer :: preposts, head, tail
532 integer :: rcvid(max_gather_block_size)
533 integer :: status(MP_STATUS_SIZE)
534 integer :: ier ! MPI error code
535
536 signal = 1.0
537 if ( present(flow_cntl) ) then
538    if (flow_cntl >= 0) then
539       gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
540       fc_gather = .true.
541    else
542       fc_gather = .false.
543    endif
544 else
545    gather_block_size = max(1,max_gather_block_size)
546    fc_gather = .true.
547 endif
548
549 if (fc_gather) then
550 
551    call mpi_comm_rank (comm, mytid, ier)
552    call mpi_comm_size (comm, mysize, ier)
553    mtag = 0
554    if (root .eq. mytid) then
555
556! prepost gather_block_size irecvs, and start receiving data
557       preposts = min(mysize-1, gather_block_size)
558       head = 0
559       count = 0
560       do p=0, mysize-1
561          if (p .ne. root) then
562             q = p+1
563             if (recvcnts(q) > 0) then
564                count = count + 1
565                if (count > preposts) then
566                   tail = mod(head,preposts) + 1
567                   call mpi_wait (rcvid(tail), status, ier)
568                end if
569                head = mod(head,preposts) + 1
570                call mpi_irecv ( recvbuf(displs(q)+1), recvcnts(q), &
571                                 recvtype, p, mtag, comm, rcvid(head), &
572                                 ier )
573                call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
574             end if
575          end if
576       end do
577
578! copy local data
579       q = mytid+1
580       do i=1,sendcnt
581          recvbuf(displs(q)+i) = sendbuf(i)
582       enddo
583
584! wait for final data
585       do i=1,min(count,preposts)
586          call mpi_wait (rcvid(i), status, ier)
587       enddo
588
589    else
590
591       if (sendcnt > 0) then
592          call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
593                          status, ier )
594          call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
595                           comm, ier )
596       end if
597
598    endif
599    if (ier /= 0) then
600       call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
601    end if
602
603 else
604 
605    call mpi_gatherv (sendbuf, sendcnt, sendtype, &
606                      recvbuf, recvcnts, displs, recvtype, &
607                      root, comm, ier)
608    if (ier /= 0) then
609       call MP_perr_die(myname_,':: MPI_GATHERV',ier)
610    end if
611
612 endif
613
614 return
615 end subroutine fc_gatherv_fp
616
617!BOP -------------------------------------------------------------------
618!
619! !IROUTINE: get_fcblocksize - return max_gather_block_size
620!
621! !DESCRIPTION:
622! This function returns the current value of max_gather_block_size
623!
624! !INTERFACE:
625
626 function get_fcblocksize()
627
628! !USES:
629!
630! No external modules are used by this function.
631
632     implicit none
633
634! !INPUT PARAMETERS:
635!
636
637! !OUTPUT PARAMETERS:
638!
639    integer           :: get_fcblocksize
640
641! !REVISION HISTORY:
642!       03Mar09 - R. Jacob (jacob@mcs.anl.gov) -- intial version
643!EOP ___________________________________________________________________
644
645  character(len=*),parameter :: myname_=myname//'::get_fcblocksize'
646
647  get_fcblocksize = max_gather_block_size
648
649 end function  get_fcblocksize
650
651!BOP -------------------------------------------------------------------
652!
653! !IROUTINE: set_fcblocksize - set max_gather_block_size
654!
655! !DESCRIPTION:
656! This function sets the current value of max_gather_block_size
657!
658! !INTERFACE:
659
660 subroutine set_fcblocksize(gather_block_size)
661
662! !USES:
663!
664! No external modules are used by this function.
665
666     implicit none
667
668! !INPUT PARAMETERS:
669!
670    integer           :: gather_block_size
671
672! !OUTPUT PARAMETERS:
673!
674
675! !REVISION HISTORY:
676!       03Mar09 - R. Jacob (jacob@mcs.anl.gov) -- intial version
677!EOP ___________________________________________________________________
678
679  character(len=*),parameter :: myname_=myname//':: set_fcblocksize'
680
681  max_gather_block_size = gather_block_size
682
683 end subroutine  set_fcblocksize
684
685 end module m_FcComms
Note: See TracBrowser for help on using the repository browser.