1*f6ced4a3SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 2*f6ced4a3SJed Brown #include <stddef.h> 3*f6ced4a3SJed Brown 4*f6ced4a3SJed Brown typedef enum {BUILDTWOSIDED_NOTSET = -1, 5*f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 6*f6ced4a3SJed Brown BUILDTWOSIDED_IBARRIER, 7*f6ced4a3SJed Brown #endif 8*f6ced4a3SJed Brown BUILDTWOSIDED_ALLREDUCE} BuildTwoSidedType; 9*f6ced4a3SJed Brown 10*f6ced4a3SJed Brown static const char *const BuildTwoSidedTypes[] = { 11*f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 12*f6ced4a3SJed Brown "IBARRIER", 13*f6ced4a3SJed Brown #endif 14*f6ced4a3SJed Brown "ALLREDUCE", 15*f6ced4a3SJed Brown "BuildTwoSidedType", 16*f6ced4a3SJed Brown "BUILDTWOSIDED_", 17*f6ced4a3SJed Brown 0 18*f6ced4a3SJed Brown }; 19*f6ced4a3SJed Brown 20*f6ced4a3SJed Brown static BuildTwoSidedType _twosided_type = BUILDTWOSIDED_NOTSET; 21*f6ced4a3SJed Brown 22*f6ced4a3SJed Brown #undef __FUNCT__ 23*f6ced4a3SJed Brown #define __FUNCT__ "PetscBuildTwoSidedGetType" 24*f6ced4a3SJed Brown static PetscErrorCode PetscBuildTwoSidedGetType(BuildTwoSidedType *twosided) 25*f6ced4a3SJed Brown { 26*f6ced4a3SJed Brown PetscErrorCode ierr; 27*f6ced4a3SJed Brown 28*f6ced4a3SJed Brown PetscFunctionBegin; 29*f6ced4a3SJed Brown if (_twosided_type == BUILDTWOSIDED_NOTSET) { 30*f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 31*f6ced4a3SJed Brown _twosided_type = BUILDTWOSIDED_IBARRIER; 32*f6ced4a3SJed Brown #else 33*f6ced4a3SJed Brown _twosided_type = BUILDTWOSIDED_ALLREDUCE; 34*f6ced4a3SJed Brown #endif 35*f6ced4a3SJed Brown ierr = PetscOptionsGetEnum(PETSC_NULL,"-build_twosided",BuildTwoSidedTypes,(PetscEnum*)&_twosided_type,PETSC_NULL);CHKERRQ(ierr); 36*f6ced4a3SJed Brown } 37*f6ced4a3SJed Brown *twosided = _twosided_type; 38*f6ced4a3SJed Brown PetscFunctionReturn(0); 39*f6ced4a3SJed Brown } 40*f6ced4a3SJed Brown 41*f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 42*f6ced4a3SJed Brown /* Segmented (extendable) array implementation */ 43*f6ced4a3SJed Brown typedef struct _SegArray *SegArray; 44*f6ced4a3SJed Brown struct _SegArray { 45*f6ced4a3SJed Brown PetscInt unitbytes; 46*f6ced4a3SJed Brown PetscInt alloc; 47*f6ced4a3SJed Brown PetscInt used; 48*f6ced4a3SJed Brown PetscInt tailused; 49*f6ced4a3SJed Brown SegArray tail; 50*f6ced4a3SJed Brown union { /* Dummy types to ensure alignment */ 51*f6ced4a3SJed Brown PetscReal dummy_real; 52*f6ced4a3SJed Brown PetscInt dummy_int; 53*f6ced4a3SJed Brown char array[1]; 54*f6ced4a3SJed Brown } u; 55*f6ced4a3SJed Brown }; 56*f6ced4a3SJed Brown 57*f6ced4a3SJed Brown #undef __FUNCT__ 58*f6ced4a3SJed Brown #define __FUNCT__ "SegArrayCreate" 59*f6ced4a3SJed Brown static PetscErrorCode SegArrayCreate(PetscInt unitbytes,PetscInt expected,SegArray *seg) 60*f6ced4a3SJed Brown { 61*f6ced4a3SJed Brown PetscErrorCode ierr; 62*f6ced4a3SJed Brown 63*f6ced4a3SJed Brown PetscFunctionBegin; 64*f6ced4a3SJed Brown ierr = PetscMalloc(offsetof(struct _SegArray,u)+expected*unitbytes,seg);CHKERRQ(ierr); 65*f6ced4a3SJed Brown ierr = PetscMemzero(*seg,offsetof(struct _SegArray,u));CHKERRQ(ierr); 66*f6ced4a3SJed Brown (*seg)->unitbytes = unitbytes; 67*f6ced4a3SJed Brown (*seg)->alloc = expected; 68*f6ced4a3SJed Brown PetscFunctionReturn(0); 69*f6ced4a3SJed Brown } 70*f6ced4a3SJed Brown 71*f6ced4a3SJed Brown #undef __FUNCT__ 72*f6ced4a3SJed Brown #define __FUNCT__ "SegArrayAlloc_Private" 73*f6ced4a3SJed Brown static PetscErrorCode SegArrayAlloc_Private(SegArray *seg,PetscInt count) 74*f6ced4a3SJed Brown { 75*f6ced4a3SJed Brown PetscErrorCode ierr; 76*f6ced4a3SJed Brown SegArray newseg,s; 77*f6ced4a3SJed Brown PetscInt alloc; 78*f6ced4a3SJed Brown 79*f6ced4a3SJed Brown PetscFunctionBegin; 80*f6ced4a3SJed Brown s = *seg; 81*f6ced4a3SJed Brown /* Grow at least fast enough to hold next item, like Fibonacci otherwise (up to 1MB chunks) */ 82*f6ced4a3SJed Brown alloc = PetscMax(s->used+count,PetscMin(1000000/s->unitbytes+1,s->alloc+s->tailused)); 83*f6ced4a3SJed Brown ierr = PetscMalloc(offsetof(struct _SegArray,u)+alloc*s->unitbytes,&newseg);CHKERRQ(ierr); 84*f6ced4a3SJed Brown ierr = PetscMemzero(newseg,offsetof(struct _SegArray,u));CHKERRQ(ierr); 85*f6ced4a3SJed Brown newseg->unitbytes = s->unitbytes; 86*f6ced4a3SJed Brown newseg->tailused = s->used + s->tailused; 87*f6ced4a3SJed Brown newseg->tail = s; 88*f6ced4a3SJed Brown newseg->alloc = alloc; 89*f6ced4a3SJed Brown *seg = newseg; 90*f6ced4a3SJed Brown PetscFunctionReturn(0); 91*f6ced4a3SJed Brown } 92*f6ced4a3SJed Brown 93*f6ced4a3SJed Brown #undef __FUNCT__ 94*f6ced4a3SJed Brown #define __FUNCT__ "SegArrayGet" 95*f6ced4a3SJed Brown static PetscErrorCode SegArrayGet(SegArray *seg,PetscInt count,void *array) 96*f6ced4a3SJed Brown { 97*f6ced4a3SJed Brown PetscErrorCode ierr; 98*f6ced4a3SJed Brown SegArray s; 99*f6ced4a3SJed Brown 100*f6ced4a3SJed Brown PetscFunctionBegin; 101*f6ced4a3SJed Brown s = *seg; 102*f6ced4a3SJed Brown if (PetscUnlikely(s->used + count > s->alloc)) {ierr = SegArrayAlloc_Private(seg,count);CHKERRQ(ierr);} 103*f6ced4a3SJed Brown s = *seg; 104*f6ced4a3SJed Brown *(char**)array = &s->u.array[s->used*s->unitbytes]; 105*f6ced4a3SJed Brown s->used += count; 106*f6ced4a3SJed Brown PetscFunctionReturn(0); 107*f6ced4a3SJed Brown } 108*f6ced4a3SJed Brown 109*f6ced4a3SJed Brown #undef __FUNCT__ 110*f6ced4a3SJed Brown #define __FUNCT__ "SegArrayDestroy" 111*f6ced4a3SJed Brown static PetscErrorCode SegArrayDestroy(SegArray *seg) 112*f6ced4a3SJed Brown { 113*f6ced4a3SJed Brown PetscErrorCode ierr; 114*f6ced4a3SJed Brown SegArray s; 115*f6ced4a3SJed Brown 116*f6ced4a3SJed Brown PetscFunctionBegin; 117*f6ced4a3SJed Brown for (s=*seg; s; ) { 118*f6ced4a3SJed Brown SegArray tail = s->tail; 119*f6ced4a3SJed Brown ierr = PetscFree(s);CHKERRQ(ierr); 120*f6ced4a3SJed Brown s = tail; 121*f6ced4a3SJed Brown } 122*f6ced4a3SJed Brown *seg = PETSC_NULL; 123*f6ced4a3SJed Brown PetscFunctionReturn(0); 124*f6ced4a3SJed Brown } 125*f6ced4a3SJed Brown 126*f6ced4a3SJed Brown #undef __FUNCT__ 127*f6ced4a3SJed Brown #define __FUNCT__ "SegArrayExtract" 128*f6ced4a3SJed Brown /* Extracts contiguous data and resets segarray */ 129*f6ced4a3SJed Brown static PetscErrorCode SegArrayExtract(SegArray *seg,void *contiguous) 130*f6ced4a3SJed Brown { 131*f6ced4a3SJed Brown PetscErrorCode ierr; 132*f6ced4a3SJed Brown PetscInt unitbytes; 133*f6ced4a3SJed Brown SegArray s,t; 134*f6ced4a3SJed Brown char *contig,*ptr; 135*f6ced4a3SJed Brown 136*f6ced4a3SJed Brown PetscFunctionBegin; 137*f6ced4a3SJed Brown s = *seg; 138*f6ced4a3SJed Brown unitbytes = s->unitbytes; 139*f6ced4a3SJed Brown ierr = PetscMalloc((s->used+s->tailused)*unitbytes,&contig);CHKERRQ(ierr); 140*f6ced4a3SJed Brown ptr = contig + s->tailused*unitbytes; 141*f6ced4a3SJed Brown ierr = PetscMemcpy(ptr,s->u.array,s->used*unitbytes);CHKERRQ(ierr); 142*f6ced4a3SJed Brown for (t=s->tail; t; ) { 143*f6ced4a3SJed Brown SegArray tail = t->tail; 144*f6ced4a3SJed Brown ptr -= t->used*unitbytes; 145*f6ced4a3SJed Brown ierr = PetscMemcpy(ptr,t->u.array,t->used*unitbytes);CHKERRQ(ierr); 146*f6ced4a3SJed Brown ierr = PetscFree(t);CHKERRQ(ierr); 147*f6ced4a3SJed Brown t = tail; 148*f6ced4a3SJed Brown } 149*f6ced4a3SJed Brown if (ptr != contig) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Tail count does not match"); 150*f6ced4a3SJed Brown s->tailused = 0; 151*f6ced4a3SJed Brown s->tail = PETSC_NULL; 152*f6ced4a3SJed Brown *(char**)contiguous = contig; 153*f6ced4a3SJed Brown PetscFunctionReturn(0); 154*f6ced4a3SJed Brown } 155*f6ced4a3SJed Brown 156*f6ced4a3SJed Brown #undef __FUNCT__ 157*f6ced4a3SJed Brown #define __FUNCT__ "BuildTwoSided_Ibarrier" 158*f6ced4a3SJed Brown static PetscErrorCode BuildTwoSided_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata) 159*f6ced4a3SJed Brown { 160*f6ced4a3SJed Brown PetscErrorCode ierr; 161*f6ced4a3SJed Brown PetscMPIInt nrecvs,tag,unitbytes,done; 162*f6ced4a3SJed Brown PetscInt i; 163*f6ced4a3SJed Brown char *tdata; 164*f6ced4a3SJed Brown MPI_Request *sendreqs,barrier; 165*f6ced4a3SJed Brown SegArray segrank,segdata; 166*f6ced4a3SJed Brown 167*f6ced4a3SJed Brown PetscFunctionBegin; 168*f6ced4a3SJed Brown ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 169*f6ced4a3SJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 170*f6ced4a3SJed Brown tdata = (char*)todata; 171*f6ced4a3SJed Brown ierr = PetscMalloc(nto*sizeof(MPI_Request),&sendreqs);CHKERRQ(ierr); 172*f6ced4a3SJed Brown for (i=0; i<nto; i++) { 173*f6ced4a3SJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 174*f6ced4a3SJed Brown } 175*f6ced4a3SJed Brown ierr = SegArrayCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 176*f6ced4a3SJed Brown ierr = SegArrayCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 177*f6ced4a3SJed Brown 178*f6ced4a3SJed Brown nrecvs = 0; 179*f6ced4a3SJed Brown barrier = MPI_REQUEST_NULL; 180*f6ced4a3SJed Brown for (done=0; !done; ) { 181*f6ced4a3SJed Brown PetscMPIInt flag; 182*f6ced4a3SJed Brown MPI_Status status; 183*f6ced4a3SJed Brown ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr); 184*f6ced4a3SJed Brown if (flag) { /* incoming message */ 185*f6ced4a3SJed Brown PetscMPIInt *recvrank; 186*f6ced4a3SJed Brown void *buf; 187*f6ced4a3SJed Brown ierr = SegArrayGet(&segrank,1,&recvrank);CHKERRQ(ierr); 188*f6ced4a3SJed Brown ierr = SegArrayGet(&segdata,count,&buf);CHKERRQ(ierr); 189*f6ced4a3SJed Brown *recvrank = status.MPI_SOURCE; 190*f6ced4a3SJed Brown ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 191*f6ced4a3SJed Brown nrecvs++; 192*f6ced4a3SJed Brown } 193*f6ced4a3SJed Brown if (barrier == MPI_REQUEST_NULL) { 194*f6ced4a3SJed Brown PetscMPIInt sent,nsends = PetscMPIIntCast(nto); 195*f6ced4a3SJed Brown ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 196*f6ced4a3SJed Brown if (sent) { 197*f6ced4a3SJed Brown ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr); 198*f6ced4a3SJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 199*f6ced4a3SJed Brown } 200*f6ced4a3SJed Brown } else { 201*f6ced4a3SJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 202*f6ced4a3SJed Brown } 203*f6ced4a3SJed Brown } 204*f6ced4a3SJed Brown *nfrom = nrecvs; 205*f6ced4a3SJed Brown ierr = SegArrayExtract(&segrank,fromranks);CHKERRQ(ierr); 206*f6ced4a3SJed Brown ierr = SegArrayDestroy(&segrank);CHKERRQ(ierr); 207*f6ced4a3SJed Brown ierr = SegArrayExtract(&segdata,fromdata);CHKERRQ(ierr); 208*f6ced4a3SJed Brown ierr = SegArrayDestroy(&segdata);CHKERRQ(ierr); 209*f6ced4a3SJed Brown PetscFunctionReturn(0); 210*f6ced4a3SJed Brown } 211*f6ced4a3SJed Brown #endif 212*f6ced4a3SJed Brown 213*f6ced4a3SJed Brown #undef __FUNCT__ 214*f6ced4a3SJed Brown #define __FUNCT__ "BuildTwoSided_Allreduce" 215*f6ced4a3SJed Brown static PetscErrorCode BuildTwoSided_Allreduce(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscInt nto,const PetscMPIInt *toranks,const void *todata,PetscInt *nfrom,PetscMPIInt **fromranks,void *fromdata) 216*f6ced4a3SJed Brown { 217*f6ced4a3SJed Brown PetscErrorCode ierr; 218*f6ced4a3SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,unitbytes,*franks; 219*f6ced4a3SJed Brown PetscInt i; 220*f6ced4a3SJed Brown char *tdata,*fdata; 221*f6ced4a3SJed Brown MPI_Request *reqs,*sendreqs; 222*f6ced4a3SJed Brown MPI_Status *statuses; 223*f6ced4a3SJed Brown 224*f6ced4a3SJed Brown PetscFunctionBegin; 225*f6ced4a3SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 226*f6ced4a3SJed Brown ierr = PetscMalloc(size*sizeof(*iflags),&iflags);CHKERRQ(ierr); 227*f6ced4a3SJed Brown ierr = PetscMemzero(iflags,size*sizeof(*iflags));CHKERRQ(ierr); 228*f6ced4a3SJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 229*f6ced4a3SJed Brown ierr = PetscGatherNumberOfMessages(comm,iflags,PETSC_NULL,&nrecvs);CHKERRQ(ierr); 230*f6ced4a3SJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 231*f6ced4a3SJed Brown 232*f6ced4a3SJed Brown ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 233*f6ced4a3SJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 234*f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 235*f6ced4a3SJed Brown tdata = (char*)todata; 236*f6ced4a3SJed Brown ierr = PetscMalloc2(nto+nrecvs,MPI_Request,&reqs,nto+nrecvs,MPI_Status,&statuses);CHKERRQ(ierr); 237*f6ced4a3SJed Brown sendreqs = reqs + nrecvs; 238*f6ced4a3SJed Brown for (i=0; i<nrecvs; i++) { 239*f6ced4a3SJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 240*f6ced4a3SJed Brown } 241*f6ced4a3SJed Brown for (i=0; i<nto; i++) { 242*f6ced4a3SJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 243*f6ced4a3SJed Brown } 244*f6ced4a3SJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 245*f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*sizeof(PetscMPIInt),&franks);CHKERRQ(ierr); 246*f6ced4a3SJed Brown for (i=0; i<nrecvs; i++) { 247*f6ced4a3SJed Brown franks[i] = statuses[i].MPI_SOURCE; 248*f6ced4a3SJed Brown } 249*f6ced4a3SJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 250*f6ced4a3SJed Brown 251*f6ced4a3SJed Brown *nfrom = nrecvs; 252*f6ced4a3SJed Brown *fromranks = franks; 253*f6ced4a3SJed Brown *(void**)fromdata = fdata; 254*f6ced4a3SJed Brown PetscFunctionReturn(0); 255*f6ced4a3SJed Brown } 256*f6ced4a3SJed Brown 257*f6ced4a3SJed Brown #undef __FUNCT__ 258*f6ced4a3SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided" 259*f6ced4a3SJed Brown /*@C 260*f6ced4a3SJed Brown PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths) 261*f6ced4a3SJed Brown 262*f6ced4a3SJed Brown Collective on MPI_Comm 263*f6ced4a3SJed Brown 264*f6ced4a3SJed Brown Input Arguments: 265*f6ced4a3SJed Brown + comm - communicator 266*f6ced4a3SJed Brown . count - number of entries to send/receive (must match on all ranks) 267*f6ced4a3SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 268*f6ced4a3SJed Brown . nto - number of ranks to send data to 269*f6ced4a3SJed Brown . toranks - ranks to send to (array of length nto) 270*f6ced4a3SJed Brown - todata - data to send to each rank (packed) 271*f6ced4a3SJed Brown 272*f6ced4a3SJed Brown Output Arguments: 273*f6ced4a3SJed Brown + nfrom - number of ranks receiving messages from 274*f6ced4a3SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 275*f6ced4a3SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 276*f6ced4a3SJed Brown 277*f6ced4a3SJed Brown Level: developer 278*f6ced4a3SJed Brown 279*f6ced4a3SJed Brown Notes: 280*f6ced4a3SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 281*f6ced4a3SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data. 282*f6ced4a3SJed Brown 283*f6ced4a3SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 284*f6ced4a3SJed Brown 285*f6ced4a3SJed Brown References: 286*f6ced4a3SJed Brown The MPI_Ibarrier implementation uses the algorithm in 287*f6ced4a3SJed Brown Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010. 288*f6ced4a3SJed Brown 289*f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 290*f6ced4a3SJed Brown @*/ 291*f6ced4a3SJed 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) 292*f6ced4a3SJed Brown { 293*f6ced4a3SJed Brown PetscErrorCode ierr; 294*f6ced4a3SJed Brown BuildTwoSidedType buildtype; 295*f6ced4a3SJed Brown 296*f6ced4a3SJed Brown PetscFunctionBegin; 297*f6ced4a3SJed Brown ierr = PetscBuildTwoSidedGetType(&buildtype);CHKERRQ(ierr); 298*f6ced4a3SJed Brown switch (buildtype) { 299*f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 300*f6ced4a3SJed Brown case BUILDTWOSIDED_IBARRIER: 301*f6ced4a3SJed Brown ierr = BuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 302*f6ced4a3SJed Brown break; 303*f6ced4a3SJed Brown #endif 304*f6ced4a3SJed Brown case BUILDTWOSIDED_ALLREDUCE: 305*f6ced4a3SJed Brown ierr = BuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 306*f6ced4a3SJed Brown break; 307*f6ced4a3SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 308*f6ced4a3SJed Brown } 309*f6ced4a3SJed Brown PetscFunctionReturn(0); 310*f6ced4a3SJed Brown } 311