1f6ced4a3SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 2f6ced4a3SJed Brown #include <stddef.h> 3f6ced4a3SJed Brown 46145cd65SJed Brown const char *const PetscBuildTwoSidedTypes[] = { 5f6ced4a3SJed Brown "ALLREDUCE", 66145cd65SJed Brown "IBARRIER", 76145cd65SJed Brown "PetscBuildTwoSidedType", 86145cd65SJed Brown "PETSC_BUILDTWOSIDED_", 9f6ced4a3SJed Brown 0 10f6ced4a3SJed Brown }; 11f6ced4a3SJed Brown 126145cd65SJed Brown static PetscBuildTwoSidedType _twosided_type = PETSC_BUILDTWOSIDED_NOTSET; 13f6ced4a3SJed Brown 14f6ced4a3SJed Brown #undef __FUNCT__ 156145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedSetType" 166145cd65SJed Brown /*@ 176145cd65SJed Brown PetscCommBuildTwoSidedSetType - set algorithm to use when building two-sided communication 186145cd65SJed Brown 196145cd65SJed Brown Logically Collective 206145cd65SJed Brown 216145cd65SJed Brown Input Arguments: 226145cd65SJed Brown + comm - PETSC_COMM_WORLD 236145cd65SJed Brown - twosided - algorithm to use in subsequent calls to PetscCommBuildTwoSided() 246145cd65SJed Brown 256145cd65SJed Brown Level: developer 266145cd65SJed Brown 276145cd65SJed Brown Note: 286145cd65SJed Brown This option is currently global, but could be made per-communicator. 296145cd65SJed Brown 306145cd65SJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedGetType() 316145cd65SJed Brown @*/ 326145cd65SJed Brown PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm comm,PetscBuildTwoSidedType twosided) 336145cd65SJed Brown { 346145cd65SJed Brown PetscFunctionBegin; 356145cd65SJed Brown #if defined(PETSC_USE_DEBUG) 366145cd65SJed Brown { /* We don't have a PetscObject so can't use PetscValidLogicalCollectiveEnum */ 376145cd65SJed Brown PetscMPIInt ierr; 386145cd65SJed Brown PetscMPIInt b1[2],b2[2]; 396145cd65SJed Brown b1[0] = -(PetscMPIInt)twosided; 406145cd65SJed Brown b1[1] = (PetscMPIInt)twosided; 416145cd65SJed Brown ierr = MPI_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); 426145cd65SJed Brown if (-b2[0] != b2[1]) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes"); 436145cd65SJed Brown } 446145cd65SJed Brown #endif 456145cd65SJed Brown _twosided_type = twosided; 466145cd65SJed Brown PetscFunctionReturn(0); 476145cd65SJed Brown } 486145cd65SJed Brown 496145cd65SJed Brown #undef __FUNCT__ 506145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSidedGetType" 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; 67f6ced4a3SJed Brown 68f6ced4a3SJed Brown PetscFunctionBegin; 696145cd65SJed Brown *twosided = PETSC_BUILDTWOSIDED_NOTSET; 706145cd65SJed Brown if (_twosided_type == PETSC_BUILDTWOSIDED_NOTSET) { 71f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 726145cd65SJed Brown # if defined(PETSC_HAVE_MPICH_CH3_SOCK) && !defined(PETSC_HAVE_MPICH_CH3_SOCK_FIXED_NBC_PROGRESS) 736145cd65SJed Brown /* Deadlock in Ibarrier: http://trac.mpich.org/projects/mpich/ticket/1785 */ 746145cd65SJed Brown _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; 75f6ced4a3SJed Brown # else 766145cd65SJed Brown _twosided_type = PETSC_BUILDTWOSIDED_IBARRIER; 77f6ced4a3SJed Brown # endif 786145cd65SJed Brown #else 796145cd65SJed Brown _twosided_type = PETSC_BUILDTWOSIDED_ALLREDUCE; 806145cd65SJed Brown #endif 81*0298fd71SBarry Smith ierr = PetscOptionsGetEnum(NULL,"-build_twosided",PetscBuildTwoSidedTypes,(PetscEnum*)&_twosided_type,NULL);CHKERRQ(ierr); 82f6ced4a3SJed Brown } 83f6ced4a3SJed Brown *twosided = _twosided_type; 84f6ced4a3SJed Brown PetscFunctionReturn(0); 85f6ced4a3SJed Brown } 86f6ced4a3SJed Brown 87f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 88f6ced4a3SJed Brown /* Segmented (extendable) array implementation */ 89f6ced4a3SJed Brown typedef struct _SegArray *SegArray; 90f6ced4a3SJed Brown struct _SegArray { 91f6ced4a3SJed Brown PetscInt unitbytes; 92f6ced4a3SJed Brown PetscInt alloc; 93f6ced4a3SJed Brown PetscInt used; 94f6ced4a3SJed Brown PetscInt tailused; 95f6ced4a3SJed Brown SegArray tail; 96f6ced4a3SJed Brown union { /* Dummy types to ensure alignment */ 97f6ced4a3SJed Brown PetscReal dummy_real; 98f6ced4a3SJed Brown PetscInt dummy_int; 99f6ced4a3SJed Brown char array[1]; 100f6ced4a3SJed Brown } u; 101f6ced4a3SJed Brown }; 102f6ced4a3SJed Brown 103f6ced4a3SJed Brown #undef __FUNCT__ 104f6ced4a3SJed Brown #define __FUNCT__ "SegArrayCreate" 105f6ced4a3SJed Brown static PetscErrorCode SegArrayCreate(PetscInt unitbytes,PetscInt expected,SegArray *seg) 106f6ced4a3SJed Brown { 107f6ced4a3SJed Brown PetscErrorCode ierr; 108f6ced4a3SJed Brown 109f6ced4a3SJed Brown PetscFunctionBegin; 110f6ced4a3SJed Brown ierr = PetscMalloc(offsetof(struct _SegArray,u)+expected*unitbytes,seg);CHKERRQ(ierr); 111f6ced4a3SJed Brown ierr = PetscMemzero(*seg,offsetof(struct _SegArray,u));CHKERRQ(ierr); 112a297a907SKarl Rupp 113f6ced4a3SJed Brown (*seg)->unitbytes = unitbytes; 114f6ced4a3SJed Brown (*seg)->alloc = expected; 115f6ced4a3SJed Brown PetscFunctionReturn(0); 116f6ced4a3SJed Brown } 117f6ced4a3SJed Brown 118f6ced4a3SJed Brown #undef __FUNCT__ 119f6ced4a3SJed Brown #define __FUNCT__ "SegArrayAlloc_Private" 120f6ced4a3SJed Brown static PetscErrorCode SegArrayAlloc_Private(SegArray *seg,PetscInt count) 121f6ced4a3SJed Brown { 122f6ced4a3SJed Brown PetscErrorCode ierr; 123f6ced4a3SJed Brown SegArray newseg,s; 124f6ced4a3SJed Brown PetscInt alloc; 125f6ced4a3SJed Brown 126f6ced4a3SJed Brown PetscFunctionBegin; 127f6ced4a3SJed Brown s = *seg; 128f6ced4a3SJed Brown /* Grow at least fast enough to hold next item, like Fibonacci otherwise (up to 1MB chunks) */ 129f6ced4a3SJed Brown alloc = PetscMax(s->used+count,PetscMin(1000000/s->unitbytes+1,s->alloc+s->tailused)); 130f6ced4a3SJed Brown ierr = PetscMalloc(offsetof(struct _SegArray,u)+alloc*s->unitbytes,&newseg);CHKERRQ(ierr); 131f6ced4a3SJed Brown ierr = PetscMemzero(newseg,offsetof(struct _SegArray,u));CHKERRQ(ierr); 132a297a907SKarl Rupp 133f6ced4a3SJed Brown newseg->unitbytes = s->unitbytes; 134f6ced4a3SJed Brown newseg->tailused = s->used + s->tailused; 135f6ced4a3SJed Brown newseg->tail = s; 136f6ced4a3SJed Brown newseg->alloc = alloc; 137f6ced4a3SJed Brown *seg = newseg; 138f6ced4a3SJed Brown PetscFunctionReturn(0); 139f6ced4a3SJed Brown } 140f6ced4a3SJed Brown 141f6ced4a3SJed Brown #undef __FUNCT__ 142f6ced4a3SJed Brown #define __FUNCT__ "SegArrayGet" 143f6ced4a3SJed Brown static PetscErrorCode SegArrayGet(SegArray *seg,PetscInt count,void *array) 144f6ced4a3SJed Brown { 145f6ced4a3SJed Brown PetscErrorCode ierr; 146f6ced4a3SJed Brown SegArray s; 147f6ced4a3SJed Brown 148f6ced4a3SJed Brown PetscFunctionBegin; 149f6ced4a3SJed Brown s = *seg; 150f6ced4a3SJed Brown if (PetscUnlikely(s->used + count > s->alloc)) {ierr = SegArrayAlloc_Private(seg,count);CHKERRQ(ierr);} 151f6ced4a3SJed Brown s = *seg; 152f6ced4a3SJed Brown *(char**)array = &s->u.array[s->used*s->unitbytes]; 153f6ced4a3SJed Brown s->used += count; 154f6ced4a3SJed Brown PetscFunctionReturn(0); 155f6ced4a3SJed Brown } 156f6ced4a3SJed Brown 157f6ced4a3SJed Brown #undef __FUNCT__ 158f6ced4a3SJed Brown #define __FUNCT__ "SegArrayDestroy" 159f6ced4a3SJed Brown static PetscErrorCode SegArrayDestroy(SegArray *seg) 160f6ced4a3SJed Brown { 161f6ced4a3SJed Brown PetscErrorCode ierr; 162f6ced4a3SJed Brown SegArray s; 163f6ced4a3SJed Brown 164f6ced4a3SJed Brown PetscFunctionBegin; 165f6ced4a3SJed Brown for (s=*seg; s;) { 166f6ced4a3SJed Brown SegArray tail = s->tail; 167f6ced4a3SJed Brown ierr = PetscFree(s);CHKERRQ(ierr); 168f6ced4a3SJed Brown s = tail; 169f6ced4a3SJed Brown } 170*0298fd71SBarry Smith *seg = NULL; 171f6ced4a3SJed Brown PetscFunctionReturn(0); 172f6ced4a3SJed Brown } 173f6ced4a3SJed Brown 174f6ced4a3SJed Brown #undef __FUNCT__ 175f6ced4a3SJed Brown #define __FUNCT__ "SegArrayExtract" 176f6ced4a3SJed Brown /* Extracts contiguous data and resets segarray */ 177f6ced4a3SJed Brown static PetscErrorCode SegArrayExtract(SegArray *seg,void *contiguous) 178f6ced4a3SJed Brown { 179f6ced4a3SJed Brown PetscErrorCode ierr; 180f6ced4a3SJed Brown PetscInt unitbytes; 181f6ced4a3SJed Brown SegArray s,t; 182f6ced4a3SJed Brown char *contig,*ptr; 183f6ced4a3SJed Brown 184f6ced4a3SJed Brown PetscFunctionBegin; 185f6ced4a3SJed Brown s = *seg; 186a297a907SKarl Rupp 187f6ced4a3SJed Brown unitbytes = s->unitbytes; 188a297a907SKarl Rupp 189f6ced4a3SJed Brown ierr = PetscMalloc((s->used+s->tailused)*unitbytes,&contig);CHKERRQ(ierr); 190f6ced4a3SJed Brown ptr = contig + s->tailused*unitbytes; 191f6ced4a3SJed Brown ierr = PetscMemcpy(ptr,s->u.array,s->used*unitbytes);CHKERRQ(ierr); 192f6ced4a3SJed Brown for (t=s->tail; t;) { 193f6ced4a3SJed Brown SegArray tail = t->tail; 194f6ced4a3SJed Brown ptr -= t->used*unitbytes; 195f6ced4a3SJed Brown ierr = PetscMemcpy(ptr,t->u.array,t->used*unitbytes);CHKERRQ(ierr); 196f6ced4a3SJed Brown ierr = PetscFree(t);CHKERRQ(ierr); 197f6ced4a3SJed Brown t = tail; 198f6ced4a3SJed Brown } 199f6ced4a3SJed Brown if (ptr != contig) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Tail count does not match"); 200f6ced4a3SJed Brown s->tailused = 0; 201*0298fd71SBarry Smith s->tail = NULL; 202f6ced4a3SJed Brown *(char**)contiguous = contig; 203f6ced4a3SJed Brown PetscFunctionReturn(0); 204f6ced4a3SJed Brown } 205f6ced4a3SJed Brown 206f6ced4a3SJed Brown #undef __FUNCT__ 2076145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Ibarrier" 2086145cd65SJed 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) 209f6ced4a3SJed Brown { 210f6ced4a3SJed Brown PetscErrorCode ierr; 211f6ced4a3SJed Brown PetscMPIInt nrecvs,tag,unitbytes,done; 212f6ced4a3SJed Brown PetscInt i; 213f6ced4a3SJed Brown char *tdata; 214f6ced4a3SJed Brown MPI_Request *sendreqs,barrier; 215f6ced4a3SJed Brown SegArray segrank,segdata; 216f6ced4a3SJed Brown 217f6ced4a3SJed Brown PetscFunctionBegin; 218f6ced4a3SJed Brown ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 219f6ced4a3SJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 220f6ced4a3SJed Brown tdata = (char*)todata; 221f6ced4a3SJed Brown ierr = PetscMalloc(nto*sizeof(MPI_Request),&sendreqs);CHKERRQ(ierr); 222f6ced4a3SJed Brown for (i=0; i<nto; i++) { 223f6ced4a3SJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 224f6ced4a3SJed Brown } 225f6ced4a3SJed Brown ierr = SegArrayCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 226f6ced4a3SJed Brown ierr = SegArrayCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 227f6ced4a3SJed Brown 228f6ced4a3SJed Brown nrecvs = 0; 229f6ced4a3SJed Brown barrier = MPI_REQUEST_NULL; 230f6ced4a3SJed Brown for (done=0; !done; ) { 231f6ced4a3SJed Brown PetscMPIInt flag; 232f6ced4a3SJed Brown MPI_Status status; 233f6ced4a3SJed Brown ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr); 234f6ced4a3SJed Brown if (flag) { /* incoming message */ 235f6ced4a3SJed Brown PetscMPIInt *recvrank; 236f6ced4a3SJed Brown void *buf; 237f6ced4a3SJed Brown ierr = SegArrayGet(&segrank,1,&recvrank);CHKERRQ(ierr); 238f6ced4a3SJed Brown ierr = SegArrayGet(&segdata,count,&buf);CHKERRQ(ierr); 239f6ced4a3SJed Brown *recvrank = status.MPI_SOURCE; 240f6ced4a3SJed Brown ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 241f6ced4a3SJed Brown nrecvs++; 242f6ced4a3SJed Brown } 243f6ced4a3SJed Brown if (barrier == MPI_REQUEST_NULL) { 2444dc2109aSBarry Smith PetscMPIInt sent,nsends; 2454dc2109aSBarry Smith ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr); 246f6ced4a3SJed Brown ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 247f6ced4a3SJed Brown if (sent) { 248f6ced4a3SJed Brown ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr); 249f6ced4a3SJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 250f6ced4a3SJed Brown } 251f6ced4a3SJed Brown } else { 252f6ced4a3SJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 253f6ced4a3SJed Brown } 254f6ced4a3SJed Brown } 255f6ced4a3SJed Brown *nfrom = nrecvs; 256f6ced4a3SJed Brown ierr = SegArrayExtract(&segrank,fromranks);CHKERRQ(ierr); 257f6ced4a3SJed Brown ierr = SegArrayDestroy(&segrank);CHKERRQ(ierr); 258f6ced4a3SJed Brown ierr = SegArrayExtract(&segdata,fromdata);CHKERRQ(ierr); 259f6ced4a3SJed Brown ierr = SegArrayDestroy(&segdata);CHKERRQ(ierr); 260f6ced4a3SJed Brown PetscFunctionReturn(0); 261f6ced4a3SJed Brown } 262f6ced4a3SJed Brown #endif 263f6ced4a3SJed Brown 264f6ced4a3SJed Brown #undef __FUNCT__ 2656145cd65SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided_Allreduce" 2666145cd65SJed 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) 267f6ced4a3SJed Brown { 268f6ced4a3SJed Brown PetscErrorCode ierr; 269f6ced4a3SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,unitbytes,*franks; 270f6ced4a3SJed Brown PetscInt i; 271f6ced4a3SJed Brown char *tdata,*fdata; 272f6ced4a3SJed Brown MPI_Request *reqs,*sendreqs; 273f6ced4a3SJed Brown MPI_Status *statuses; 274f6ced4a3SJed Brown 275f6ced4a3SJed Brown PetscFunctionBegin; 276f6ced4a3SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 277f6ced4a3SJed Brown ierr = PetscMalloc(size*sizeof(*iflags),&iflags);CHKERRQ(ierr); 278f6ced4a3SJed Brown ierr = PetscMemzero(iflags,size*sizeof(*iflags));CHKERRQ(ierr); 279f6ced4a3SJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 280*0298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,iflags,NULL,&nrecvs);CHKERRQ(ierr); 281f6ced4a3SJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 282f6ced4a3SJed Brown 283f6ced4a3SJed Brown ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 284f6ced4a3SJed Brown ierr = MPI_Type_size(dtype,&unitbytes);CHKERRQ(ierr); 285f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 286f6ced4a3SJed Brown tdata = (char*)todata; 287f6ced4a3SJed Brown ierr = PetscMalloc2(nto+nrecvs,MPI_Request,&reqs,nto+nrecvs,MPI_Status,&statuses);CHKERRQ(ierr); 288f6ced4a3SJed Brown sendreqs = reqs + nrecvs; 289f6ced4a3SJed Brown for (i=0; i<nrecvs; i++) { 290f6ced4a3SJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 291f6ced4a3SJed Brown } 292f6ced4a3SJed Brown for (i=0; i<nto; i++) { 293f6ced4a3SJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 294f6ced4a3SJed Brown } 295f6ced4a3SJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 296f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*sizeof(PetscMPIInt),&franks);CHKERRQ(ierr); 297a297a907SKarl Rupp for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 298f6ced4a3SJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 299f6ced4a3SJed Brown 300f6ced4a3SJed Brown *nfrom = nrecvs; 301f6ced4a3SJed Brown *fromranks = franks; 302f6ced4a3SJed Brown *(void**)fromdata = fdata; 303f6ced4a3SJed Brown PetscFunctionReturn(0); 304f6ced4a3SJed Brown } 305f6ced4a3SJed Brown 306f6ced4a3SJed Brown #undef __FUNCT__ 307f6ced4a3SJed Brown #define __FUNCT__ "PetscCommBuildTwoSided" 308f6ced4a3SJed Brown /*@C 309f6ced4a3SJed Brown PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths) 310f6ced4a3SJed Brown 311f6ced4a3SJed Brown Collective on MPI_Comm 312f6ced4a3SJed Brown 313f6ced4a3SJed Brown Input Arguments: 314f6ced4a3SJed Brown + comm - communicator 315f6ced4a3SJed Brown . count - number of entries to send/receive (must match on all ranks) 316f6ced4a3SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 317f6ced4a3SJed Brown . nto - number of ranks to send data to 318f6ced4a3SJed Brown . toranks - ranks to send to (array of length nto) 319f6ced4a3SJed Brown - todata - data to send to each rank (packed) 320f6ced4a3SJed Brown 321f6ced4a3SJed Brown Output Arguments: 322f6ced4a3SJed Brown + nfrom - number of ranks receiving messages from 323f6ced4a3SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 324f6ced4a3SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 325f6ced4a3SJed Brown 326f6ced4a3SJed Brown Level: developer 327f6ced4a3SJed Brown 328f6ced4a3SJed Brown Notes: 329f6ced4a3SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 330f6ced4a3SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data. 331f6ced4a3SJed Brown 332f6ced4a3SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 333f6ced4a3SJed Brown 334f6ced4a3SJed Brown References: 335f6ced4a3SJed Brown The MPI_Ibarrier implementation uses the algorithm in 336f6ced4a3SJed Brown Hoefler, Siebert and Lumsdaine, Scalable communication protocols for dynamic sparse data exchange, 2010. 337f6ced4a3SJed Brown 338f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 339f6ced4a3SJed Brown @*/ 340f6ced4a3SJed 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) 341f6ced4a3SJed Brown { 342f6ced4a3SJed Brown PetscErrorCode ierr; 3436145cd65SJed Brown PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET; 344f6ced4a3SJed Brown 345f6ced4a3SJed Brown PetscFunctionBegin; 3466145cd65SJed Brown ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr); 347f6ced4a3SJed Brown switch (buildtype) { 3486145cd65SJed Brown case PETSC_BUILDTWOSIDED_IBARRIER: 349f6ced4a3SJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 3506145cd65SJed Brown ierr = PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 3516145cd65SJed Brown #else 3526145cd65SJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)"); 353f6ced4a3SJed Brown #endif 3546145cd65SJed Brown break; 3556145cd65SJed Brown case PETSC_BUILDTWOSIDED_ALLREDUCE: 3566145cd65SJed Brown ierr = PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 357f6ced4a3SJed Brown break; 358f6ced4a3SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 359f6ced4a3SJed Brown } 360f6ced4a3SJed Brown PetscFunctionReturn(0); 361f6ced4a3SJed Brown } 362