1f6ced4a3SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 2f6ced4a3SJed Brown 3*27a32ea5SJed Brown PetscLogEvent PETSC_BuildTwoSided,PETSC_BuildTwoSidedF; 43b3561c8SJed Brown 56145cd65SJed Brown const char *const PetscBuildTwoSidedTypes[] = { 6f6ced4a3SJed Brown "ALLREDUCE", 76145cd65SJed Brown "IBARRIER", 81654bf6bSJed Brown "REDSCATTER", 96145cd65SJed Brown "PetscBuildTwoSidedType", 106145cd65SJed Brown "PETSC_BUILDTWOSIDED_", 11f6ced4a3SJed Brown 0 12f6ced4a3SJed Brown }; 13f6ced4a3SJed Brown 146145cd65SJed Brown static PetscBuildTwoSidedType _twosided_type = PETSC_BUILDTWOSIDED_NOTSET; 15f6ced4a3SJed Brown 16f6ced4a3SJed Brown #undef __FUNCT__ 176145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedSetType" 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; 436145cd65SJed Brown ierr = MPI_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 #undef __FUNCT__ 526145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedGetType" 536145cd65SJed Brown /*@ 546145cd65SJed Brown PetscCommBuildTwoSidedGetType - set algorithm to use when building two-sided communication 556145cd65SJed Brown 566145cd65SJed Brown Logically Collective 576145cd65SJed Brown 586145cd65SJed Brown Output Arguments: 596145cd65SJed Brown + comm - communicator on which to query algorithm 606145cd65SJed Brown - twosided - algorithm to use for PetscCommBuildTwoSided() 616145cd65SJed Brown 626145cd65SJed Brown Level: developer 636145cd65SJed Brown 646145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedSetType() 656145cd65SJed Brown @*/ 666145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedGetType(MPI_Comm comm,PetscBuildTwoSidedType *twosided) 67f6ced4a3SJed Brown { 68f6ced4a3SJed Brown PetscErrorCode ierr; 69f6ced4a3SJed Brown 70f6ced4a3SJed Brown PetscFunctionBegin; 716145cd65SJed Brown *twosided = PETSC_BUILDTWOSIDED_NOTSET; 726145cd65SJed Brown if (_twosided_type == PETSC_BUILDTWOSIDED_NOTSET) { 73f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 746145cd65SJed Brown # if defined(PETSC_HAVE_MPICH_CH3_SOCK) && !defined(PETSC_HAVE_MPICH_CH3_SOCK_FIXED_NBC_PROGRESS) 756145cd65SJed Brown /* Deadlock in Ibarrier: http://trac.mpich.org/projects/mpich/ticket/1785 */ 766145cd65SJed Brown _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; 77f6ced4a3SJed Brown # else 786145cd65SJed Brown _twosided_type = PETSC_BUILDTWOSIDED_IBARRIER; 79f6ced4a3SJed Brown # endif 806145cd65SJed Brown #else 816145cd65SJed Brown _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; 826145cd65SJed Brown #endif 830298fd71SBarry Smith ierr = PetscOptionsGetEnum(NULL,"-build_twosided",PetscBuildTwoSidedTypes,(PetscEnum*)&_twosided_type,NULL);CHKERRQ(ierr); 84f6ced4a3SJed Brown } 85f6ced4a3SJed Brown *twosided = _twosided_type; 86f6ced4a3SJed Brown PetscFunctionReturn(0); 87f6ced4a3SJed Brown } 88f6ced4a3SJed Brown 89b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 90f6ced4a3SJed Brown 91f6ced4a3SJed Brown #undef __FUNCT__ 926145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Ibarrier" 93cf4b5b4fSJed 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) 94f6ced4a3SJed Brown { 95f6ced4a3SJed Brown PetscErrorCode ierr; 96c3a59b84SJed Brown PetscMPIInt nrecvs,tag,done,i; 97c3a59b84SJed Brown MPI_Aint lb,unitbytes; 98f6ced4a3SJed Brown char *tdata; 99f6ced4a3SJed Brown MPI_Request *sendreqs,barrier; 1000f453b92SJed Brown PetscSegBuffer segrank,segdata; 101f6ced4a3SJed Brown 102f6ced4a3SJed Brown PetscFunctionBegin; 103a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 104c3a59b84SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 105c3a59b84SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 106f6ced4a3SJed Brown tdata = (char*)todata; 107785e854fSJed Brown ierr = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr); 108f6ced4a3SJed Brown for (i=0; i<nto; i++) { 109f6ced4a3SJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 110f6ced4a3SJed Brown } 1110f453b92SJed Brown ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 1120f453b92SJed Brown ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 113f6ced4a3SJed Brown 114f6ced4a3SJed Brown nrecvs = 0; 115f6ced4a3SJed Brown barrier = MPI_REQUEST_NULL; 116f6ced4a3SJed Brown for (done=0; !done; ) { 117f6ced4a3SJed Brown PetscMPIInt flag; 118f6ced4a3SJed Brown MPI_Status status; 119f6ced4a3SJed Brown ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr); 120f6ced4a3SJed Brown if (flag) { /* incoming message */ 121f6ced4a3SJed Brown PetscMPIInt *recvrank; 122f6ced4a3SJed Brown void *buf; 123137cf7b6SJed Brown ierr = PetscSegBufferGet(segrank,1,&recvrank);CHKERRQ(ierr); 124137cf7b6SJed Brown ierr = PetscSegBufferGet(segdata,count,&buf);CHKERRQ(ierr); 125f6ced4a3SJed Brown *recvrank = status.MPI_SOURCE; 126f6ced4a3SJed Brown ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 127f6ced4a3SJed Brown nrecvs++; 128f6ced4a3SJed Brown } 129f6ced4a3SJed Brown if (barrier == MPI_REQUEST_NULL) { 1304dc2109aSBarry Smith PetscMPIInt sent,nsends; 1314dc2109aSBarry Smith ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr); 132f6ced4a3SJed Brown ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 133f6ced4a3SJed Brown if (sent) { 134b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 135f6ced4a3SJed Brown ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr); 136b5a01e9fSJed Brown #elif defined(PETSC_HAVE_MPIX_IBARRIER) 137b5a01e9fSJed Brown ierr = MPIX_Ibarrier(comm,&barrier);CHKERRQ(ierr); 138b5a01e9fSJed Brown #endif 139f6ced4a3SJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 140f6ced4a3SJed Brown } 141f6ced4a3SJed Brown } else { 142f6ced4a3SJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 143f6ced4a3SJed Brown } 144f6ced4a3SJed Brown } 145f6ced4a3SJed Brown *nfrom = nrecvs; 146137cf7b6SJed Brown ierr = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr); 1470f453b92SJed Brown ierr = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr); 148137cf7b6SJed Brown ierr = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr); 1490f453b92SJed Brown ierr = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr); 150a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 151f6ced4a3SJed Brown PetscFunctionReturn(0); 152f6ced4a3SJed Brown } 153f6ced4a3SJed Brown #endif 154f6ced4a3SJed Brown 155f6ced4a3SJed Brown #undef __FUNCT__ 1566145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Allreduce" 157cf4b5b4fSJed 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) 158f6ced4a3SJed Brown { 159f6ced4a3SJed Brown PetscErrorCode ierr; 160c3a59b84SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,*franks,i; 161c3a59b84SJed Brown MPI_Aint lb,unitbytes; 162f6ced4a3SJed Brown char *tdata,*fdata; 163f6ced4a3SJed Brown MPI_Request *reqs,*sendreqs; 164f6ced4a3SJed Brown MPI_Status *statuses; 165f6ced4a3SJed Brown 166f6ced4a3SJed Brown PetscFunctionBegin; 167f6ced4a3SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1681795a4d1SJed Brown ierr = PetscCalloc1(size,&iflags);CHKERRQ(ierr); 169f6ced4a3SJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 1700298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,iflags,NULL,&nrecvs);CHKERRQ(ierr); 171f6ced4a3SJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 172f6ced4a3SJed Brown 173a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 174c3a59b84SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 175c3a59b84SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 176f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 177f6ced4a3SJed Brown tdata = (char*)todata; 178dcca6d9dSJed Brown ierr = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr); 17906926cafSJed Brown sendreqs = reqs + nrecvs; 180f6ced4a3SJed Brown for (i=0; i<nrecvs; i++) { 181f6ced4a3SJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 182f6ced4a3SJed Brown } 183f6ced4a3SJed Brown for (i=0; i<nto; i++) { 184f6ced4a3SJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 185f6ced4a3SJed Brown } 186f6ced4a3SJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 187785e854fSJed Brown ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr); 188a297a907SKarl Rupp for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 18906926cafSJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 190a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 191f6ced4a3SJed Brown 192f6ced4a3SJed Brown *nfrom = nrecvs; 193f6ced4a3SJed Brown *fromranks = franks; 194f6ced4a3SJed Brown *(void**)fromdata = fdata; 195f6ced4a3SJed Brown PetscFunctionReturn(0); 196f6ced4a3SJed Brown } 197f6ced4a3SJed Brown 1981654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 1991654bf6bSJed Brown #undef __FUNCT__ 2001654bf6bSJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_RedScatter" 201cf4b5b4fSJed 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) 2021654bf6bSJed Brown { 2031654bf6bSJed Brown PetscErrorCode ierr; 204c3a59b84SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,*franks,i; 205c3a59b84SJed Brown MPI_Aint lb,unitbytes; 2061654bf6bSJed Brown char *tdata,*fdata; 2071654bf6bSJed Brown MPI_Request *reqs,*sendreqs; 2081654bf6bSJed Brown MPI_Status *statuses; 2091654bf6bSJed Brown 2101654bf6bSJed Brown PetscFunctionBegin; 2111654bf6bSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2121654bf6bSJed Brown ierr = PetscMalloc1(size,&iflags);CHKERRQ(ierr); 2131654bf6bSJed Brown ierr = PetscMemzero(iflags,size*sizeof(*iflags));CHKERRQ(ierr); 2141654bf6bSJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 2151654bf6bSJed Brown ierr = MPI_Reduce_scatter_block(iflags,&nrecvs,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2161654bf6bSJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 2171654bf6bSJed Brown 218a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 219c3a59b84SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 220c3a59b84SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 2211654bf6bSJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 2221654bf6bSJed Brown tdata = (char*)todata; 2231654bf6bSJed Brown ierr = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr); 2241654bf6bSJed Brown sendreqs = reqs + nrecvs; 2251654bf6bSJed Brown for (i=0; i<nrecvs; i++) { 2261654bf6bSJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 2271654bf6bSJed Brown } 2281654bf6bSJed Brown for (i=0; i<nto; i++) { 2291654bf6bSJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 2301654bf6bSJed Brown } 2311654bf6bSJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 2321654bf6bSJed Brown ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr); 2331654bf6bSJed Brown for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 2341654bf6bSJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 235a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 2361654bf6bSJed Brown 2371654bf6bSJed Brown *nfrom = nrecvs; 2381654bf6bSJed Brown *fromranks = franks; 2391654bf6bSJed Brown *(void**)fromdata = fdata; 2401654bf6bSJed Brown PetscFunctionReturn(0); 2411654bf6bSJed Brown } 2421654bf6bSJed Brown #endif 2431654bf6bSJed Brown 244f6ced4a3SJed Brown #undef __FUNCT__ 245f6ced4a3SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided" 246f6ced4a3SJed Brown /*@C 247f6ced4a3SJed Brown PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths) 248f6ced4a3SJed Brown 249f6ced4a3SJed Brown Collective on MPI_Comm 250f6ced4a3SJed Brown 251f6ced4a3SJed Brown Input Arguments: 252f6ced4a3SJed Brown + comm - communicator 253f6ced4a3SJed Brown . count - number of entries to send/receive (must match on all ranks) 254f6ced4a3SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 255f6ced4a3SJed Brown . nto - number of ranks to send data to 256f6ced4a3SJed Brown . toranks - ranks to send to (array of length nto) 257f6ced4a3SJed Brown - todata - data to send to each rank (packed) 258f6ced4a3SJed Brown 259f6ced4a3SJed Brown Output Arguments: 260f6ced4a3SJed Brown + nfrom - number of ranks receiving messages from 261f6ced4a3SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 262f6ced4a3SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 263f6ced4a3SJed Brown 264f6ced4a3SJed Brown Level: developer 265f6ced4a3SJed Brown 2661654bf6bSJed Brown Options Database Keys: 2671654bf6bSJed Brown . -build_twosided <allreduce|ibarrier|redscatter> - algorithm to set up two-sided communication 2681654bf6bSJed Brown 269f6ced4a3SJed Brown Notes: 270f6ced4a3SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 271f6ced4a3SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data. 272f6ced4a3SJed Brown 273f6ced4a3SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 274f6ced4a3SJed Brown 275f6ced4a3SJed Brown References: 276f6ced4a3SJed Brown The MPI_Ibarrier implementation uses the algorithm in 277f6ced4a3SJed Brown Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010. 278f6ced4a3SJed Brown 279f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 280f6ced4a3SJed Brown @*/ 281cf4b5b4fSJed 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) 282f6ced4a3SJed Brown { 283f6ced4a3SJed Brown PetscErrorCode ierr; 2846145cd65SJed Brown PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET; 285f6ced4a3SJed Brown 286f6ced4a3SJed Brown PetscFunctionBegin; 2873b3561c8SJed Brown ierr = PetscSysInitializePackage();CHKERRQ(ierr); 2883b3561c8SJed Brown ierr = PetscLogEventBegin(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr); 2896145cd65SJed Brown ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr); 290f6ced4a3SJed Brown switch (buildtype) { 2916145cd65SJed Brown case PETSC_BUILDTWOSIDED_IBARRIER: 292b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 2936145cd65SJed Brown ierr = PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 2946145cd65SJed Brown #else 2956145cd65SJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)"); 296f6ced4a3SJed Brown #endif 2976145cd65SJed Brown break; 2986145cd65SJed Brown case PETSC_BUILDTWOSIDED_ALLREDUCE: 2996145cd65SJed Brown ierr = PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 300f6ced4a3SJed Brown break; 3011654bf6bSJed Brown case PETSC_BUILDTWOSIDED_REDSCATTER: 3021654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 3031654bf6bSJed Brown ierr = PetscCommBuildTwoSided_RedScatter(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 3041654bf6bSJed Brown #else 3051654bf6bSJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Reduce_scatter_block (part of MPI-2.2)"); 3061654bf6bSJed Brown #endif 3071654bf6bSJed Brown break; 308f6ced4a3SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 309f6ced4a3SJed Brown } 3103b3561c8SJed Brown ierr = PetscLogEventEnd(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr); 311f6ced4a3SJed Brown PetscFunctionReturn(0); 312f6ced4a3SJed Brown } 313d815da10SJed Brown 314d815da10SJed Brown #undef __FUNCT__ 3156afbe7ddSJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedFReq_Reference" 3166afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Reference(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata, 3176afbe7ddSJed Brown PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs, 318d815da10SJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 319d815da10SJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 320d815da10SJed Brown { 321d815da10SJed Brown PetscErrorCode ierr; 322d815da10SJed Brown PetscMPIInt i,*tag; 323d815da10SJed Brown MPI_Aint lb,unitbytes; 324d815da10SJed Brown MPI_Request *sendreq,*recvreq; 325d815da10SJed Brown 326d815da10SJed Brown PetscFunctionBegin; 327d815da10SJed Brown ierr = PetscMalloc1(ntags,&tag);CHKERRQ(ierr); 328d815da10SJed Brown if (ntags > 0) { 329d815da10SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag[0]);CHKERRQ(ierr); 330d815da10SJed Brown } 331d815da10SJed Brown for (i=1; i<ntags; i++) { 332d815da10SJed Brown ierr = PetscCommGetNewTag(comm,&tag[i]);CHKERRQ(ierr); 333d815da10SJed Brown } 334d815da10SJed Brown 335d815da10SJed Brown /* Perform complete initial rendezvous */ 336d815da10SJed Brown ierr = PetscCommBuildTwoSided(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 337d815da10SJed Brown 3386afbe7ddSJed Brown ierr = PetscMalloc1(nto*ntags,&sendreq);CHKERRQ(ierr); 3396afbe7ddSJed Brown ierr = PetscMalloc1(*nfrom*ntags,&recvreq);CHKERRQ(ierr); 340d815da10SJed Brown 341d815da10SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 342d815da10SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 343d815da10SJed Brown for (i=0; i<nto; i++) { 34489157a57SJed Brown PetscMPIInt k; 34589157a57SJed Brown for (k=0; k<ntags; k++) sendreq[i*ntags+k] = MPI_REQUEST_NULL; 346d815da10SJed Brown ierr = (*send)(comm,tag,i,toranks[i],((char*)todata)+count*unitbytes*i,sendreq+i*ntags,ctx);CHKERRQ(ierr); 347d815da10SJed Brown } 348d815da10SJed Brown for (i=0; i<*nfrom; i++) { 349d815da10SJed Brown void *header = (*(char**)fromdata) + count*unitbytes*i; 35089157a57SJed Brown PetscMPIInt k; 35189157a57SJed Brown for (k=0; k<ntags; k++) recvreq[i*ntags+k] = MPI_REQUEST_NULL; 352d815da10SJed Brown ierr = (*recv)(comm,tag,(*fromranks)[i],header,recvreq+i*ntags,ctx);CHKERRQ(ierr); 353d815da10SJed Brown } 354d815da10SJed Brown ierr = PetscFree(tag);CHKERRQ(ierr); 355d815da10SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 3566afbe7ddSJed Brown *toreqs = sendreq; 3576afbe7ddSJed Brown *fromreqs = recvreq; 358d815da10SJed Brown PetscFunctionReturn(0); 359d815da10SJed Brown } 360d815da10SJed Brown 361b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 362afb254cdSJed Brown 363afb254cdSJed Brown #undef __FUNCT__ 3646afbe7ddSJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedFReq_Ibarrier" 3656afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata, 3666afbe7ddSJed Brown PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs, 367afb254cdSJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 368afb254cdSJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 369afb254cdSJed Brown { 370afb254cdSJed Brown PetscErrorCode ierr; 371afb254cdSJed Brown PetscMPIInt nrecvs,tag,*tags,done,i; 372afb254cdSJed Brown MPI_Aint lb,unitbytes; 373afb254cdSJed Brown char *tdata; 374afb254cdSJed Brown MPI_Request *sendreqs,*usendreqs,*req,barrier; 375afb254cdSJed Brown PetscSegBuffer segrank,segdata,segreq; 376afb254cdSJed Brown 377afb254cdSJed Brown PetscFunctionBegin; 378afb254cdSJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 379afb254cdSJed Brown ierr = PetscMalloc1(ntags,&tags);CHKERRQ(ierr); 380afb254cdSJed Brown for (i=0; i<ntags; i++) { 381afb254cdSJed Brown ierr = PetscCommGetNewTag(comm,&tags[i]);CHKERRQ(ierr); 382afb254cdSJed Brown } 383afb254cdSJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 384afb254cdSJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 385afb254cdSJed Brown tdata = (char*)todata; 3866afbe7ddSJed Brown ierr = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr); 3876afbe7ddSJed Brown ierr = PetscMalloc1(nto*ntags,&usendreqs);CHKERRQ(ierr); 388afb254cdSJed Brown /* Post synchronous sends */ 389afb254cdSJed Brown for (i=0; i<nto; i++) { 390afb254cdSJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 391afb254cdSJed Brown } 392afb254cdSJed Brown /* Post actual payloads. These are typically larger messages. Hopefully sending these later does not slow down the 393afb254cdSJed Brown * synchronous messages above. */ 394afb254cdSJed Brown for (i=0; i<nto; i++) { 39589157a57SJed Brown PetscMPIInt k; 39689157a57SJed Brown for (k=0; k<ntags; k++) usendreqs[i*ntags+k] = MPI_REQUEST_NULL; 397afb254cdSJed Brown ierr = (*send)(comm,tags,i,toranks[i],tdata+count*unitbytes*i,usendreqs+i*ntags,ctx);CHKERRQ(ierr); 398afb254cdSJed Brown } 399afb254cdSJed Brown 400afb254cdSJed Brown ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 401afb254cdSJed Brown ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 402afb254cdSJed Brown ierr = PetscSegBufferCreate(sizeof(MPI_Request),4,&segreq);CHKERRQ(ierr); 403afb254cdSJed Brown 404afb254cdSJed Brown nrecvs = 0; 405afb254cdSJed Brown barrier = MPI_REQUEST_NULL; 406afb254cdSJed Brown for (done=0; !done; ) { 407afb254cdSJed Brown PetscMPIInt flag; 408afb254cdSJed Brown MPI_Status status; 409afb254cdSJed Brown ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr); 410afb254cdSJed Brown if (flag) { /* incoming message */ 41189157a57SJed Brown PetscMPIInt *recvrank,k; 412afb254cdSJed Brown void *buf; 413afb254cdSJed Brown ierr = PetscSegBufferGet(segrank,1,&recvrank);CHKERRQ(ierr); 414afb254cdSJed Brown ierr = PetscSegBufferGet(segdata,count,&buf);CHKERRQ(ierr); 415afb254cdSJed Brown *recvrank = status.MPI_SOURCE; 416afb254cdSJed Brown ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 417afb254cdSJed Brown ierr = PetscSegBufferGet(segreq,ntags,&req);CHKERRQ(ierr); 41889157a57SJed Brown for (k=0; k<ntags; k++) req[k] = MPI_REQUEST_NULL; 419afb254cdSJed Brown ierr = (*recv)(comm,tags,status.MPI_SOURCE,buf,req,ctx);CHKERRQ(ierr); 420afb254cdSJed Brown nrecvs++; 421afb254cdSJed Brown } 422afb254cdSJed Brown if (barrier == MPI_REQUEST_NULL) { 423afb254cdSJed Brown PetscMPIInt sent,nsends; 424afb254cdSJed Brown ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr); 425afb254cdSJed Brown ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 426afb254cdSJed Brown if (sent) { 427b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 428afb254cdSJed Brown ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr); 429b5a01e9fSJed Brown #elif defined(PETSC_HAVE_MPIX_IBARRIER) 430b5a01e9fSJed Brown ierr = MPIX_Ibarrier(comm,&barrier);CHKERRQ(ierr); 431b5a01e9fSJed Brown #endif 432afb254cdSJed Brown } 433afb254cdSJed Brown } else { 434afb254cdSJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 435afb254cdSJed Brown } 436afb254cdSJed Brown } 437afb254cdSJed Brown *nfrom = nrecvs; 438afb254cdSJed Brown ierr = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr); 439afb254cdSJed Brown ierr = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr); 440afb254cdSJed Brown ierr = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr); 441afb254cdSJed Brown ierr = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr); 4426afbe7ddSJed Brown *toreqs = usendreqs; 4436afbe7ddSJed Brown ierr = PetscSegBufferExtractAlloc(segreq,fromreqs);CHKERRQ(ierr); 444afb254cdSJed Brown ierr = PetscSegBufferDestroy(&segreq);CHKERRQ(ierr); 4456afbe7ddSJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 446afb254cdSJed Brown ierr = PetscFree(tags);CHKERRQ(ierr); 447afb254cdSJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 448afb254cdSJed Brown PetscFunctionReturn(0); 449afb254cdSJed Brown } 450afb254cdSJed Brown #endif 451afb254cdSJed Brown 452d815da10SJed Brown #undef __FUNCT__ 453d815da10SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedF" 454d815da10SJed Brown /*@C 455d815da10SJed Brown PetscCommBuildTwoSidedF - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous 456d815da10SJed Brown 457d815da10SJed Brown Collective on MPI_Comm 458d815da10SJed Brown 459d815da10SJed Brown Input Arguments: 460d815da10SJed Brown + comm - communicator 461d815da10SJed Brown . count - number of entries to send/receive in initial rendezvous (must match on all ranks) 462d815da10SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 463d815da10SJed Brown . nto - number of ranks to send data to 464d815da10SJed Brown . toranks - ranks to send to (array of length nto) 465d815da10SJed Brown . todata - data to send to each rank (packed) 466d815da10SJed Brown . ntags - number of tags needed by send/recv callbacks 467d815da10SJed Brown . send - callback invoked on sending process when ready to send primary payload 468d815da10SJed Brown . recv - callback invoked on receiving process after delivery of rendezvous message 469d815da10SJed Brown - ctx - context for callbacks 470d815da10SJed Brown 471d815da10SJed Brown Output Arguments: 472d815da10SJed Brown + nfrom - number of ranks receiving messages from 473d815da10SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 474d815da10SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 475d815da10SJed Brown 476d815da10SJed Brown Level: developer 477d815da10SJed Brown 478d815da10SJed Brown Notes: 479d815da10SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 480d815da10SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data. 481d815da10SJed Brown 482d815da10SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 483d815da10SJed Brown 484d815da10SJed Brown References: 485d815da10SJed Brown The MPI_Ibarrier implementation uses the algorithm in 486d815da10SJed Brown Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010. 487d815da10SJed Brown 4886afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedFReq(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 489d815da10SJed Brown @*/ 490d815da10SJed 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, 491d815da10SJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 492d815da10SJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 493d815da10SJed Brown { 4946afbe7ddSJed Brown PetscErrorCode ierr; 4956afbe7ddSJed Brown MPI_Request *toreqs,*fromreqs; 4966afbe7ddSJed Brown 4976afbe7ddSJed Brown PetscFunctionBegin; 4986afbe7ddSJed Brown ierr = PetscCommBuildTwoSidedFReq(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,&toreqs,&fromreqs,send,recv,ctx);CHKERRQ(ierr); 4996afbe7ddSJed Brown ierr = MPI_Waitall(nto*ntags,toreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 5006afbe7ddSJed Brown ierr = MPI_Waitall(*nfrom*ntags,fromreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 5016afbe7ddSJed Brown ierr = PetscFree(toreqs);CHKERRQ(ierr); 5026afbe7ddSJed Brown ierr = PetscFree(fromreqs);CHKERRQ(ierr); 5036afbe7ddSJed Brown PetscFunctionReturn(0); 5046afbe7ddSJed Brown } 5056afbe7ddSJed Brown 5066afbe7ddSJed Brown #undef __FUNCT__ 5076afbe7ddSJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedFReq" 5086afbe7ddSJed Brown /*@C 5096afbe7ddSJed Brown PetscCommBuildTwoSidedFReq - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous, returns requests 5106afbe7ddSJed Brown 5116afbe7ddSJed Brown Collective on MPI_Comm 5126afbe7ddSJed Brown 5136afbe7ddSJed Brown Input Arguments: 5146afbe7ddSJed Brown + comm - communicator 5156afbe7ddSJed Brown . count - number of entries to send/receive in initial rendezvous (must match on all ranks) 5166afbe7ddSJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 5176afbe7ddSJed Brown . nto - number of ranks to send data to 5186afbe7ddSJed Brown . toranks - ranks to send to (array of length nto) 5196afbe7ddSJed Brown . todata - data to send to each rank (packed) 5206afbe7ddSJed Brown . ntags - number of tags needed by send/recv callbacks 5216afbe7ddSJed Brown . send - callback invoked on sending process when ready to send primary payload 5226afbe7ddSJed Brown . recv - callback invoked on receiving process after delivery of rendezvous message 5236afbe7ddSJed Brown - ctx - context for callbacks 5246afbe7ddSJed Brown 5256afbe7ddSJed Brown Output Arguments: 5266afbe7ddSJed Brown + nfrom - number of ranks receiving messages from 5276afbe7ddSJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 5286afbe7ddSJed Brown . fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 5296afbe7ddSJed Brown . toreqs - array of nto*ntags sender requests (caller must wait on these, then PetscFree()) 5306afbe7ddSJed Brown - fromreqs - array of nfrom*ntags receiver requests (caller must wait on these, then PetscFree()) 5316afbe7ddSJed Brown 5326afbe7ddSJed Brown Level: developer 5336afbe7ddSJed Brown 5346afbe7ddSJed Brown Notes: 5356afbe7ddSJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 5366afbe7ddSJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data. 5376afbe7ddSJed Brown 5386afbe7ddSJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 5396afbe7ddSJed Brown 5406afbe7ddSJed Brown References: 5416afbe7ddSJed Brown The MPI_Ibarrier implementation uses the algorithm in 5426afbe7ddSJed Brown Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010. 5436afbe7ddSJed Brown 5446afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedF(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 5456afbe7ddSJed Brown @*/ 5466afbe7ddSJed Brown PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata, 5476afbe7ddSJed Brown PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs, 5486afbe7ddSJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 5496afbe7ddSJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 5506afbe7ddSJed Brown { 5516afbe7ddSJed Brown PetscErrorCode ierr,(*f)(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*, 5526afbe7ddSJed Brown PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,MPI_Request**,MPI_Request**, 553d815da10SJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 554d815da10SJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx); 555d815da10SJed Brown PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET; 556d815da10SJed Brown 557d815da10SJed Brown PetscFunctionBegin; 558d815da10SJed Brown ierr = PetscSysInitializePackage();CHKERRQ(ierr); 559*27a32ea5SJed Brown ierr = PetscLogEventBegin(PETSC_BuildTwoSidedF,0,0,0,0);CHKERRQ(ierr); 560d815da10SJed Brown ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr); 561d815da10SJed Brown switch (buildtype) { 562d815da10SJed Brown case PETSC_BUILDTWOSIDED_IBARRIER: 563b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 5646afbe7ddSJed Brown f = PetscCommBuildTwoSidedFReq_Ibarrier; 565afb254cdSJed Brown #else 566afb254cdSJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)"); 567afb254cdSJed Brown #endif 568afb254cdSJed Brown break; 569d815da10SJed Brown case PETSC_BUILDTWOSIDED_ALLREDUCE: 570d815da10SJed Brown case PETSC_BUILDTWOSIDED_REDSCATTER: 5716afbe7ddSJed Brown f = PetscCommBuildTwoSidedFReq_Reference; 572d815da10SJed Brown break; 573d815da10SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 574d815da10SJed Brown } 5756afbe7ddSJed Brown ierr = (*f)(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,toreqs,fromreqs,send,recv,ctx);CHKERRQ(ierr); 576*27a32ea5SJed Brown ierr = PetscLogEventEnd(PETSC_BuildTwoSidedF,0,0,0,0);CHKERRQ(ierr); 577d815da10SJed Brown PetscFunctionReturn(0); 578d815da10SJed Brown } 579