123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685 |
- !BOP -------------------------------------------------------------------
- !
- ! !MODULE: m_FcComms - MPI collective communication operators
- ! with explict flow control
- !
- ! !DESCRIPTION:
- !
- ! This module includes implementations of MPI collective operators that
- ! have proven problematic on certain systems when run at scale. By
- ! introducing additonal flow control, these problems (exhausting internal
- ! system resources) can be avoided. These routines were ported from
- ! the Community Atmosphere Model's spmd_utils.F90.
- !
- ! !INTERFACE:
- !
- ! Workaround for performance issue with rsend on cray systems with
- ! gemini interconnect
- !
- #ifdef _NO_MPI_RSEND
- #define MPI_RSEND MPI_SEND
- #define mpi_rsend mpi_send
- #define MPI_IRSEND MPI_ISEND
- #define mpi_irsend mpi_isend
- #endif
- module m_FcComms
- implicit none
- private ! except
- public :: fc_gather_int ! flow control version of mpi_gather for integer vectors
- public :: fc_gather_fp ! flow control version of mpi_gather for FP vectors
- public :: fc_gatherv_int ! flow control version of mpi_gatherv for integer vectors
- public :: fc_gatherv_fp ! flow control version of mpi_gatherv for integer vectors
- public :: get_fcblocksize ! get current value of max_gather_block_size
- public :: set_fcblocksize ! set current value of max_gather_block_size
- ! !REVISION HISTORY:
- ! 30Jan09 - P.H. Worley <worleyph@ornl.gov> - imported routines
- ! from CAM's spmd_utils to create this module.
- integer, public :: max_gather_block_size = 64
- character(len=*),parameter :: myname='MCT(MPEU)::m_FcComms'
- contains
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: fc_gather_int - Gather an array of type integer
- !
- ! !DESCRIPTION:
- ! This routine gathers a {\em distributed} array of type {\em integer}
- ! to the {\tt root} process. Explicit handshaking messages are used
- ! to control the number of processes communicating with the root
- ! at any one time.
- !
- ! If flow_cntl optional parameter
- ! < 0 : use MPI_Gather
- ! >= 0: use point-to-point with handshaking messages and
- ! preposting receive requests up to
- ! min(max(1,flow_cntl),max_gather_block_size)
- ! ahead if optional flow_cntl parameter is present.
- ! Otherwise, max_gather_block_size is used in its place.
- ! Default value is max_gather_block_size.
- ! !INTERFACE:
- !
- subroutine fc_gather_int (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnt, recvtype, &
- root, comm, flow_cntl )
- !
- ! !USES:
- !
- use m_die
- use m_mpif90
- !
- ! !INPUT PARAMETERS:
- !
- integer, intent(in) :: sendbuf(*)
- integer, intent(in) :: sendcnt
- integer, intent(in) :: sendtype
- integer, intent(in) :: recvcnt
- integer, intent(in) :: recvtype
- integer, intent(in) :: root
- integer, intent(in) :: comm
- integer, optional, intent(in) :: flow_cntl
- ! !OUTPUT PARAMETERS:
- !
- integer, intent(out) :: recvbuf(*)
- ! !REVISION HISTORY:
- ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::fc_gather_int'
- integer :: signal
- logical fc_gather ! use explicit flow control?
- integer gather_block_size ! number of preposted receive requests
- integer :: mytid, mysize, mtag, p, i, count, displs
- integer :: preposts, head, tail
- integer :: rcvid(max_gather_block_size)
- integer :: status(MP_STATUS_SIZE)
- integer :: ier ! MPI error code
- signal = 1
- if ( present(flow_cntl) ) then
- if (flow_cntl >= 0) then
- gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
- fc_gather = .true.
- else
- fc_gather = .false.
- endif
- else
- gather_block_size = max(1,max_gather_block_size)
- fc_gather = .true.
- endif
- if (fc_gather) then
-
- call mpi_comm_rank (comm, mytid, ier)
- call mpi_comm_size (comm, mysize, ier)
- mtag = 0
- if (root .eq. mytid) then
- ! prepost gather_block_size irecvs, and start receiving data
- preposts = min(mysize-1, gather_block_size)
- head = 0
- count = 0
- do p=0, mysize-1
- if (p .ne. root) then
- if (recvcnt > 0) then
- count = count + 1
- if (count > preposts) then
- tail = mod(head,preposts) + 1
- call mpi_wait (rcvid(tail), status, ier)
- end if
- head = mod(head,preposts) + 1
- displs = p*recvcnt
- call mpi_irecv ( recvbuf(displs+1), recvcnt, &
- recvtype, p, mtag, comm, rcvid(head), &
- ier )
- call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
- end if
- end if
- end do
- ! copy local data
- displs = mytid*recvcnt
- do i=1,sendcnt
- recvbuf(displs+i) = sendbuf(i)
- enddo
- ! wait for final data
- do i=1,min(count,preposts)
- call mpi_wait (rcvid(i), status, ier)
- enddo
- else
- if (sendcnt > 0) then
- call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
- status, ier )
- call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
- comm, ier )
- end if
- endif
- if (ier /= 0) then
- call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
- end if
- else
-
- call mpi_gather (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnt, recvtype, &
- root, comm, ier)
- if (ier /= 0) then
- call MP_perr_die(myname_,':: MPI_GATHER',ier)
- end if
- endif
- return
- end subroutine fc_gather_int
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: fc_gather_fp - Gather an array of type FP
- !
- ! !DESCRIPTION:
- ! This routine gathers a {\em distributed} array of type {\em FP} to
- ! the {\tt root} process. Explicit handshaking messages are used
- ! to control the number of processes communicating with the root
- ! at any one time.
- !
- ! If flow_cntl optional parameter
- ! < 0 : use MPI_Gather
- ! >= 0: use point-to-point with handshaking messages and
- ! preposting receive requests up to
- ! min(max(1,flow_cntl),max_gather_block_size)
- ! ahead if optional flow_cntl parameter is present.
- ! Otherwise, max_gather_block_size is used in its place.
- ! Default value is max_gather_block_size.
- ! !INTERFACE:
- !
- subroutine fc_gather_fp (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnt, recvtype, &
- root, comm, flow_cntl )
- !
- ! !USES:
- !
- use m_realkinds, only : FP
- use m_die
- use m_mpif90
- !
- ! !INPUT PARAMETERS:
- !
- real (FP), intent(in) :: sendbuf(*)
- integer, intent(in) :: sendcnt
- integer, intent(in) :: sendtype
- integer, intent(in) :: recvcnt
- integer, intent(in) :: recvtype
- integer, intent(in) :: root
- integer, intent(in) :: comm
- integer, optional, intent(in) :: flow_cntl
- ! !OUTPUT PARAMETERS:
- !
- real (FP), intent(out) :: recvbuf(*)
- ! !REVISION HISTORY:
- ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::fc_gather_fp'
- real (FP) :: signal
- logical fc_gather ! use explicit flow control?
- integer gather_block_size ! number of preposted receive requests
- integer :: mytid, mysize, mtag, p, i, count, displs
- integer :: preposts, head, tail
- integer :: rcvid(max_gather_block_size)
- integer :: status(MP_STATUS_SIZE)
- integer :: ier ! MPI error code
- signal = 1.0
- if ( present(flow_cntl) ) then
- if (flow_cntl >= 0) then
- gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
- fc_gather = .true.
- else
- fc_gather = .false.
- endif
- else
- gather_block_size = max(1,max_gather_block_size)
- fc_gather = .true.
- endif
- if (fc_gather) then
-
- call mpi_comm_rank (comm, mytid, ier)
- call mpi_comm_size (comm, mysize, ier)
- mtag = 0
- if (root .eq. mytid) then
- ! prepost gather_block_size irecvs, and start receiving data
- preposts = min(mysize-1, gather_block_size)
- head = 0
- count = 0
- do p=0, mysize-1
- if (p .ne. root) then
- if (recvcnt > 0) then
- count = count + 1
- if (count > preposts) then
- tail = mod(head,preposts) + 1
- call mpi_wait (rcvid(tail), status, ier)
- end if
- head = mod(head,preposts) + 1
- displs = p*recvcnt
- call mpi_irecv ( recvbuf(displs+1), recvcnt, &
- recvtype, p, mtag, comm, rcvid(head), &
- ier )
- call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
- end if
- end if
- end do
- ! copy local data
- displs = mytid*recvcnt
- do i=1,sendcnt
- recvbuf(displs+i) = sendbuf(i)
- enddo
- ! wait for final data
- do i=1,min(count,preposts)
- call mpi_wait (rcvid(i), status, ier)
- enddo
- else
- if (sendcnt > 0) then
- call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
- status, ier )
- call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
- comm, ier )
- end if
- endif
- if (ier /= 0) then
- call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
- end if
- else
-
- call mpi_gather (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnt, recvtype, &
- root, comm, ier)
- if (ier /= 0) then
- call MP_perr_die(myname_,':: MPI_GATHER',ier)
- end if
- endif
- return
- end subroutine fc_gather_fp
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: fc_gatherv_int - Gather an array of type integer
- !
- ! !DESCRIPTION:
- ! This routine gathers a {\em distributed} array of type {\em integer}
- ! to the {\tt root} process. Explicit handshaking messages are used
- ! to control the number of processes communicating with the root
- ! at any one time.
- !
- ! If flow_cntl optional parameter
- ! < 0 : use MPI_Gatherv
- ! >= 0: use point-to-point with handshaking messages and
- ! preposting receive requests up to
- ! min(max(1,flow_cntl),max_gather_block_size)
- ! ahead if optional flow_cntl parameter is present.
- ! Otherwise, max_gather_block_size is used in its place.
- ! Default value is max_gather_block_size.
- ! !INTERFACE:
- !
- subroutine fc_gatherv_int (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnts, displs, recvtype, &
- root, comm, flow_cntl )
- !
- ! !USES:
- !
- use m_die
- use m_mpif90
- !
- ! !INPUT PARAMETERS:
- !
- integer, intent(in) :: sendbuf(*)
- integer, intent(in) :: sendcnt
- integer, intent(in) :: sendtype
- integer, intent(in) :: recvcnts(*)
- integer, intent(in) :: displs(*)
- integer, intent(in) :: recvtype
- integer, intent(in) :: root
- integer, intent(in) :: comm
- integer, optional, intent(in) :: flow_cntl
- ! !OUTPUT PARAMETERS:
- !
- integer, intent(out) :: recvbuf(*)
- ! !REVISION HISTORY:
- ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::fc_gatherv_int'
- integer :: signal
- logical fc_gather ! use explicit flow control?
- integer gather_block_size ! number of preposted receive requests
- integer :: mytid, mysize, mtag, p, q, i, count
- integer :: preposts, head, tail
- integer :: rcvid(max_gather_block_size)
- integer :: status(MP_STATUS_SIZE)
- integer :: ier ! MPI error code
- signal = 1
- if ( present(flow_cntl) ) then
- if (flow_cntl >= 0) then
- gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
- fc_gather = .true.
- else
- fc_gather = .false.
- endif
- else
- gather_block_size = max(1,max_gather_block_size)
- fc_gather = .true.
- endif
- if (fc_gather) then
-
- call mpi_comm_rank (comm, mytid, ier)
- call mpi_comm_size (comm, mysize, ier)
- mtag = 0
- if (root .eq. mytid) then
- ! prepost gather_block_size irecvs, and start receiving data
- preposts = min(mysize-1, gather_block_size)
- head = 0
- count = 0
- do p=0, mysize-1
- if (p .ne. root) then
- q = p+1
- if (recvcnts(q) > 0) then
- count = count + 1
- if (count > preposts) then
- tail = mod(head,preposts) + 1
- call mpi_wait (rcvid(tail), status, ier)
- end if
- head = mod(head,preposts) + 1
- call mpi_irecv ( recvbuf(displs(q)+1), recvcnts(q), &
- recvtype, p, mtag, comm, rcvid(head), &
- ier )
- call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
- end if
- end if
- end do
- ! copy local data
- q = mytid+1
- do i=1,sendcnt
- recvbuf(displs(q)+i) = sendbuf(i)
- enddo
- ! wait for final data
- do i=1,min(count,preposts)
- call mpi_wait (rcvid(i), status, ier)
- enddo
- else
- if (sendcnt > 0) then
- call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
- status, ier )
- call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
- comm, ier )
- end if
- endif
- if (ier /= 0) then
- call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
- end if
- else
-
- call mpi_gatherv (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnts, displs, recvtype, &
- root, comm, ier)
- if (ier /= 0) then
- call MP_perr_die(myname_,':: MPI_GATHERV',ier)
- end if
- endif
- return
- end subroutine fc_gatherv_int
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: fc_gatherv_fp - Gather an array of type FP
- !
- ! !DESCRIPTION:
- ! This routine gathers a {\em distributed} array of type {\em FP} to
- ! the {\tt root} process. Explicit handshaking messages are used
- ! to control the number of processes communicating with the root
- ! at any one time.
- !
- ! If flow_cntl optional parameter
- ! < 0 : use MPI_Gatherv
- ! >= 0: use point-to-point with handshaking messages and
- ! preposting receive requests up to
- ! min(max(1,flow_cntl),max_gather_block_size)
- ! ahead if optional flow_cntl parameter is present.
- ! Otherwise, max_gather_block_size is used in its place.
- ! Default value is max_gather_block_size.
- ! !INTERFACE:
- !
- subroutine fc_gatherv_fp (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnts, displs, recvtype, &
- root, comm, flow_cntl )
- !
- ! !USES:
- !
- use m_realkinds, only : FP
- use m_die
- use m_mpif90
- !
- ! !INPUT PARAMETERS:
- !
- real (FP), intent(in) :: sendbuf(*)
- integer, intent(in) :: sendcnt
- integer, intent(in) :: sendtype
- integer, intent(in) :: recvcnts(*)
- integer, intent(in) :: displs(*)
- integer, intent(in) :: recvtype
- integer, intent(in) :: root
- integer, intent(in) :: comm
- integer, optional, intent(in) :: flow_cntl
- ! !OUTPUT PARAMETERS:
- !
- real (FP), intent(out) :: recvbuf(*)
- ! !REVISION HISTORY:
- ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::fc_gatherv_fp'
- real (FP) :: signal
- logical fc_gather ! use explicit flow control?
- integer gather_block_size ! number of preposted receive requests
- integer :: mytid, mysize, mtag, p, q, i, count
- integer :: preposts, head, tail
- integer :: rcvid(max_gather_block_size)
- integer :: status(MP_STATUS_SIZE)
- integer :: ier ! MPI error code
- signal = 1.0
- if ( present(flow_cntl) ) then
- if (flow_cntl >= 0) then
- gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
- fc_gather = .true.
- else
- fc_gather = .false.
- endif
- else
- gather_block_size = max(1,max_gather_block_size)
- fc_gather = .true.
- endif
- if (fc_gather) then
-
- call mpi_comm_rank (comm, mytid, ier)
- call mpi_comm_size (comm, mysize, ier)
- mtag = 0
- if (root .eq. mytid) then
- ! prepost gather_block_size irecvs, and start receiving data
- preposts = min(mysize-1, gather_block_size)
- head = 0
- count = 0
- do p=0, mysize-1
- if (p .ne. root) then
- q = p+1
- if (recvcnts(q) > 0) then
- count = count + 1
- if (count > preposts) then
- tail = mod(head,preposts) + 1
- call mpi_wait (rcvid(tail), status, ier)
- end if
- head = mod(head,preposts) + 1
- call mpi_irecv ( recvbuf(displs(q)+1), recvcnts(q), &
- recvtype, p, mtag, comm, rcvid(head), &
- ier )
- call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
- end if
- end if
- end do
- ! copy local data
- q = mytid+1
- do i=1,sendcnt
- recvbuf(displs(q)+i) = sendbuf(i)
- enddo
- ! wait for final data
- do i=1,min(count,preposts)
- call mpi_wait (rcvid(i), status, ier)
- enddo
- else
- if (sendcnt > 0) then
- call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
- status, ier )
- call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
- comm, ier )
- end if
- endif
- if (ier /= 0) then
- call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
- end if
- else
-
- call mpi_gatherv (sendbuf, sendcnt, sendtype, &
- recvbuf, recvcnts, displs, recvtype, &
- root, comm, ier)
- if (ier /= 0) then
- call MP_perr_die(myname_,':: MPI_GATHERV',ier)
- end if
- endif
- return
- end subroutine fc_gatherv_fp
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: get_fcblocksize - return max_gather_block_size
- !
- ! !DESCRIPTION:
- ! This function returns the current value of max_gather_block_size
- !
- ! !INTERFACE:
- function get_fcblocksize()
- ! !USES:
- !
- ! No external modules are used by this function.
- implicit none
- ! !INPUT PARAMETERS:
- !
- ! !OUTPUT PARAMETERS:
- !
- integer :: get_fcblocksize
- ! !REVISION HISTORY:
- ! 03Mar09 - R. Jacob (jacob@mcs.anl.gov) -- intial version
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::get_fcblocksize'
- get_fcblocksize = max_gather_block_size
- end function get_fcblocksize
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: set_fcblocksize - set max_gather_block_size
- !
- ! !DESCRIPTION:
- ! This function sets the current value of max_gather_block_size
- !
- ! !INTERFACE:
- subroutine set_fcblocksize(gather_block_size)
- ! !USES:
- !
- ! No external modules are used by this function.
- implicit none
- ! !INPUT PARAMETERS:
- !
- integer :: gather_block_size
- ! !OUTPUT PARAMETERS:
- !
- ! !REVISION HISTORY:
- ! 03Mar09 - R. Jacob (jacob@mcs.anl.gov) -- intial version
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//':: set_fcblocksize'
- max_gather_block_size = gather_block_size
- end subroutine set_fcblocksize
- end module m_FcComms
|