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

Last change on this file since 6331 was 6331, checked in by aclsce, 15 months ago

Moved oasis-mct_5.0 in oasis3-mct/branches directory.

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