1f6ced4a3SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 2f6ced4a3SJed Brown 33b3561c8SJed Brown PetscLogEvent PETSC_BuildTwoSided; 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 89f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 90f6ced4a3SJed Brown 91f6ced4a3SJed Brown #undef __FUNCT__ 926145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Ibarrier" 936145cd65SJed 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) 94f6ced4a3SJed Brown { 95f6ced4a3SJed Brown PetscErrorCode ierr; 96f6ced4a3SJed Brown PetscMPIInt nrecvs,tag,unitbytes,done; 97f6ced4a3SJed Brown PetscInt i; 98f6ced4a3SJed Brown char *tdata; 99f6ced4a3SJed Brown MPI_Request *sendreqs,barrier; 1000f453b92SJed Brown PetscSegBuffer segrank,segdata; 101f6ced4a3SJed Brown 102f6ced4a3SJed Brown PetscFunctionBegin; 103*a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 104f6ced4a3SJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 105f6ced4a3SJed Brown tdata = (char*)todata; 106785e854fSJed Brown ierr = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr); 107f6ced4a3SJed Brown for (i=0; i<nto; i++) { 108f6ced4a3SJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 109f6ced4a3SJed Brown } 1100f453b92SJed Brown ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 1110f453b92SJed Brown ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 112f6ced4a3SJed Brown 113f6ced4a3SJed Brown nrecvs = 0; 114f6ced4a3SJed Brown barrier = MPI_REQUEST_NULL; 115f6ced4a3SJed Brown for (done=0; !done; ) { 116f6ced4a3SJed Brown PetscMPIInt flag; 117f6ced4a3SJed Brown MPI_Status status; 118f6ced4a3SJed Brown ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr); 119f6ced4a3SJed Brown if (flag) { /* incoming message */ 120f6ced4a3SJed Brown PetscMPIInt *recvrank; 121f6ced4a3SJed Brown void *buf; 122137cf7b6SJed Brown ierr = PetscSegBufferGet(segrank,1,&recvrank);CHKERRQ(ierr); 123137cf7b6SJed Brown ierr = PetscSegBufferGet(segdata,count,&buf);CHKERRQ(ierr); 124f6ced4a3SJed Brown *recvrank = status.MPI_SOURCE; 125f6ced4a3SJed Brown ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 126f6ced4a3SJed Brown nrecvs++; 127f6ced4a3SJed Brown } 128f6ced4a3SJed Brown if (barrier == MPI_REQUEST_NULL) { 1294dc2109aSBarry Smith PetscMPIInt sent,nsends; 1304dc2109aSBarry Smith ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr); 131f6ced4a3SJed Brown ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 132f6ced4a3SJed Brown if (sent) { 133f6ced4a3SJed Brown ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr); 134f6ced4a3SJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 135f6ced4a3SJed Brown } 136f6ced4a3SJed Brown } else { 137f6ced4a3SJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 138f6ced4a3SJed Brown } 139f6ced4a3SJed Brown } 140f6ced4a3SJed Brown *nfrom = nrecvs; 141137cf7b6SJed Brown ierr = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr); 1420f453b92SJed Brown ierr = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr); 143137cf7b6SJed Brown ierr = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr); 1440f453b92SJed Brown ierr = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr); 145*a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 146f6ced4a3SJed Brown PetscFunctionReturn(0); 147f6ced4a3SJed Brown } 148f6ced4a3SJed Brown #endif 149f6ced4a3SJed Brown 150f6ced4a3SJed Brown #undef __FUNCT__ 1516145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Allreduce" 1526145cd65SJed 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) 153f6ced4a3SJed Brown { 154f6ced4a3SJed Brown PetscErrorCode ierr; 155f6ced4a3SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,unitbytes,*franks; 156f6ced4a3SJed Brown PetscInt i; 157f6ced4a3SJed Brown char *tdata,*fdata; 158f6ced4a3SJed Brown MPI_Request *reqs,*sendreqs; 159f6ced4a3SJed Brown MPI_Status *statuses; 160f6ced4a3SJed Brown 161f6ced4a3SJed Brown PetscFunctionBegin; 162f6ced4a3SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1631795a4d1SJed Brown ierr = PetscCalloc1(size,&iflags);CHKERRQ(ierr); 164f6ced4a3SJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 1650298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,iflags,NULL,&nrecvs);CHKERRQ(ierr); 166f6ced4a3SJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 167f6ced4a3SJed Brown 168*a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 169f6ced4a3SJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 170f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 171f6ced4a3SJed Brown tdata = (char*)todata; 172dcca6d9dSJed Brown ierr = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr); 17306926cafSJed Brown sendreqs = reqs + nrecvs; 174f6ced4a3SJed Brown for (i=0; i<nrecvs; i++) { 175f6ced4a3SJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 176f6ced4a3SJed Brown } 177f6ced4a3SJed Brown for (i=0; i<nto; i++) { 178f6ced4a3SJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 179f6ced4a3SJed Brown } 180f6ced4a3SJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 181785e854fSJed Brown ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr); 182a297a907SKarl Rupp for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 18306926cafSJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 184*a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 185f6ced4a3SJed Brown 186f6ced4a3SJed Brown *nfrom = nrecvs; 187f6ced4a3SJed Brown *fromranks = franks; 188f6ced4a3SJed Brown *(void**)fromdata = fdata; 189f6ced4a3SJed Brown PetscFunctionReturn(0); 190f6ced4a3SJed Brown } 191f6ced4a3SJed Brown 1921654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 1931654bf6bSJed Brown #undef __FUNCT__ 1941654bf6bSJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_RedScatter" 1951654bf6bSJed Brown static PetscErrorCode PetscCommBuildTwoSided_RedScatter(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata) 1961654bf6bSJed Brown { 1971654bf6bSJed Brown PetscErrorCode ierr; 1981654bf6bSJed Brown PetscMPIInt size,*iflags,nrecvs,tag,unitbytes,*franks; 1991654bf6bSJed Brown PetscInt i; 2001654bf6bSJed Brown char *tdata,*fdata; 2011654bf6bSJed Brown MPI_Request *reqs,*sendreqs; 2021654bf6bSJed Brown MPI_Status *statuses; 2031654bf6bSJed Brown 2041654bf6bSJed Brown PetscFunctionBegin; 2051654bf6bSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2061654bf6bSJed Brown ierr = PetscMalloc1(size,&iflags);CHKERRQ(ierr); 2071654bf6bSJed Brown ierr = PetscMemzero(iflags,size*sizeof(*iflags));CHKERRQ(ierr); 2081654bf6bSJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 2091654bf6bSJed Brown ierr = MPI_Reduce_scatter_block(iflags,&nrecvs,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2101654bf6bSJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 2111654bf6bSJed Brown 212*a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 2131654bf6bSJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 2141654bf6bSJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 2151654bf6bSJed Brown tdata = (char*)todata; 2161654bf6bSJed Brown ierr = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr); 2171654bf6bSJed Brown sendreqs = reqs + nrecvs; 2181654bf6bSJed Brown for (i=0; i<nrecvs; i++) { 2191654bf6bSJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 2201654bf6bSJed Brown } 2211654bf6bSJed Brown for (i=0; i<nto; i++) { 2221654bf6bSJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 2231654bf6bSJed Brown } 2241654bf6bSJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 2251654bf6bSJed Brown ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr); 2261654bf6bSJed Brown for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 2271654bf6bSJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 228*a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 2291654bf6bSJed Brown 2301654bf6bSJed Brown *nfrom = nrecvs; 2311654bf6bSJed Brown *fromranks = franks; 2321654bf6bSJed Brown *(void**)fromdata = fdata; 2331654bf6bSJed Brown PetscFunctionReturn(0); 2341654bf6bSJed Brown } 2351654bf6bSJed Brown #endif 2361654bf6bSJed Brown 237f6ced4a3SJed Brown #undef __FUNCT__ 238f6ced4a3SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided" 239f6ced4a3SJed Brown /*@C 240f6ced4a3SJed Brown PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths) 241f6ced4a3SJed Brown 242f6ced4a3SJed Brown Collective on MPI_Comm 243f6ced4a3SJed Brown 244f6ced4a3SJed Brown Input Arguments: 245f6ced4a3SJed Brown + comm - communicator 246f6ced4a3SJed Brown . count - number of entries to send/receive (must match on all ranks) 247f6ced4a3SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 248f6ced4a3SJed Brown . nto - number of ranks to send data to 249f6ced4a3SJed Brown . toranks - ranks to send to (array of length nto) 250f6ced4a3SJed Brown - todata - data to send to each rank (packed) 251f6ced4a3SJed Brown 252f6ced4a3SJed Brown Output Arguments: 253f6ced4a3SJed Brown + nfrom - number of ranks receiving messages from 254f6ced4a3SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 255f6ced4a3SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 256f6ced4a3SJed Brown 257f6ced4a3SJed Brown Level: developer 258f6ced4a3SJed Brown 2591654bf6bSJed Brown Options Database Keys: 2601654bf6bSJed Brown . -build_twosided <allreduce|ibarrier|redscatter> - algorithm to set up two-sided communication 2611654bf6bSJed Brown 262f6ced4a3SJed Brown Notes: 263f6ced4a3SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 264f6ced4a3SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data. 265f6ced4a3SJed Brown 266f6ced4a3SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 267f6ced4a3SJed Brown 268f6ced4a3SJed Brown References: 269f6ced4a3SJed Brown The MPI_Ibarrier implementation uses the algorithm in 270f6ced4a3SJed Brown Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010. 271f6ced4a3SJed Brown 272f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 273f6ced4a3SJed Brown @*/ 274f6ced4a3SJed 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) 275f6ced4a3SJed Brown { 276f6ced4a3SJed Brown PetscErrorCode ierr; 2776145cd65SJed Brown PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET; 278f6ced4a3SJed Brown 279f6ced4a3SJed Brown PetscFunctionBegin; 2803b3561c8SJed Brown ierr = PetscSysInitializePackage();CHKERRQ(ierr); 2813b3561c8SJed Brown ierr = PetscLogEventBegin(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr); 2826145cd65SJed Brown ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr); 283f6ced4a3SJed Brown switch (buildtype) { 2846145cd65SJed Brown case PETSC_BUILDTWOSIDED_IBARRIER: 285f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 2866145cd65SJed Brown ierr = PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 2876145cd65SJed Brown #else 2886145cd65SJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)"); 289f6ced4a3SJed Brown #endif 2906145cd65SJed Brown break; 2916145cd65SJed Brown case PETSC_BUILDTWOSIDED_ALLREDUCE: 2926145cd65SJed Brown ierr = PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 293f6ced4a3SJed Brown break; 2941654bf6bSJed Brown case PETSC_BUILDTWOSIDED_REDSCATTER: 2951654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 2961654bf6bSJed Brown ierr = PetscCommBuildTwoSided_RedScatter(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 2971654bf6bSJed Brown #else 2981654bf6bSJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Reduce_scatter_block (part of MPI-2.2)"); 2991654bf6bSJed Brown #endif 3001654bf6bSJed Brown break; 301f6ced4a3SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 302f6ced4a3SJed Brown } 3033b3561c8SJed Brown ierr = PetscLogEventEnd(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr); 304f6ced4a3SJed Brown PetscFunctionReturn(0); 305f6ced4a3SJed Brown } 306