xref: /petsc/src/sys/utils/mpits.c (revision 05342407404a810c6afec3612a6850394c3ef881)
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 
236145cd65SJed Brown    Input Arguments:
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;
376145cd65SJed Brown #if defined(PETSC_USE_DEBUG)
386145cd65SJed Brown   {                             /* We don't have a PetscObject so can't use PetscValidLogicalCollectiveEnum */
396145cd65SJed Brown     PetscMPIInt ierr;
406145cd65SJed Brown     PetscMPIInt b1[2],b2[2];
416145cd65SJed Brown     b1[0] = -(PetscMPIInt)twosided;
426145cd65SJed Brown     b1[1] = (PetscMPIInt)twosided;
43b2566f29SBarry Smith     ierr  = MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
446145cd65SJed Brown     if (-b2[0] != b2[1]) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes");
456145cd65SJed Brown   }
466145cd65SJed Brown #endif
476145cd65SJed Brown   _twosided_type = twosided;
486145cd65SJed Brown   PetscFunctionReturn(0);
496145cd65SJed Brown }
506145cd65SJed Brown 
516145cd65SJed Brown /*@
526145cd65SJed Brown    PetscCommBuildTwoSidedGetType - set algorithm to use when building two-sided communication
536145cd65SJed Brown 
546145cd65SJed Brown    Logically Collective
556145cd65SJed Brown 
566145cd65SJed Brown    Output Arguments:
576145cd65SJed Brown +  comm - communicator on which to query algorithm
586145cd65SJed Brown -  twosided - algorithm to use for PetscCommBuildTwoSided()
596145cd65SJed Brown 
606145cd65SJed Brown    Level: developer
616145cd65SJed Brown 
626145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType()
636145cd65SJed Brown @*/
646145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm comm,PetscBuildTwoSidedType *twosided)
65f6ced4a3SJed Brown {
66f6ced4a3SJed Brown   PetscErrorCode ierr;
67*05342407SJunchao Zhang   PetscMPIInt    size;
68f6ced4a3SJed Brown 
69f6ced4a3SJed Brown   PetscFunctionBegin;
706145cd65SJed Brown   *twosided = PETSC_BUILDTWOSIDED_NOTSET;
716145cd65SJed Brown   if (_twosided_type == PETSC_BUILDTWOSIDED_NOTSET) {
72*05342407SJunchao Zhang     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
73*05342407SJunchao Zhang     _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; /* default for small comms, see https://gitlab.com/petsc/petsc/-/merge_requests/2611 */
74f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER)
75*05342407SJunchao Zhang     if (size > 1024) _twosided_type = PETSC_BUILDTWOSIDED_IBARRIER;
766145cd65SJed Brown #endif
77c5929fdfSBarry Smith     ierr = PetscOptionsGetEnum(NULL,NULL,"-build_twosided",PetscBuildTwoSidedTypes,(PetscEnum*)&_twosided_type,NULL);CHKERRQ(ierr);
78f6ced4a3SJed Brown   }
79f6ced4a3SJed Brown   *twosided = _twosided_type;
80f6ced4a3SJed Brown   PetscFunctionReturn(0);
81f6ced4a3SJed Brown }
82f6ced4a3SJed Brown 
83b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER)
84f6ced4a3SJed Brown 
85cf4b5b4fSJed 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)
86f6ced4a3SJed Brown {
87f6ced4a3SJed Brown   PetscErrorCode ierr;
88c3a59b84SJed Brown   PetscMPIInt    nrecvs,tag,done,i;
89c3a59b84SJed Brown   MPI_Aint       lb,unitbytes;
90f6ced4a3SJed Brown   char           *tdata;
91f6ced4a3SJed Brown   MPI_Request    *sendreqs,barrier;
920f453b92SJed Brown   PetscSegBuffer segrank,segdata;
9346e50bb6SJed Brown   PetscBool      barrier_started;
94f6ced4a3SJed Brown 
95f6ced4a3SJed Brown   PetscFunctionBegin;
96a502f807SJed Brown   ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr);
97c3a59b84SJed Brown   ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr);
98c3a59b84SJed Brown   if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
99f6ced4a3SJed Brown   tdata = (char*)todata;
100785e854fSJed Brown   ierr  = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr);
101f6ced4a3SJed Brown   for (i=0; i<nto; i++) {
102f6ced4a3SJed Brown     ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr);
103f6ced4a3SJed Brown   }
1040f453b92SJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr);
1050f453b92SJed Brown   ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr);
106f6ced4a3SJed Brown 
107f6ced4a3SJed Brown   nrecvs  = 0;
108f6ced4a3SJed Brown   barrier = MPI_REQUEST_NULL;
10946e50bb6SJed Brown   /* MPICH-3.2 sometimes does not create a request in some "optimized" cases.  This is arguably a standard violation,
11046e50bb6SJed Brown    * but we need to work around it. */
11146e50bb6SJed Brown   barrier_started = PETSC_FALSE;
112f6ced4a3SJed Brown   for (done=0; !done; ) {
113f6ced4a3SJed Brown     PetscMPIInt flag;
114f6ced4a3SJed Brown     MPI_Status  status;
115f6ced4a3SJed Brown     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr);
116f6ced4a3SJed Brown     if (flag) {                 /* incoming message */
117f6ced4a3SJed Brown       PetscMPIInt *recvrank;
118f6ced4a3SJed Brown       void        *buf;
119137cf7b6SJed Brown       ierr      = PetscSegBufferGet(segrank,1,&recvrank);CHKERRQ(ierr);
120137cf7b6SJed Brown       ierr      = PetscSegBufferGet(segdata,count,&buf);CHKERRQ(ierr);
121f6ced4a3SJed Brown       *recvrank = status.MPI_SOURCE;
122f6ced4a3SJed Brown       ierr      = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr);
123f6ced4a3SJed Brown       nrecvs++;
124f6ced4a3SJed Brown     }
12546e50bb6SJed Brown     if (!barrier_started) {
1264dc2109aSBarry Smith       PetscMPIInt sent,nsends;
1274dc2109aSBarry Smith       ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr);
128f6ced4a3SJed Brown       ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
129f6ced4a3SJed Brown       if (sent) {
130b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER)
131f6ced4a3SJed Brown         ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr);
132b5a01e9fSJed Brown #elif defined(PETSC_HAVE_MPIX_IBARRIER)
133b5a01e9fSJed Brown         ierr = MPIX_Ibarrier(comm,&barrier);CHKERRQ(ierr);
134b5a01e9fSJed Brown #endif
13546e50bb6SJed Brown         barrier_started = PETSC_TRUE;
136f6ced4a3SJed Brown         ierr = PetscFree(sendreqs);CHKERRQ(ierr);
137f6ced4a3SJed Brown       }
138f6ced4a3SJed Brown     } else {
139f6ced4a3SJed Brown       ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr);
140f6ced4a3SJed Brown     }
141f6ced4a3SJed Brown   }
142f6ced4a3SJed Brown   *nfrom = nrecvs;
143137cf7b6SJed Brown   ierr   = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr);
1440f453b92SJed Brown   ierr   = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr);
145137cf7b6SJed Brown   ierr   = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr);
1460f453b92SJed Brown   ierr   = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr);
147a502f807SJed Brown   ierr   = PetscCommDestroy(&comm);CHKERRQ(ierr);
148f6ced4a3SJed Brown   PetscFunctionReturn(0);
149f6ced4a3SJed Brown }
150f6ced4a3SJed Brown #endif
151f6ced4a3SJed Brown 
152cf4b5b4fSJed 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)
153f6ced4a3SJed Brown {
154f6ced4a3SJed Brown   PetscErrorCode   ierr;
155*05342407SJunchao Zhang   PetscMPIInt      size,rank,*iflags,nrecvs,tag,*franks,i,flg;
156c3a59b84SJed Brown   MPI_Aint         lb,unitbytes;
157f6ced4a3SJed Brown   char             *tdata,*fdata;
158f6ced4a3SJed Brown   MPI_Request      *reqs,*sendreqs;
159f6ced4a3SJed Brown   MPI_Status       *statuses;
160*05342407SJunchao Zhang   PetscCommCounter *counter;
161f6ced4a3SJed Brown 
162f6ced4a3SJed Brown   PetscFunctionBegin;
163f6ced4a3SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
164*05342407SJunchao Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
165a502f807SJed Brown   ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr);
166*05342407SJunchao Zhang   ierr = MPI_Comm_get_attr(comm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
167*05342407SJunchao Zhang   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner PETSc communicator does not have its tag/name counter attribute set");
168*05342407SJunchao Zhang   if (!counter->iflags) {
169*05342407SJunchao Zhang     ierr   = PetscCalloc1(size,&counter->iflags);CHKERRQ(ierr);
170*05342407SJunchao Zhang     iflags = counter->iflags;
171*05342407SJunchao Zhang   } else {
172*05342407SJunchao Zhang     iflags = counter->iflags;
173*05342407SJunchao Zhang     ierr   = PetscArrayzero(iflags,size);CHKERRQ(ierr);
174*05342407SJunchao Zhang   }
175*05342407SJunchao Zhang   for (i=0; i<nto; i++) iflags[toranks[i]] = 1;
176*05342407SJunchao Zhang   ierr     = MPIU_Allreduce(MPI_IN_PLACE,iflags,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
177*05342407SJunchao Zhang   nrecvs   = iflags[rank];
178c3a59b84SJed Brown   ierr     = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr);
179c3a59b84SJed Brown   if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
180f6ced4a3SJed Brown   ierr     = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr);
181f6ced4a3SJed Brown   tdata    = (char*)todata;
182dcca6d9dSJed Brown   ierr     = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr);
18306926cafSJed Brown   sendreqs = reqs + nrecvs;
184f6ced4a3SJed Brown   for (i=0; i<nrecvs; i++) {
185f6ced4a3SJed Brown     ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr);
186f6ced4a3SJed Brown   }
187f6ced4a3SJed Brown   for (i=0; i<nto; i++) {
188f6ced4a3SJed Brown     ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr);
189f6ced4a3SJed Brown   }
190f6ced4a3SJed Brown   ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr);
191785e854fSJed Brown   ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr);
192a297a907SKarl Rupp   for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE;
19306926cafSJed Brown   ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr);
194a502f807SJed Brown   ierr = PetscCommDestroy(&comm);CHKERRQ(ierr);
195f6ced4a3SJed Brown 
196f6ced4a3SJed Brown   *nfrom            = nrecvs;
197f6ced4a3SJed Brown   *fromranks        = franks;
198f6ced4a3SJed Brown   *(void**)fromdata = fdata;
199f6ced4a3SJed Brown   PetscFunctionReturn(0);
200f6ced4a3SJed Brown }
201f6ced4a3SJed Brown 
2021654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
203cf4b5b4fSJed 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)
2041654bf6bSJed Brown {
2051654bf6bSJed Brown   PetscErrorCode ierr;
206*05342407SJunchao Zhang   PetscMPIInt    size,*iflags,nrecvs,tag,*franks,i,flg;
207c3a59b84SJed Brown   MPI_Aint       lb,unitbytes;
2081654bf6bSJed Brown   char           *tdata,*fdata;
2091654bf6bSJed Brown   MPI_Request    *reqs,*sendreqs;
2101654bf6bSJed Brown   MPI_Status     *statuses;
211*05342407SJunchao Zhang   PetscCommCounter *counter;
2121654bf6bSJed Brown 
2131654bf6bSJed Brown   PetscFunctionBegin;
2141654bf6bSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
215*05342407SJunchao Zhang   ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr);
216*05342407SJunchao Zhang   ierr = MPI_Comm_get_attr(comm,Petsc_Counter_keyval,&counter,&flg);CHKERRQ(ierr);
217*05342407SJunchao Zhang   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner PETSc communicator does not have its tag/name counter attribute set");
218*05342407SJunchao Zhang   if (!counter->iflags) {
219*05342407SJunchao Zhang     ierr   = PetscCalloc1(size,&counter->iflags);CHKERRQ(ierr);
220*05342407SJunchao Zhang     iflags = counter->iflags;
221*05342407SJunchao Zhang   } else {
222*05342407SJunchao Zhang     iflags = counter->iflags;
223580bdb30SBarry Smith     ierr   = PetscArrayzero(iflags,size);CHKERRQ(ierr);
224*05342407SJunchao Zhang   }
2251654bf6bSJed Brown   for (i=0; i<nto; i++) iflags[toranks[i]] = 1;
2261654bf6bSJed Brown   ierr     = MPI_Reduce_scatter_block(iflags,&nrecvs,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
227c3a59b84SJed Brown   ierr     = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr);
228c3a59b84SJed Brown   if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
2291654bf6bSJed Brown   ierr     = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr);
2301654bf6bSJed Brown   tdata    = (char*)todata;
2311654bf6bSJed Brown   ierr     = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr);
2321654bf6bSJed Brown   sendreqs = reqs + nrecvs;
2331654bf6bSJed Brown   for (i=0; i<nrecvs; i++) {
2341654bf6bSJed Brown     ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr);
2351654bf6bSJed Brown   }
2361654bf6bSJed Brown   for (i=0; i<nto; i++) {
2371654bf6bSJed Brown     ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr);
2381654bf6bSJed Brown   }
2391654bf6bSJed Brown   ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr);
2401654bf6bSJed Brown   ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr);
2411654bf6bSJed Brown   for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE;
2421654bf6bSJed Brown   ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr);
243a502f807SJed Brown   ierr = PetscCommDestroy(&comm);CHKERRQ(ierr);
2441654bf6bSJed Brown 
2451654bf6bSJed Brown   *nfrom            = nrecvs;
2461654bf6bSJed Brown   *fromranks        = franks;
2471654bf6bSJed Brown   *(void**)fromdata = fdata;
2481654bf6bSJed Brown   PetscFunctionReturn(0);
2491654bf6bSJed Brown }
2501654bf6bSJed Brown #endif
2511654bf6bSJed Brown 
252f6ced4a3SJed Brown /*@C
253f6ced4a3SJed Brown    PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths)
254f6ced4a3SJed Brown 
255d083f849SBarry Smith    Collective
256f6ced4a3SJed Brown 
257f6ced4a3SJed Brown    Input Arguments:
258f6ced4a3SJed Brown +  comm - communicator
259f6ced4a3SJed Brown .  count - number of entries to send/receive (must match on all ranks)
260f6ced4a3SJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
261f6ced4a3SJed Brown .  nto - number of ranks to send data to
262f6ced4a3SJed Brown .  toranks - ranks to send to (array of length nto)
263f6ced4a3SJed Brown -  todata - data to send to each rank (packed)
264f6ced4a3SJed Brown 
265f6ced4a3SJed Brown    Output Arguments:
266f6ced4a3SJed Brown +  nfrom - number of ranks receiving messages from
267f6ced4a3SJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
268f6ced4a3SJed Brown -  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
269f6ced4a3SJed Brown 
270f6ced4a3SJed Brown    Level: developer
271f6ced4a3SJed Brown 
2721654bf6bSJed Brown    Options Database Keys:
273*05342407SJunchao Zhang .  -build_twosided <allreduce|ibarrier|redscatter> - algorithm to set up two-sided communication. Default is allreduce for communicators with <= 1024 ranks, otherwise ibarrier.
2741654bf6bSJed Brown 
275f6ced4a3SJed Brown    Notes:
276f6ced4a3SJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
277f6ced4a3SJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data.
278f6ced4a3SJed Brown 
279f6ced4a3SJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
280f6ced4a3SJed Brown 
281f6ced4a3SJed Brown    References:
28296a0c994SBarry Smith .  1. - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in
2834f02bc6aSBarry Smith    Scalable communication protocols for dynamic sparse data exchange, 2010.
284f6ced4a3SJed Brown 
285f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
286f6ced4a3SJed Brown @*/
287cf4b5b4fSJed 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)
288f6ced4a3SJed Brown {
289f6ced4a3SJed Brown   PetscErrorCode         ierr;
2906145cd65SJed Brown   PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET;
291f6ced4a3SJed Brown 
292f6ced4a3SJed Brown   PetscFunctionBegin;
2933b3561c8SJed Brown   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
29462872c28SLisandro Dalcin   ierr = PetscLogEventSync(PETSC_BuildTwoSided,comm);CHKERRQ(ierr);
2953b3561c8SJed Brown   ierr = PetscLogEventBegin(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr);
2966145cd65SJed Brown   ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr);
297f6ced4a3SJed Brown   switch (buildtype) {
2986145cd65SJed Brown   case PETSC_BUILDTWOSIDED_IBARRIER:
299b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER)
3006145cd65SJed Brown     ierr = PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr);
3016145cd65SJed Brown #else
3026145cd65SJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)");
303f6ced4a3SJed Brown #endif
3046145cd65SJed Brown     break;
3056145cd65SJed Brown   case PETSC_BUILDTWOSIDED_ALLREDUCE:
3066145cd65SJed Brown     ierr = PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr);
307f6ced4a3SJed Brown     break;
3081654bf6bSJed Brown   case PETSC_BUILDTWOSIDED_REDSCATTER:
3091654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK)
3101654bf6bSJed Brown     ierr = PetscCommBuildTwoSided_RedScatter(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr);
3111654bf6bSJed Brown #else
3121654bf6bSJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Reduce_scatter_block (part of MPI-2.2)");
3131654bf6bSJed Brown #endif
3141654bf6bSJed Brown     break;
315f6ced4a3SJed Brown   default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication");
316f6ced4a3SJed Brown   }
3173b3561c8SJed Brown   ierr = PetscLogEventEnd(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr);
318f6ced4a3SJed Brown   PetscFunctionReturn(0);
319f6ced4a3SJed Brown }
320d815da10SJed Brown 
3216afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Reference(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,
3226afbe7ddSJed Brown                                                            PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs,
323d815da10SJed Brown                                                            PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
324d815da10SJed Brown                                                            PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
325d815da10SJed Brown {
326d815da10SJed Brown   PetscErrorCode ierr;
327d815da10SJed Brown   PetscMPIInt i,*tag;
328d815da10SJed Brown   MPI_Aint    lb,unitbytes;
329d815da10SJed Brown   MPI_Request *sendreq,*recvreq;
330d815da10SJed Brown 
331d815da10SJed Brown   PetscFunctionBegin;
332d815da10SJed Brown   ierr = PetscMalloc1(ntags,&tag);CHKERRQ(ierr);
333d815da10SJed Brown   if (ntags > 0) {
334d815da10SJed Brown     ierr = PetscCommDuplicate(comm,&comm,&tag[0]);CHKERRQ(ierr);
335d815da10SJed Brown   }
336d815da10SJed Brown   for (i=1; i<ntags; i++) {
337d815da10SJed Brown     ierr = PetscCommGetNewTag(comm,&tag[i]);CHKERRQ(ierr);
338d815da10SJed Brown   }
339d815da10SJed Brown 
340d815da10SJed Brown   /* Perform complete initial rendezvous */
341d815da10SJed Brown   ierr = PetscCommBuildTwoSided(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr);
342d815da10SJed Brown 
3436afbe7ddSJed Brown   ierr = PetscMalloc1(nto*ntags,&sendreq);CHKERRQ(ierr);
3446afbe7ddSJed Brown   ierr = PetscMalloc1(*nfrom*ntags,&recvreq);CHKERRQ(ierr);
345d815da10SJed Brown 
346d815da10SJed Brown   ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr);
347d815da10SJed Brown   if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
348d815da10SJed Brown   for (i=0; i<nto; i++) {
34989157a57SJed Brown     PetscMPIInt k;
35089157a57SJed Brown     for (k=0; k<ntags; k++) sendreq[i*ntags+k] = MPI_REQUEST_NULL;
351d815da10SJed Brown     ierr = (*send)(comm,tag,i,toranks[i],((char*)todata)+count*unitbytes*i,sendreq+i*ntags,ctx);CHKERRQ(ierr);
352d815da10SJed Brown   }
353d815da10SJed Brown   for (i=0; i<*nfrom; i++) {
354d815da10SJed Brown     void *header = (*(char**)fromdata) + count*unitbytes*i;
35589157a57SJed Brown     PetscMPIInt k;
35689157a57SJed Brown     for (k=0; k<ntags; k++) recvreq[i*ntags+k] = MPI_REQUEST_NULL;
357d815da10SJed Brown     ierr = (*recv)(comm,tag,(*fromranks)[i],header,recvreq+i*ntags,ctx);CHKERRQ(ierr);
358d815da10SJed Brown   }
359d815da10SJed Brown   ierr = PetscFree(tag);CHKERRQ(ierr);
360d815da10SJed Brown   ierr = PetscCommDestroy(&comm);CHKERRQ(ierr);
3616afbe7ddSJed Brown   *toreqs = sendreq;
3626afbe7ddSJed Brown   *fromreqs = recvreq;
363d815da10SJed Brown   PetscFunctionReturn(0);
364d815da10SJed Brown }
365d815da10SJed Brown 
366b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER)
367afb254cdSJed Brown 
3686afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,
3696afbe7ddSJed Brown                                                           PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs,
370afb254cdSJed Brown                                                           PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
371afb254cdSJed Brown                                                           PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
372afb254cdSJed Brown {
373afb254cdSJed Brown   PetscErrorCode ierr;
374afb254cdSJed Brown   PetscMPIInt    nrecvs,tag,*tags,done,i;
375afb254cdSJed Brown   MPI_Aint       lb,unitbytes;
376afb254cdSJed Brown   char           *tdata;
377afb254cdSJed Brown   MPI_Request    *sendreqs,*usendreqs,*req,barrier;
378afb254cdSJed Brown   PetscSegBuffer segrank,segdata,segreq;
37946e50bb6SJed Brown   PetscBool      barrier_started;
380afb254cdSJed Brown 
381afb254cdSJed Brown   PetscFunctionBegin;
382afb254cdSJed Brown   ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr);
383afb254cdSJed Brown   ierr = PetscMalloc1(ntags,&tags);CHKERRQ(ierr);
384afb254cdSJed Brown   for (i=0; i<ntags; i++) {
385afb254cdSJed Brown     ierr = PetscCommGetNewTag(comm,&tags[i]);CHKERRQ(ierr);
386afb254cdSJed Brown   }
387afb254cdSJed Brown   ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr);
388afb254cdSJed Brown   if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb);
389afb254cdSJed Brown   tdata = (char*)todata;
3906afbe7ddSJed Brown   ierr = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr);
3916afbe7ddSJed Brown   ierr = PetscMalloc1(nto*ntags,&usendreqs);CHKERRQ(ierr);
392afb254cdSJed Brown   /* Post synchronous sends */
393afb254cdSJed Brown   for (i=0; i<nto; i++) {
394afb254cdSJed Brown     ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr);
395afb254cdSJed Brown   }
396afb254cdSJed Brown   /* Post actual payloads.  These are typically larger messages.  Hopefully sending these later does not slow down the
397afb254cdSJed Brown    * synchronous messages above. */
398afb254cdSJed Brown   for (i=0; i<nto; i++) {
39989157a57SJed Brown     PetscMPIInt k;
40089157a57SJed Brown     for (k=0; k<ntags; k++) usendreqs[i*ntags+k] = MPI_REQUEST_NULL;
401afb254cdSJed Brown     ierr = (*send)(comm,tags,i,toranks[i],tdata+count*unitbytes*i,usendreqs+i*ntags,ctx);CHKERRQ(ierr);
402afb254cdSJed Brown   }
403afb254cdSJed Brown 
404afb254cdSJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr);
405afb254cdSJed Brown   ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr);
406afb254cdSJed Brown   ierr = PetscSegBufferCreate(sizeof(MPI_Request),4,&segreq);CHKERRQ(ierr);
407afb254cdSJed Brown 
408afb254cdSJed Brown   nrecvs  = 0;
409afb254cdSJed Brown   barrier = MPI_REQUEST_NULL;
41046e50bb6SJed Brown   /* MPICH-3.2 sometimes does not create a request in some "optimized" cases.  This is arguably a standard violation,
41146e50bb6SJed Brown    * but we need to work around it. */
41246e50bb6SJed Brown   barrier_started = PETSC_FALSE;
413afb254cdSJed Brown   for (done=0; !done; ) {
414afb254cdSJed Brown     PetscMPIInt flag;
415afb254cdSJed Brown     MPI_Status  status;
416afb254cdSJed Brown     ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr);
417afb254cdSJed Brown     if (flag) {                 /* incoming message */
41889157a57SJed Brown       PetscMPIInt *recvrank,k;
419afb254cdSJed Brown       void        *buf;
420afb254cdSJed Brown       ierr = PetscSegBufferGet(segrank,1,&recvrank);CHKERRQ(ierr);
421afb254cdSJed Brown       ierr = PetscSegBufferGet(segdata,count,&buf);CHKERRQ(ierr);
422afb254cdSJed Brown       *recvrank = status.MPI_SOURCE;
423afb254cdSJed Brown       ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr);
424afb254cdSJed Brown       ierr = PetscSegBufferGet(segreq,ntags,&req);CHKERRQ(ierr);
42589157a57SJed Brown       for (k=0; k<ntags; k++) req[k] = MPI_REQUEST_NULL;
426afb254cdSJed Brown       ierr = (*recv)(comm,tags,status.MPI_SOURCE,buf,req,ctx);CHKERRQ(ierr);
427afb254cdSJed Brown       nrecvs++;
428afb254cdSJed Brown     }
42946e50bb6SJed Brown     if (!barrier_started) {
430afb254cdSJed Brown       PetscMPIInt sent,nsends;
431afb254cdSJed Brown       ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr);
432afb254cdSJed Brown       ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
433afb254cdSJed Brown       if (sent) {
434b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER)
435afb254cdSJed Brown         ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr);
436b5a01e9fSJed Brown #elif defined(PETSC_HAVE_MPIX_IBARRIER)
437b5a01e9fSJed Brown         ierr = MPIX_Ibarrier(comm,&barrier);CHKERRQ(ierr);
438b5a01e9fSJed Brown #endif
43946e50bb6SJed Brown         barrier_started = PETSC_TRUE;
440afb254cdSJed Brown       }
441afb254cdSJed Brown     } else {
442afb254cdSJed Brown       ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr);
443afb254cdSJed Brown     }
444afb254cdSJed Brown   }
445afb254cdSJed Brown   *nfrom = nrecvs;
446afb254cdSJed Brown   ierr = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr);
447afb254cdSJed Brown   ierr = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr);
448afb254cdSJed Brown   ierr = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr);
449afb254cdSJed Brown   ierr = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr);
4506afbe7ddSJed Brown   *toreqs = usendreqs;
4516afbe7ddSJed Brown   ierr = PetscSegBufferExtractAlloc(segreq,fromreqs);CHKERRQ(ierr);
452afb254cdSJed Brown   ierr = PetscSegBufferDestroy(&segreq);CHKERRQ(ierr);
4536afbe7ddSJed Brown   ierr = PetscFree(sendreqs);CHKERRQ(ierr);
454afb254cdSJed Brown   ierr = PetscFree(tags);CHKERRQ(ierr);
455afb254cdSJed Brown   ierr = PetscCommDestroy(&comm);CHKERRQ(ierr);
456afb254cdSJed Brown   PetscFunctionReturn(0);
457afb254cdSJed Brown }
458afb254cdSJed Brown #endif
459afb254cdSJed Brown 
460d815da10SJed Brown /*@C
461d815da10SJed Brown    PetscCommBuildTwoSidedF - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous
462d815da10SJed Brown 
463d083f849SBarry Smith    Collective
464d815da10SJed Brown 
465d815da10SJed Brown    Input Arguments:
466d815da10SJed Brown +  comm - communicator
467d815da10SJed Brown .  count - number of entries to send/receive in initial rendezvous (must match on all ranks)
468d815da10SJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
469d815da10SJed Brown .  nto - number of ranks to send data to
470d815da10SJed Brown .  toranks - ranks to send to (array of length nto)
471d815da10SJed Brown .  todata - data to send to each rank (packed)
472d815da10SJed Brown .  ntags - number of tags needed by send/recv callbacks
473d815da10SJed Brown .  send - callback invoked on sending process when ready to send primary payload
474d815da10SJed Brown .  recv - callback invoked on receiving process after delivery of rendezvous message
475d815da10SJed Brown -  ctx - context for callbacks
476d815da10SJed Brown 
477d815da10SJed Brown    Output Arguments:
478d815da10SJed Brown +  nfrom - number of ranks receiving messages from
479d815da10SJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
480d815da10SJed Brown -  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
481d815da10SJed Brown 
482d815da10SJed Brown    Level: developer
483d815da10SJed Brown 
484d815da10SJed Brown    Notes:
485d815da10SJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
486d815da10SJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data.
487d815da10SJed Brown 
488d815da10SJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
489d815da10SJed Brown 
490d815da10SJed Brown    References:
49196a0c994SBarry Smith .  1. - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in
4924f02bc6aSBarry Smith    Scalable communication protocols for dynamic sparse data exchange, 2010.
493d815da10SJed Brown 
4946afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedFReq(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
495d815da10SJed Brown @*/
496d815da10SJed 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,
497d815da10SJed Brown                                        PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
498d815da10SJed Brown                                        PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
499d815da10SJed Brown {
5006afbe7ddSJed Brown   PetscErrorCode ierr;
5016afbe7ddSJed Brown   MPI_Request    *toreqs,*fromreqs;
5026afbe7ddSJed Brown 
5036afbe7ddSJed Brown   PetscFunctionBegin;
5046afbe7ddSJed Brown   ierr = PetscCommBuildTwoSidedFReq(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,&toreqs,&fromreqs,send,recv,ctx);CHKERRQ(ierr);
5056afbe7ddSJed Brown   ierr = MPI_Waitall(nto*ntags,toreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
5066afbe7ddSJed Brown   ierr = MPI_Waitall(*nfrom*ntags,fromreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
5076afbe7ddSJed Brown   ierr = PetscFree(toreqs);CHKERRQ(ierr);
5086afbe7ddSJed Brown   ierr = PetscFree(fromreqs);CHKERRQ(ierr);
5096afbe7ddSJed Brown   PetscFunctionReturn(0);
5106afbe7ddSJed Brown }
5116afbe7ddSJed Brown 
5126afbe7ddSJed Brown /*@C
5136afbe7ddSJed Brown    PetscCommBuildTwoSidedFReq - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous, returns requests
5146afbe7ddSJed Brown 
515d083f849SBarry Smith    Collective
5166afbe7ddSJed Brown 
5176afbe7ddSJed Brown    Input Arguments:
5186afbe7ddSJed Brown +  comm - communicator
5196afbe7ddSJed Brown .  count - number of entries to send/receive in initial rendezvous (must match on all ranks)
5206afbe7ddSJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
5216afbe7ddSJed Brown .  nto - number of ranks to send data to
5226afbe7ddSJed Brown .  toranks - ranks to send to (array of length nto)
5236afbe7ddSJed Brown .  todata - data to send to each rank (packed)
5246afbe7ddSJed Brown .  ntags - number of tags needed by send/recv callbacks
5256afbe7ddSJed Brown .  send - callback invoked on sending process when ready to send primary payload
5266afbe7ddSJed Brown .  recv - callback invoked on receiving process after delivery of rendezvous message
5276afbe7ddSJed Brown -  ctx - context for callbacks
5286afbe7ddSJed Brown 
5296afbe7ddSJed Brown    Output Arguments:
5306afbe7ddSJed Brown +  nfrom - number of ranks receiving messages from
5316afbe7ddSJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
5326afbe7ddSJed Brown .  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
5336afbe7ddSJed Brown .  toreqs - array of nto*ntags sender requests (caller must wait on these, then PetscFree())
5346afbe7ddSJed Brown -  fromreqs - array of nfrom*ntags receiver requests (caller must wait on these, then PetscFree())
5356afbe7ddSJed Brown 
5366afbe7ddSJed Brown    Level: developer
5376afbe7ddSJed Brown 
5386afbe7ddSJed Brown    Notes:
5396afbe7ddSJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
5406afbe7ddSJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data.
5416afbe7ddSJed Brown 
5426afbe7ddSJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
5436afbe7ddSJed Brown 
5446afbe7ddSJed Brown    References:
54596a0c994SBarry Smith .  1. - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in
5464f02bc6aSBarry Smith    Scalable communication protocols for dynamic sparse data exchange, 2010.
5476afbe7ddSJed Brown 
5486afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedF(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
5496afbe7ddSJed Brown @*/
5506afbe7ddSJed Brown PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata,
5516afbe7ddSJed Brown                                           PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs,
5526afbe7ddSJed Brown                                           PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
5536afbe7ddSJed Brown                                           PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx)
5546afbe7ddSJed Brown {
5556afbe7ddSJed Brown   PetscErrorCode         ierr,(*f)(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*,
5566afbe7ddSJed Brown                                    PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,MPI_Request**,MPI_Request**,
557d815da10SJed Brown                                    PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*),
558d815da10SJed Brown                                    PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx);
559d815da10SJed Brown   PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET;
5603461502dSJed Brown   PetscMPIInt i,size;
561d815da10SJed Brown 
562d815da10SJed Brown   PetscFunctionBegin;
563d815da10SJed Brown   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
5643461502dSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
5653461502dSJed Brown   for (i=0; i<nto; i++) {
5663461502dSJed Brown     if (toranks[i] < 0 || size <= toranks[i]) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"toranks[%d] %d not in comm size %d",i,toranks[i],size);
5673461502dSJed Brown   }
56862872c28SLisandro Dalcin   ierr = PetscLogEventSync(PETSC_BuildTwoSidedF,comm);CHKERRQ(ierr);
56927a32ea5SJed Brown   ierr = PetscLogEventBegin(PETSC_BuildTwoSidedF,0,0,0,0);CHKERRQ(ierr);
570d815da10SJed Brown   ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr);
571d815da10SJed Brown   switch (buildtype) {
572d815da10SJed Brown   case PETSC_BUILDTWOSIDED_IBARRIER:
573b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER)
5746afbe7ddSJed Brown     f = PetscCommBuildTwoSidedFReq_Ibarrier;
575afb254cdSJed Brown #else
576afb254cdSJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)");
577afb254cdSJed Brown #endif
578afb254cdSJed Brown     break;
579d815da10SJed Brown   case PETSC_BUILDTWOSIDED_ALLREDUCE:
580d815da10SJed Brown   case PETSC_BUILDTWOSIDED_REDSCATTER:
5816afbe7ddSJed Brown     f = PetscCommBuildTwoSidedFReq_Reference;
582d815da10SJed Brown     break;
583d815da10SJed Brown   default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication");
584d815da10SJed Brown   }
5856afbe7ddSJed Brown   ierr = (*f)(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,toreqs,fromreqs,send,recv,ctx);CHKERRQ(ierr);
58627a32ea5SJed Brown   ierr = PetscLogEventEnd(PETSC_BuildTwoSidedF,0,0,0,0);CHKERRQ(ierr);
587d815da10SJed Brown   PetscFunctionReturn(0);
588d815da10SJed Brown }
589