xref: /petsc/src/sys/utils/mpits.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1f6ced4a3SJed Brown #include <petscsys.h>        /*I  "petscsys.h"  I*/
295c0884eSLisandro Dalcin #include <petsc/private/petscimpl.h>
3f6ced4a3SJed Brown 
495c0884eSLisandro Dalcin PetscLogEvent PETSC_BuildTwoSided;
595c0884eSLisandro Dalcin PetscLogEvent PETSC_BuildTwoSidedF;
63b3561c8SJed Brown 
76145cd65SJed Brown const char *const PetscBuildTwoSidedTypes[] = {
8f6ced4a3SJed Brown   "ALLREDUCE",
96145cd65SJed Brown   "IBARRIER",
101654bf6bSJed Brown   "REDSCATTER",
116145cd65SJed Brown   "PetscBuildTwoSidedType",
126145cd65SJed Brown   "PETSC_BUILDTWOSIDED_",
1302c9f0b5SLisandro Dalcin   NULL
14f6ced4a3SJed Brown };
15f6ced4a3SJed Brown 
166145cd65SJed Brown static PetscBuildTwoSidedType _twosided_type = PETSC_BUILDTWOSIDED_NOTSET;
17f6ced4a3SJed Brown 
186145cd65SJed Brown /*@
196145cd65SJed Brown    PetscCommBuildTwoSidedSetType - set algorithm to use when building two-sided communication
206145cd65SJed Brown 
216145cd65SJed Brown    Logically Collective
226145cd65SJed Brown 
234165533cSJose E. Roman    Input Parameters:
246145cd65SJed Brown +  comm - PETSC_COMM_WORLD
256145cd65SJed Brown -  twosided - algorithm to use in subsequent calls to PetscCommBuildTwoSided()
266145cd65SJed Brown 
276145cd65SJed Brown    Level: developer
286145cd65SJed Brown 
296145cd65SJed Brown    Note:
306145cd65SJed Brown    This option is currently global, but could be made per-communicator.
316145cd65SJed Brown 
326145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedGetType()
336145cd65SJed Brown @*/
346145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm comm,PetscBuildTwoSidedType twosided)
356145cd65SJed Brown {
366145cd65SJed Brown   PetscFunctionBegin;
3776bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {                             /* We don't have a PetscObject so can't use PetscValidLogicalCollectiveEnum */
386145cd65SJed Brown     PetscMPIInt b1[2],b2[2];
396145cd65SJed Brown     b1[0] = -(PetscMPIInt)twosided;
406145cd65SJed Brown     b1[1] = (PetscMPIInt)twosided;
415f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,comm));
422c71b3e2SJacob Faibussowitsch     PetscCheckFalse(-b2[0] != b2[1],comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes");
436145cd65SJed Brown   }
446145cd65SJed Brown   _twosided_type = twosided;
456145cd65SJed Brown   PetscFunctionReturn(0);
466145cd65SJed Brown }
476145cd65SJed Brown 
486145cd65SJed Brown /*@
496145cd65SJed Brown    PetscCommBuildTwoSidedGetType - set algorithm to use when building two-sided communication
506145cd65SJed Brown 
516145cd65SJed Brown    Logically Collective
526145cd65SJed Brown 
534165533cSJose E. Roman    Output Parameters:
546145cd65SJed Brown +  comm - communicator on which to query algorithm
556145cd65SJed Brown -  twosided - algorithm to use for PetscCommBuildTwoSided()
566145cd65SJed Brown 
576145cd65SJed Brown    Level: developer
586145cd65SJed Brown 
596145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType()
606145cd65SJed Brown @*/
616145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm comm,PetscBuildTwoSidedType *twosided)
62f6ced4a3SJed Brown {
6305342407SJunchao Zhang   PetscMPIInt    size;
64f6ced4a3SJed Brown 
65f6ced4a3SJed Brown   PetscFunctionBegin;
666145cd65SJed Brown   *twosided = PETSC_BUILDTWOSIDED_NOTSET;
676145cd65SJed Brown   if (_twosided_type == PETSC_BUILDTWOSIDED_NOTSET) {
685f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_size(comm,&size));
6905342407SJunchao Zhang     _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; /* default for small comms, see https://gitlab.com/petsc/petsc/-/merge_requests/2611 */
708ee3e7ecSJunchao Zhang #if defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
7105342407SJunchao Zhang     if (size > 1024) _twosided_type = PETSC_BUILDTWOSIDED_IBARRIER;
726145cd65SJed Brown #endif
735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetEnum(NULL,NULL,"-build_twosided",PetscBuildTwoSidedTypes,(PetscEnum*)&_twosided_type,NULL));
74f6ced4a3SJed Brown   }
75f6ced4a3SJed Brown   *twosided = _twosided_type;
76f6ced4a3SJed Brown   PetscFunctionReturn(0);
77f6ced4a3SJed Brown }
78f6ced4a3SJed Brown 
798ee3e7ecSJunchao Zhang #if defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
80cf4b5b4fSJed Brown static PetscErrorCode PetscCommBuildTwoSided_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
81f6ced4a3SJed Brown {
82c3a59b84SJed Brown   PetscMPIInt    nrecvs,tag,done,i;
83c3a59b84SJed Brown   MPI_Aint       lb,unitbytes;
84f6ced4a3SJed Brown   char           *tdata;
85f6ced4a3SJed Brown   MPI_Request    *sendreqs,barrier;
860f453b92SJed Brown   PetscSegBuffer segrank,segdata;
8746e50bb6SJed Brown   PetscBool      barrier_started;
88f6ced4a3SJed Brown 
89f6ced4a3SJed Brown   PetscFunctionBegin;
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDuplicate(comm,&comm,&tag));
915f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_get_extent(dtype,&lb,&unitbytes));
922c71b3e2SJacob Faibussowitsch   PetscCheckFalse(lb != 0,comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld",(long)lb);
93f6ced4a3SJed Brown   tdata = (char*)todata;
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nto,&sendreqs));
95f6ced4a3SJed Brown   for (i=0; i<nto; i++) {
965f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i));
97f6ced4a3SJed Brown   }
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferCreate(unitbytes,4*count,&segdata));
100f6ced4a3SJed Brown 
101f6ced4a3SJed Brown   nrecvs  = 0;
102f6ced4a3SJed Brown   barrier = MPI_REQUEST_NULL;
10346e50bb6SJed Brown   /* MPICH-3.2 sometimes does not create a request in some "optimized" cases.  This is arguably a standard violation,
10446e50bb6SJed Brown    * but we need to work around it. */
10546e50bb6SJed Brown   barrier_started = PETSC_FALSE;
106f6ced4a3SJed Brown   for (done=0; !done;) {
107f6ced4a3SJed Brown     PetscMPIInt flag;
108f6ced4a3SJed Brown     MPI_Status  status;
1095f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status));
110f6ced4a3SJed Brown     if (flag) {                 /* incoming message */
111f6ced4a3SJed Brown       PetscMPIInt *recvrank;
112f6ced4a3SJed Brown       void        *buf;
1135f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSegBufferGet(segrank,1,&recvrank));
1145f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSegBufferGet(segdata,count,&buf));
115f6ced4a3SJed Brown       *recvrank = status.MPI_SOURCE;
1165f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE));
117f6ced4a3SJed Brown       nrecvs++;
118f6ced4a3SJed Brown     }
11946e50bb6SJed Brown     if (!barrier_started) {
1204dc2109aSBarry Smith       PetscMPIInt sent,nsends;
1215f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMPIIntCast(nto,&nsends));
1225f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE));
123f6ced4a3SJed Brown       if (sent) {
1245f80ce2aSJacob Faibussowitsch         CHKERRMPI(MPI_Ibarrier(comm,&barrier));
12546e50bb6SJed Brown         barrier_started = PETSC_TRUE;
1265f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(sendreqs));
127f6ced4a3SJed Brown       }
128f6ced4a3SJed Brown     } else {
1295f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Test(&barrier,&done,MPI_STATUS_IGNORE));
130f6ced4a3SJed Brown     }
131f6ced4a3SJed Brown   }
132f6ced4a3SJed Brown   *nfrom = nrecvs;
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferExtractAlloc(segrank,fromranks));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferDestroy(&segrank));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferExtractAlloc(segdata,fromdata));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferDestroy(&segdata));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDestroy(&comm));
138f6ced4a3SJed Brown   PetscFunctionReturn(0);
139f6ced4a3SJed Brown }
140f6ced4a3SJed Brown #endif
141f6ced4a3SJed Brown 
142cf4b5b4fSJed Brown static PetscErrorCode PetscCommBuildTwoSided_Allreduce(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
143f6ced4a3SJed Brown {
14405342407SJunchao Zhang   PetscMPIInt      size,rank,*iflags,nrecvs,tag,*franks,i,flg;
145c3a59b84SJed Brown   MPI_Aint         lb,unitbytes;
146f6ced4a3SJed Brown   char             *tdata,*fdata;
147f6ced4a3SJed Brown   MPI_Request      *reqs,*sendreqs;
148f6ced4a3SJed Brown   MPI_Status       *statuses;
14905342407SJunchao Zhang   PetscCommCounter *counter;
150f6ced4a3SJed Brown 
151f6ced4a3SJed Brown   PetscFunctionBegin;
1525f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
1535f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDuplicate(comm,&comm,&tag));
1555f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_get_attr(comm,Petsc_Counter_keyval,&counter,&flg));
156*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner PETSc communicator does not have its tag/name counter attribute set");
15705342407SJunchao Zhang   if (!counter->iflags) {
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(size,&counter->iflags));
15905342407SJunchao Zhang     iflags = counter->iflags;
16005342407SJunchao Zhang   } else {
16105342407SJunchao Zhang     iflags = counter->iflags;
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(iflags,size));
16305342407SJunchao Zhang   }
16405342407SJunchao Zhang   for (i=0; i<nto; i++) iflags[toranks[i]] = 1;
1655f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,iflags,size,MPI_INT,MPI_SUM,comm));
16605342407SJunchao Zhang   nrecvs   = iflags[rank];
1675f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_get_extent(dtype,&lb,&unitbytes));
1685f80ce2aSJacob Faibussowitsch   PetscCheck(lb == 0,comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld",(long)lb);
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc(nrecvs*count*unitbytes,&fdata));
170f6ced4a3SJed Brown   tdata    = (char*)todata;
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses));
17206926cafSJed Brown   sendreqs = reqs + nrecvs;
173f6ced4a3SJed Brown   for (i=0; i<nrecvs; i++) {
1745f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i));
175f6ced4a3SJed Brown   }
176f6ced4a3SJed Brown   for (i=0; i<nto; i++) {
1775f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i));
178f6ced4a3SJed Brown   }
1795f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(nto+nrecvs,reqs,statuses));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,&franks));
181a297a907SKarl Rupp   for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE;
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(reqs,statuses));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDestroy(&comm));
184f6ced4a3SJed Brown 
185f6ced4a3SJed Brown   *nfrom            = nrecvs;
186f6ced4a3SJed Brown   *fromranks        = franks;
187f6ced4a3SJed Brown   *(void**)fromdata = fdata;
188f6ced4a3SJed Brown   PetscFunctionReturn(0);
189f6ced4a3SJed Brown }
190f6ced4a3SJed Brown 
1911654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
192cf4b5b4fSJed Brown static PetscErrorCode PetscCommBuildTwoSided_RedScatter(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
1931654bf6bSJed Brown {
19405342407SJunchao Zhang   PetscMPIInt    size,*iflags,nrecvs,tag,*franks,i,flg;
195c3a59b84SJed Brown   MPI_Aint       lb,unitbytes;
1961654bf6bSJed Brown   char           *tdata,*fdata;
1971654bf6bSJed Brown   MPI_Request    *reqs,*sendreqs;
1981654bf6bSJed Brown   MPI_Status     *statuses;
19905342407SJunchao Zhang   PetscCommCounter *counter;
2001654bf6bSJed Brown 
2011654bf6bSJed Brown   PetscFunctionBegin;
2025f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDuplicate(comm,&comm,&tag));
2045f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_get_attr(comm,Petsc_Counter_keyval,&counter,&flg));
205*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner PETSc communicator does not have its tag/name counter attribute set");
20605342407SJunchao Zhang   if (!counter->iflags) {
2075f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(size,&counter->iflags));
20805342407SJunchao Zhang     iflags = counter->iflags;
20905342407SJunchao Zhang   } else {
21005342407SJunchao Zhang     iflags = counter->iflags;
2115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(iflags,size));
21205342407SJunchao Zhang   }
2131654bf6bSJed Brown   for (i=0; i<nto; i++) iflags[toranks[i]] = 1;
2145f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Reduce_scatter_block(iflags,&nrecvs,1,MPI_INT,MPI_SUM,comm));
2155f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_get_extent(dtype,&lb,&unitbytes));
2162c71b3e2SJacob Faibussowitsch   PetscCheckFalse(lb != 0,comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld",(long)lb);
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc(nrecvs*count*unitbytes,&fdata));
2181654bf6bSJed Brown   tdata    = (char*)todata;
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses));
2201654bf6bSJed Brown   sendreqs = reqs + nrecvs;
2211654bf6bSJed Brown   for (i=0; i<nrecvs; i++) {
2225f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i));
2231654bf6bSJed Brown   }
2241654bf6bSJed Brown   for (i=0; i<nto; i++) {
2255f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i));
2261654bf6bSJed Brown   }
2275f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(nto+nrecvs,reqs,statuses));
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,&franks));
2291654bf6bSJed Brown   for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE;
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(reqs,statuses));
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDestroy(&comm));
2321654bf6bSJed Brown 
2331654bf6bSJed Brown   *nfrom            = nrecvs;
2341654bf6bSJed Brown   *fromranks        = franks;
2351654bf6bSJed Brown   *(void**)fromdata = fdata;
2361654bf6bSJed Brown   PetscFunctionReturn(0);
2371654bf6bSJed Brown }
2381654bf6bSJed Brown #endif
2391654bf6bSJed Brown 
240f6ced4a3SJed Brown /*@C
241f6ced4a3SJed Brown    PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths)
242f6ced4a3SJed Brown 
243d083f849SBarry Smith    Collective
244f6ced4a3SJed Brown 
2454165533cSJose E. Roman    Input Parameters:
246f6ced4a3SJed Brown +  comm - communicator
247f6ced4a3SJed Brown .  count - number of entries to send/receive (must match on all ranks)
248f6ced4a3SJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
249f6ced4a3SJed Brown .  nto - number of ranks to send data to
250f6ced4a3SJed Brown .  toranks - ranks to send to (array of length nto)
251f6ced4a3SJed Brown -  todata - data to send to each rank (packed)
252f6ced4a3SJed Brown 
2534165533cSJose E. Roman    Output Parameters:
254f6ced4a3SJed Brown +  nfrom - number of ranks receiving messages from
255f6ced4a3SJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
256f6ced4a3SJed Brown -  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
257f6ced4a3SJed Brown 
258f6ced4a3SJed Brown    Level: developer
259f6ced4a3SJed Brown 
2601654bf6bSJed Brown    Options Database Keys:
26105342407SJunchao Zhang .  -build_twosided <allreduce|ibarrier|redscatter> - algorithm to set up two-sided communication. Default is allreduce for communicators with <= 1024 ranks, otherwise ibarrier.
2621654bf6bSJed Brown 
263f6ced4a3SJed Brown    Notes:
264f6ced4a3SJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
265f6ced4a3SJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data.
266f6ced4a3SJed Brown 
267f6ced4a3SJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
268f6ced4a3SJed Brown 
269f6ced4a3SJed Brown    References:
270606c0280SSatish Balay .  * - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in
2714f02bc6aSBarry Smith    Scalable communication protocols for dynamic sparse data exchange, 2010.
272f6ced4a3SJed Brown 
273f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
274f6ced4a3SJed Brown @*/
275cf4b5b4fSJed Brown PetscErrorCode PetscCommBuildTwoSided(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
276f6ced4a3SJed Brown {
2776145cd65SJed Brown   PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET;
278f6ced4a3SJed Brown 
279f6ced4a3SJed Brown   PetscFunctionBegin;
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSysInitializePackage());
2815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventSync(PETSC_BuildTwoSided,comm));
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(PETSC_BuildTwoSided,0,0,0,0));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommBuildTwoSidedGetType(comm,&buildtype));
284f6ced4a3SJed Brown   switch (buildtype) {
2856145cd65SJed Brown   case PETSC_BUILDTWOSIDED_IBARRIER:
2868ee3e7ecSJunchao Zhang #if defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
2875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata));
288b458e8f1SJose E. Roman     break;
2896145cd65SJed Brown #else
2906145cd65SJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)");
291f6ced4a3SJed Brown #endif
2926145cd65SJed Brown   case PETSC_BUILDTWOSIDED_ALLREDUCE:
2935f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata));
294f6ced4a3SJed Brown     break;
2951654bf6bSJed Brown   case PETSC_BUILDTWOSIDED_REDSCATTER:
2961654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
2975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommBuildTwoSided_RedScatter(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata));
298b458e8f1SJose E. Roman     break;
2991654bf6bSJed Brown #else
3001654bf6bSJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Reduce_scatter_block (part of MPI-2.2)");
3011654bf6bSJed Brown #endif
302f6ced4a3SJed Brown   default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication");
303f6ced4a3SJed Brown   }
3045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(PETSC_BuildTwoSided,0,0,0,0));
305f6ced4a3SJed Brown   PetscFunctionReturn(0);
306f6ced4a3SJed Brown }
307d815da10SJed Brown 
3086afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Reference(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,
3096afbe7ddSJed Brown                                                            PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs,
310d815da10SJed Brown                                                            PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
311d815da10SJed Brown                                                            PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
312d815da10SJed Brown {
313d815da10SJed Brown   PetscMPIInt i,*tag;
314d815da10SJed Brown   MPI_Aint    lb,unitbytes;
315d815da10SJed Brown   MPI_Request *sendreq,*recvreq;
316d815da10SJed Brown 
317d815da10SJed Brown   PetscFunctionBegin;
3185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ntags,&tag));
319d815da10SJed Brown   if (ntags > 0) {
3205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommDuplicate(comm,&comm,&tag[0]));
321d815da10SJed Brown   }
322d815da10SJed Brown   for (i=1; i<ntags; i++) {
3235f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommGetNewTag(comm,&tag[i]));
324d815da10SJed Brown   }
325d815da10SJed Brown 
326d815da10SJed Brown   /* Perform complete initial rendezvous */
3275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommBuildTwoSided(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata));
328d815da10SJed Brown 
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nto*ntags,&sendreq));
3305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(*nfrom*ntags,&recvreq));
331d815da10SJed Brown 
3325f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_get_extent(dtype,&lb,&unitbytes));
3332c71b3e2SJacob Faibussowitsch   PetscCheckFalse(lb != 0,comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld",(long)lb);
334d815da10SJed Brown   for (i=0; i<nto; i++) {
33589157a57SJed Brown     PetscMPIInt k;
33689157a57SJed Brown     for (k=0; k<ntags; k++) sendreq[i*ntags+k] = MPI_REQUEST_NULL;
3375f80ce2aSJacob Faibussowitsch     CHKERRQ((*send)(comm,tag,i,toranks[i],((char*)todata)+count*unitbytes*i,sendreq+i*ntags,ctx));
338d815da10SJed Brown   }
339d815da10SJed Brown   for (i=0; i<*nfrom; i++) {
340d815da10SJed Brown     void *header = (*(char**)fromdata) + count*unitbytes*i;
34189157a57SJed Brown     PetscMPIInt k;
34289157a57SJed Brown     for (k=0; k<ntags; k++) recvreq[i*ntags+k] = MPI_REQUEST_NULL;
3435f80ce2aSJacob Faibussowitsch     CHKERRQ((*recv)(comm,tag,(*fromranks)[i],header,recvreq+i*ntags,ctx));
344d815da10SJed Brown   }
3455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tag));
3465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDestroy(&comm));
3476afbe7ddSJed Brown   *toreqs = sendreq;
3486afbe7ddSJed Brown   *fromreqs = recvreq;
349d815da10SJed Brown   PetscFunctionReturn(0);
350d815da10SJed Brown }
351d815da10SJed Brown 
3528ee3e7ecSJunchao Zhang #if defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
353afb254cdSJed Brown 
3546afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,
3556afbe7ddSJed Brown                                                           PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs,
356afb254cdSJed Brown                                                           PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
357afb254cdSJed Brown                                                           PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
358afb254cdSJed Brown {
359afb254cdSJed Brown   PetscMPIInt    nrecvs,tag,*tags,done,i;
360afb254cdSJed Brown   MPI_Aint       lb,unitbytes;
361afb254cdSJed Brown   char           *tdata;
362afb254cdSJed Brown   MPI_Request    *sendreqs,*usendreqs,*req,barrier;
363afb254cdSJed Brown   PetscSegBuffer segrank,segdata,segreq;
36446e50bb6SJed Brown   PetscBool      barrier_started;
365afb254cdSJed Brown 
366afb254cdSJed Brown   PetscFunctionBegin;
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDuplicate(comm,&comm,&tag));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ntags,&tags));
369afb254cdSJed Brown   for (i=0; i<ntags; i++) {
3705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCommGetNewTag(comm,&tags[i]));
371afb254cdSJed Brown   }
3725f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Type_get_extent(dtype,&lb,&unitbytes));
3732c71b3e2SJacob Faibussowitsch   PetscCheckFalse(lb != 0,comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld",(long)lb);
374afb254cdSJed Brown   tdata = (char*)todata;
3755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nto,&sendreqs));
3765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nto*ntags,&usendreqs));
377afb254cdSJed Brown   /* Post synchronous sends */
378afb254cdSJed Brown   for (i=0; i<nto; i++) {
3795f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i));
380afb254cdSJed Brown   }
381afb254cdSJed Brown   /* Post actual payloads.  These are typically larger messages.  Hopefully sending these later does not slow down the
382afb254cdSJed Brown    * synchronous messages above. */
383afb254cdSJed Brown   for (i=0; i<nto; i++) {
38489157a57SJed Brown     PetscMPIInt k;
38589157a57SJed Brown     for (k=0; k<ntags; k++) usendreqs[i*ntags+k] = MPI_REQUEST_NULL;
3865f80ce2aSJacob Faibussowitsch     CHKERRQ((*send)(comm,tags,i,toranks[i],tdata+count*unitbytes*i,usendreqs+i*ntags,ctx));
387afb254cdSJed Brown   }
388afb254cdSJed Brown 
3895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank));
3905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferCreate(unitbytes,4*count,&segdata));
3915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferCreate(sizeof(MPI_Request),4,&segreq));
392afb254cdSJed Brown 
393afb254cdSJed Brown   nrecvs  = 0;
394afb254cdSJed Brown   barrier = MPI_REQUEST_NULL;
39546e50bb6SJed Brown   /* MPICH-3.2 sometimes does not create a request in some "optimized" cases.  This is arguably a standard violation,
39646e50bb6SJed Brown    * but we need to work around it. */
39746e50bb6SJed Brown   barrier_started = PETSC_FALSE;
398afb254cdSJed Brown   for (done=0; !done;) {
399afb254cdSJed Brown     PetscMPIInt flag;
400afb254cdSJed Brown     MPI_Status  status;
4015f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status));
402afb254cdSJed Brown     if (flag) {                 /* incoming message */
40389157a57SJed Brown       PetscMPIInt *recvrank,k;
404afb254cdSJed Brown       void        *buf;
4055f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSegBufferGet(segrank,1,&recvrank));
4065f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSegBufferGet(segdata,count,&buf));
407afb254cdSJed Brown       *recvrank = status.MPI_SOURCE;
4085f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE));
4095f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSegBufferGet(segreq,ntags,&req));
41089157a57SJed Brown       for (k=0; k<ntags; k++) req[k] = MPI_REQUEST_NULL;
4115f80ce2aSJacob Faibussowitsch       CHKERRQ((*recv)(comm,tags,status.MPI_SOURCE,buf,req,ctx));
412afb254cdSJed Brown       nrecvs++;
413afb254cdSJed Brown     }
41446e50bb6SJed Brown     if (!barrier_started) {
415afb254cdSJed Brown       PetscMPIInt sent,nsends;
4165f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMPIIntCast(nto,&nsends));
4175f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE));
418afb254cdSJed Brown       if (sent) {
4195f80ce2aSJacob Faibussowitsch         CHKERRMPI(MPI_Ibarrier(comm,&barrier));
42046e50bb6SJed Brown         barrier_started = PETSC_TRUE;
421afb254cdSJed Brown       }
422afb254cdSJed Brown     } else {
4235f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Test(&barrier,&done,MPI_STATUS_IGNORE));
424afb254cdSJed Brown     }
425afb254cdSJed Brown   }
426afb254cdSJed Brown   *nfrom = nrecvs;
4275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferExtractAlloc(segrank,fromranks));
4285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferDestroy(&segrank));
4295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferExtractAlloc(segdata,fromdata));
4305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferDestroy(&segdata));
4316afbe7ddSJed Brown   *toreqs = usendreqs;
4325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferExtractAlloc(segreq,fromreqs));
4335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSegBufferDestroy(&segreq));
4345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sendreqs));
4355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tags));
4365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommDestroy(&comm));
437afb254cdSJed Brown   PetscFunctionReturn(0);
438afb254cdSJed Brown }
439afb254cdSJed Brown #endif
440afb254cdSJed Brown 
441d815da10SJed Brown /*@C
442d815da10SJed Brown    PetscCommBuildTwoSidedF - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous
443d815da10SJed Brown 
444d083f849SBarry Smith    Collective
445d815da10SJed Brown 
4464165533cSJose E. Roman    Input Parameters:
447d815da10SJed Brown +  comm - communicator
448d815da10SJed Brown .  count - number of entries to send/receive in initial rendezvous (must match on all ranks)
449d815da10SJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
450d815da10SJed Brown .  nto - number of ranks to send data to
451d815da10SJed Brown .  toranks - ranks to send to (array of length nto)
452d815da10SJed Brown .  todata - data to send to each rank (packed)
453d815da10SJed Brown .  ntags - number of tags needed by send/recv callbacks
454d815da10SJed Brown .  send - callback invoked on sending process when ready to send primary payload
455d815da10SJed Brown .  recv - callback invoked on receiving process after delivery of rendezvous message
456d815da10SJed Brown -  ctx - context for callbacks
457d815da10SJed Brown 
4584165533cSJose E. Roman    Output Parameters:
459d815da10SJed Brown +  nfrom - number of ranks receiving messages from
460d815da10SJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
461d815da10SJed Brown -  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
462d815da10SJed Brown 
463d815da10SJed Brown    Level: developer
464d815da10SJed Brown 
465d815da10SJed Brown    Notes:
466d815da10SJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
467d815da10SJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data.
468d815da10SJed Brown 
469d815da10SJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
470d815da10SJed Brown 
471d815da10SJed Brown    References:
472606c0280SSatish Balay .  * - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in
4734f02bc6aSBarry Smith    Scalable communication protocols for dynamic sparse data exchange, 2010.
474d815da10SJed Brown 
4756afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedFReq(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
476d815da10SJed Brown @*/
477d815da10SJed Brown PetscErrorCode PetscCommBuildTwoSidedF(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,
478d815da10SJed Brown                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
479d815da10SJed Brown                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
480d815da10SJed Brown {
4816afbe7ddSJed Brown   MPI_Request    *toreqs,*fromreqs;
4826afbe7ddSJed Brown 
4836afbe7ddSJed Brown   PetscFunctionBegin;
4845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommBuildTwoSidedFReq(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,&toreqs,&fromreqs,send,recv,ctx));
4855f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(nto*ntags,toreqs,MPI_STATUSES_IGNORE));
4865f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Waitall(*nfrom*ntags,fromreqs,MPI_STATUSES_IGNORE));
4875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(toreqs));
4885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(fromreqs));
4896afbe7ddSJed Brown   PetscFunctionReturn(0);
4906afbe7ddSJed Brown }
4916afbe7ddSJed Brown 
4926afbe7ddSJed Brown /*@C
4936afbe7ddSJed Brown    PetscCommBuildTwoSidedFReq - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous, returns requests
4946afbe7ddSJed Brown 
495d083f849SBarry Smith    Collective
4966afbe7ddSJed Brown 
4974165533cSJose E. Roman    Input Parameters:
4986afbe7ddSJed Brown +  comm - communicator
4996afbe7ddSJed Brown .  count - number of entries to send/receive in initial rendezvous (must match on all ranks)
5006afbe7ddSJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
5016afbe7ddSJed Brown .  nto - number of ranks to send data to
5026afbe7ddSJed Brown .  toranks - ranks to send to (array of length nto)
5036afbe7ddSJed Brown .  todata - data to send to each rank (packed)
5046afbe7ddSJed Brown .  ntags - number of tags needed by send/recv callbacks
5056afbe7ddSJed Brown .  send - callback invoked on sending process when ready to send primary payload
5066afbe7ddSJed Brown .  recv - callback invoked on receiving process after delivery of rendezvous message
5076afbe7ddSJed Brown -  ctx - context for callbacks
5086afbe7ddSJed Brown 
5094165533cSJose E. Roman    Output Parameters:
5106afbe7ddSJed Brown +  nfrom - number of ranks receiving messages from
5116afbe7ddSJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
5126afbe7ddSJed Brown .  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
5136afbe7ddSJed Brown .  toreqs - array of nto*ntags sender requests (caller must wait on these, then PetscFree())
5146afbe7ddSJed Brown -  fromreqs - array of nfrom*ntags receiver requests (caller must wait on these, then PetscFree())
5156afbe7ddSJed Brown 
5166afbe7ddSJed Brown    Level: developer
5176afbe7ddSJed Brown 
5186afbe7ddSJed Brown    Notes:
5196afbe7ddSJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
5206afbe7ddSJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data.
5216afbe7ddSJed Brown 
5226afbe7ddSJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
5236afbe7ddSJed Brown 
5246afbe7ddSJed Brown    References:
525606c0280SSatish Balay .  * - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in
5264f02bc6aSBarry Smith    Scalable communication protocols for dynamic sparse data exchange, 2010.
5276afbe7ddSJed Brown 
5286afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedF(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
5296afbe7ddSJed Brown @*/
5306afbe7ddSJed Brown PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,
5316afbe7ddSJed Brown                                           PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs,
5326afbe7ddSJed Brown                                           PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
5336afbe7ddSJed Brown                                           PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
5346afbe7ddSJed Brown {
5355f80ce2aSJacob Faibussowitsch   PetscErrorCode         (*f)(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,
5366afbe7ddSJed Brown                               PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,MPI_Request**,MPI_Request**,
537d815da10SJed Brown                               PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
538d815da10SJed Brown                               PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx);
539d815da10SJed Brown   PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET;
5403461502dSJed Brown   PetscMPIInt i,size;
541d815da10SJed Brown 
542d815da10SJed Brown   PetscFunctionBegin;
5435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSysInitializePackage());
5445f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
5453461502dSJed Brown   for (i=0; i<nto; i++) {
5462c71b3e2SJacob Faibussowitsch     PetscCheckFalse(toranks[i] < 0 || size <= toranks[i],comm,PETSC_ERR_ARG_OUTOFRANGE,"toranks[%d] %d not in comm size %d",i,toranks[i],size);
5473461502dSJed Brown   }
5485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventSync(PETSC_BuildTwoSidedF,comm));
5495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(PETSC_BuildTwoSidedF,0,0,0,0));
5505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommBuildTwoSidedGetType(comm,&buildtype));
551d815da10SJed Brown   switch (buildtype) {
552d815da10SJed Brown   case PETSC_BUILDTWOSIDED_IBARRIER:
5538ee3e7ecSJunchao Zhang #if defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
5546afbe7ddSJed Brown     f = PetscCommBuildTwoSidedFReq_Ibarrier;
555b458e8f1SJose E. Roman     break;
556afb254cdSJed Brown #else
557afb254cdSJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)");
558afb254cdSJed Brown #endif
559d815da10SJed Brown   case PETSC_BUILDTWOSIDED_ALLREDUCE:
560d815da10SJed Brown   case PETSC_BUILDTWOSIDED_REDSCATTER:
5616afbe7ddSJed Brown     f = PetscCommBuildTwoSidedFReq_Reference;
562d815da10SJed Brown     break;
563d815da10SJed Brown   default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication");
564d815da10SJed Brown   }
5655f80ce2aSJacob Faibussowitsch   CHKERRQ((*f)(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,toreqs,fromreqs,send,recv,ctx));
5665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(PETSC_BuildTwoSidedF,0,0,0,0));
567d815da10SJed Brown   PetscFunctionReturn(0);
568d815da10SJed Brown }
569