xref: /petsc/src/sys/utils/mpits.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
1f6ced4a3SJed Brown #include <petscsys.h>        /*I  "petscsys.h"  I*/
2f6ced4a3SJed Brown 
36145cd65SJed Brown const char *const PetscBuildTwoSidedTypes[] = {
4f6ced4a3SJed Brown   "ALLREDUCE",
56145cd65SJed Brown   "IBARRIER",
66145cd65SJed Brown   "PetscBuildTwoSidedType",
76145cd65SJed Brown   "PETSC_BUILDTWOSIDED_",
8f6ced4a3SJed Brown   0
9f6ced4a3SJed Brown };
10f6ced4a3SJed Brown 
116145cd65SJed Brown static PetscBuildTwoSidedType _twosided_type = PETSC_BUILDTWOSIDED_NOTSET;
12f6ced4a3SJed Brown 
13f6ced4a3SJed Brown #undef __FUNCT__
146145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedSetType"
156145cd65SJed Brown /*@
166145cd65SJed Brown    PetscCommBuildTwoSidedSetType - set algorithm to use when building two-sided communication
176145cd65SJed Brown 
186145cd65SJed Brown    Logically Collective
196145cd65SJed Brown 
206145cd65SJed Brown    Input Arguments:
216145cd65SJed Brown +  comm - PETSC_COMM_WORLD
226145cd65SJed Brown -  twosided - algorithm to use in subsequent calls to PetscCommBuildTwoSided()
236145cd65SJed Brown 
246145cd65SJed Brown    Level: developer
256145cd65SJed Brown 
266145cd65SJed Brown    Note:
276145cd65SJed Brown    This option is currently global, but could be made per-communicator.
286145cd65SJed Brown 
296145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedGetType()
306145cd65SJed Brown @*/
316145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm comm,PetscBuildTwoSidedType twosided)
326145cd65SJed Brown {
336145cd65SJed Brown   PetscFunctionBegin;
346145cd65SJed Brown #if defined(PETSC_USE_DEBUG)
356145cd65SJed Brown   {                             /* We don't have a PetscObject so can't use PetscValidLogicalCollectiveEnum */
366145cd65SJed Brown     PetscMPIInt ierr;
376145cd65SJed Brown     PetscMPIInt b1[2],b2[2];
386145cd65SJed Brown     b1[0] = -(PetscMPIInt)twosided;
396145cd65SJed Brown     b1[1] = (PetscMPIInt)twosided;
406145cd65SJed Brown     ierr  = MPI_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr);
416145cd65SJed Brown     if (-b2[0] != b2[1]) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes");
426145cd65SJed Brown   }
436145cd65SJed Brown #endif
446145cd65SJed Brown   _twosided_type = twosided;
456145cd65SJed Brown   PetscFunctionReturn(0);
466145cd65SJed Brown }
476145cd65SJed Brown 
486145cd65SJed Brown #undef __FUNCT__
496145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedGetType"
506145cd65SJed Brown /*@
516145cd65SJed Brown    PetscCommBuildTwoSidedGetType - set algorithm to use when building two-sided communication
526145cd65SJed Brown 
536145cd65SJed Brown    Logically Collective
546145cd65SJed Brown 
556145cd65SJed Brown    Output Arguments:
566145cd65SJed Brown +  comm - communicator on which to query algorithm
576145cd65SJed Brown -  twosided - algorithm to use for PetscCommBuildTwoSided()
586145cd65SJed Brown 
596145cd65SJed Brown    Level: developer
606145cd65SJed Brown 
616145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType()
626145cd65SJed Brown @*/
636145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm comm,PetscBuildTwoSidedType *twosided)
64f6ced4a3SJed Brown {
65f6ced4a3SJed Brown   PetscErrorCode ierr;
66f6ced4a3SJed Brown 
67f6ced4a3SJed Brown   PetscFunctionBegin;
686145cd65SJed Brown   *twosided = PETSC_BUILDTWOSIDED_NOTSET;
696145cd65SJed Brown   if (_twosided_type == PETSC_BUILDTWOSIDED_NOTSET) {
70f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER)
716145cd65SJed Brown #  if defined(PETSC_HAVE_MPICH_CH3_SOCK) && !defined(PETSC_HAVE_MPICH_CH3_SOCK_FIXED_NBC_PROGRESS)
726145cd65SJed Brown     /* Deadlock in Ibarrier: http://trac.mpich.org/projects/mpich/ticket/1785 */
736145cd65SJed Brown     _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE;
74f6ced4a3SJed Brown #  else
756145cd65SJed Brown     _twosided_type = PETSC_BUILDTWOSIDED_IBARRIER;
76f6ced4a3SJed Brown #  endif
776145cd65SJed Brown #else
786145cd65SJed Brown     _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE;
796145cd65SJed Brown #endif
800298fd71SBarry Smith     ierr = PetscOptionsGetEnum(NULL,"-build_twosided",PetscBuildTwoSidedTypes,(PetscEnum*)&_twosided_type,NULL);CHKERRQ(ierr);
81f6ced4a3SJed Brown   }
82f6ced4a3SJed Brown   *twosided = _twosided_type;
83f6ced4a3SJed Brown   PetscFunctionReturn(0);
84f6ced4a3SJed Brown }
85f6ced4a3SJed Brown 
86f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER)
87f6ced4a3SJed Brown 
88f6ced4a3SJed Brown #undef __FUNCT__
896145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Ibarrier"
906145cd65SJed Brown static PetscErrorCode PetscCommBuildTwoSided_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
91f6ced4a3SJed Brown {
92f6ced4a3SJed Brown   PetscErrorCode ierr;
93f6ced4a3SJed Brown   PetscMPIInt    nrecvs,tag,unitbytes,done;
94f6ced4a3SJed Brown   PetscInt       i;
95f6ced4a3SJed Brown   char           *tdata;
96f6ced4a3SJed Brown   MPI_Request    *sendreqs,barrier;
970f453b92SJed Brown   PetscSegBuffer segrank,segdata;
98f6ced4a3SJed Brown 
99f6ced4a3SJed Brown   PetscFunctionBegin;
100f6ced4a3SJed Brown   ierr  = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
101f6ced4a3SJed Brown   ierr  = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr);
102f6ced4a3SJed Brown   tdata = (char*)todata;
103785e854fSJed Brown   ierr  = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr);
104f6ced4a3SJed Brown   for (i=0; i<nto; i++) {
105f6ced4a3SJed Brown     ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr);
106f6ced4a3SJed Brown   }
1070f453b92SJed Brown   ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr);
1080f453b92SJed Brown   ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr);
109f6ced4a3SJed Brown 
110f6ced4a3SJed Brown   nrecvs  = 0;
111f6ced4a3SJed Brown   barrier = MPI_REQUEST_NULL;
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     }
125f6ced4a3SJed Brown     if (barrier == MPI_REQUEST_NULL) {
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) {
130f6ced4a3SJed Brown         ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr);
131f6ced4a3SJed Brown         ierr = PetscFree(sendreqs);CHKERRQ(ierr);
132f6ced4a3SJed Brown       }
133f6ced4a3SJed Brown     } else {
134f6ced4a3SJed Brown       ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr);
135f6ced4a3SJed Brown     }
136f6ced4a3SJed Brown   }
137f6ced4a3SJed Brown   *nfrom = nrecvs;
138137cf7b6SJed Brown   ierr   = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr);
1390f453b92SJed Brown   ierr   = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr);
140137cf7b6SJed Brown   ierr   = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr);
1410f453b92SJed Brown   ierr   = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr);
142f6ced4a3SJed Brown   PetscFunctionReturn(0);
143f6ced4a3SJed Brown }
144f6ced4a3SJed Brown #endif
145f6ced4a3SJed Brown 
146f6ced4a3SJed Brown #undef __FUNCT__
1476145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Allreduce"
1486145cd65SJed Brown static PetscErrorCode PetscCommBuildTwoSided_Allreduce(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
149f6ced4a3SJed Brown {
150f6ced4a3SJed Brown   PetscErrorCode ierr;
151f6ced4a3SJed Brown   PetscMPIInt    size,*iflags,nrecvs,tag,unitbytes,*franks;
152f6ced4a3SJed Brown   PetscInt       i;
153f6ced4a3SJed Brown   char           *tdata,*fdata;
154f6ced4a3SJed Brown   MPI_Request    *reqs,*sendreqs;
155f6ced4a3SJed Brown   MPI_Status     *statuses;
156f6ced4a3SJed Brown 
157f6ced4a3SJed Brown   PetscFunctionBegin;
158f6ced4a3SJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
159*1795a4d1SJed Brown   ierr = PetscCalloc1(size,&iflags);CHKERRQ(ierr);
160f6ced4a3SJed Brown   for (i=0; i<nto; i++) iflags[toranks[i]] = 1;
1610298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,iflags,NULL,&nrecvs);CHKERRQ(ierr);
162f6ced4a3SJed Brown   ierr = PetscFree(iflags);CHKERRQ(ierr);
163f6ced4a3SJed Brown 
164f6ced4a3SJed Brown   ierr     = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
165f6ced4a3SJed Brown   ierr     = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr);
166f6ced4a3SJed Brown   ierr     = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr);
167f6ced4a3SJed Brown   tdata    = (char*)todata;
168dcca6d9dSJed Brown   ierr     = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr);
16906926cafSJed Brown   sendreqs = reqs + nrecvs;
170f6ced4a3SJed Brown   for (i=0; i<nrecvs; i++) {
171f6ced4a3SJed Brown     ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr);
172f6ced4a3SJed Brown   }
173f6ced4a3SJed Brown   for (i=0; i<nto; i++) {
174f6ced4a3SJed Brown     ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr);
175f6ced4a3SJed Brown   }
176f6ced4a3SJed Brown   ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr);
177785e854fSJed Brown   ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr);
178a297a907SKarl Rupp   for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE;
17906926cafSJed Brown   ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr);
180f6ced4a3SJed Brown 
181f6ced4a3SJed Brown   *nfrom            = nrecvs;
182f6ced4a3SJed Brown   *fromranks        = franks;
183f6ced4a3SJed Brown   *(void**)fromdata = fdata;
184f6ced4a3SJed Brown   PetscFunctionReturn(0);
185f6ced4a3SJed Brown }
186f6ced4a3SJed Brown 
187f6ced4a3SJed Brown #undef __FUNCT__
188f6ced4a3SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided"
189f6ced4a3SJed Brown /*@C
190f6ced4a3SJed Brown    PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths)
191f6ced4a3SJed Brown 
192f6ced4a3SJed Brown    Collective on MPI_Comm
193f6ced4a3SJed Brown 
194f6ced4a3SJed Brown    Input Arguments:
195f6ced4a3SJed Brown +  comm - communicator
196f6ced4a3SJed Brown .  count - number of entries to send/receive (must match on all ranks)
197f6ced4a3SJed Brown .  dtype - datatype to send/receive from each rank (must match on all ranks)
198f6ced4a3SJed Brown .  nto - number of ranks to send data to
199f6ced4a3SJed Brown .  toranks - ranks to send to (array of length nto)
200f6ced4a3SJed Brown -  todata - data to send to each rank (packed)
201f6ced4a3SJed Brown 
202f6ced4a3SJed Brown    Output Arguments:
203f6ced4a3SJed Brown +  nfrom - number of ranks receiving messages from
204f6ced4a3SJed Brown .  fromranks - ranks receiving messages from (length nfrom; caller should PetscFree())
205f6ced4a3SJed Brown -  fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree())
206f6ced4a3SJed Brown 
207f6ced4a3SJed Brown    Level: developer
208f6ced4a3SJed Brown 
209f6ced4a3SJed Brown    Notes:
210f6ced4a3SJed Brown    This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and
211f6ced4a3SJed Brown    PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data.
212f6ced4a3SJed Brown 
213f6ced4a3SJed Brown    Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not.
214f6ced4a3SJed Brown 
215f6ced4a3SJed Brown    References:
216f6ced4a3SJed Brown    The MPI_Ibarrier implementation uses the algorithm in
217f6ced4a3SJed Brown    Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010.
218f6ced4a3SJed Brown 
219f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths()
220f6ced4a3SJed Brown @*/
221f6ced4a3SJed Brown PetscErrorCode PetscCommBuildTwoSided(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata)
222f6ced4a3SJed Brown {
223f6ced4a3SJed Brown   PetscErrorCode         ierr;
2246145cd65SJed Brown   PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET;
225f6ced4a3SJed Brown 
226f6ced4a3SJed Brown   PetscFunctionBegin;
2276145cd65SJed Brown   ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr);
228f6ced4a3SJed Brown   switch (buildtype) {
2296145cd65SJed Brown   case PETSC_BUILDTWOSIDED_IBARRIER:
230f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER)
2316145cd65SJed Brown     ierr = PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr);
2326145cd65SJed Brown #else
2336145cd65SJed Brown     SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)");
234f6ced4a3SJed Brown #endif
2356145cd65SJed Brown     break;
2366145cd65SJed Brown   case PETSC_BUILDTWOSIDED_ALLREDUCE:
2376145cd65SJed Brown     ierr = PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr);
238f6ced4a3SJed Brown     break;
239f6ced4a3SJed Brown   default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication");
240f6ced4a3SJed Brown   }
241f6ced4a3SJed Brown   PetscFunctionReturn(0);
242f6ced4a3SJed Brown }
243