1f6ced4a3SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 2*95c0884eSLisandro Dalcin #include <petsc/private/petscimpl.h> 3f6ced4a3SJed Brown 4*95c0884eSLisandro Dalcin PetscLogEvent PETSC_BuildTwoSided; 5*95c0884eSLisandro Dalcin PetscLogEvent PETSC_BuildTwoSidedF; 63b3561c8SJed Brown 76145cd65SJed Brown const char *const PetscBuildTwoSidedTypes[] = { 8f6ced4a3SJed Brown "ALLREDUCE", 96145cd65SJed Brown "IBARRIER", 101654bf6bSJed Brown "REDSCATTER", 116145cd65SJed Brown "PetscBuildTwoSidedType", 126145cd65SJed Brown "PETSC_BUILDTWOSIDED_", 13f6ced4a3SJed Brown 0 14f6ced4a3SJed Brown }; 15f6ced4a3SJed Brown 166145cd65SJed Brown static PetscBuildTwoSidedType _twosided_type = PETSC_BUILDTWOSIDED_NOTSET; 17f6ced4a3SJed Brown 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; 43b2566f29SBarry Smith ierr = MPIU_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 /*@ 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 81c5929fdfSBarry Smith ierr = PetscOptionsGetEnum(NULL,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 87b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 88f6ced4a3SJed Brown 89cf4b5b4fSJed 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) 90f6ced4a3SJed Brown { 91f6ced4a3SJed Brown PetscErrorCode ierr; 92c3a59b84SJed Brown PetscMPIInt nrecvs,tag,done,i; 93c3a59b84SJed Brown MPI_Aint lb,unitbytes; 94f6ced4a3SJed Brown char *tdata; 95f6ced4a3SJed Brown MPI_Request *sendreqs,barrier; 960f453b92SJed Brown PetscSegBuffer segrank,segdata; 9746e50bb6SJed Brown PetscBool barrier_started; 98f6ced4a3SJed Brown 99f6ced4a3SJed Brown PetscFunctionBegin; 100a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 101c3a59b84SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 102c3a59b84SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 103f6ced4a3SJed Brown tdata = (char*)todata; 104785e854fSJed Brown ierr = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr); 105f6ced4a3SJed Brown for (i=0; i<nto; i++) { 106f6ced4a3SJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 107f6ced4a3SJed Brown } 1080f453b92SJed Brown ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 1090f453b92SJed Brown ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 110f6ced4a3SJed Brown 111f6ced4a3SJed Brown nrecvs = 0; 112f6ced4a3SJed Brown barrier = MPI_REQUEST_NULL; 11346e50bb6SJed Brown /* MPICH-3.2 sometimes does not create a request in some "optimized" cases. This is arguably a standard violation, 11446e50bb6SJed Brown * but we need to work around it. */ 11546e50bb6SJed Brown barrier_started = PETSC_FALSE; 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 } 12946e50bb6SJed Brown if (!barrier_started) { 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 13946e50bb6SJed Brown barrier_started = PETSC_TRUE; 140f6ced4a3SJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 141f6ced4a3SJed Brown } 142f6ced4a3SJed Brown } else { 143f6ced4a3SJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 144f6ced4a3SJed Brown } 145f6ced4a3SJed Brown } 146f6ced4a3SJed Brown *nfrom = nrecvs; 147137cf7b6SJed Brown ierr = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr); 1480f453b92SJed Brown ierr = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr); 149137cf7b6SJed Brown ierr = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr); 1500f453b92SJed Brown ierr = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr); 151a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 152f6ced4a3SJed Brown PetscFunctionReturn(0); 153f6ced4a3SJed Brown } 154f6ced4a3SJed Brown #endif 155f6ced4a3SJed Brown 156cf4b5b4fSJed 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) 157f6ced4a3SJed Brown { 158f6ced4a3SJed Brown PetscErrorCode ierr; 159c3a59b84SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,*franks,i; 160c3a59b84SJed Brown MPI_Aint lb,unitbytes; 161f6ced4a3SJed Brown char *tdata,*fdata; 162f6ced4a3SJed Brown MPI_Request *reqs,*sendreqs; 163f6ced4a3SJed Brown MPI_Status *statuses; 164f6ced4a3SJed Brown 165f6ced4a3SJed Brown PetscFunctionBegin; 166f6ced4a3SJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1671795a4d1SJed Brown ierr = PetscCalloc1(size,&iflags);CHKERRQ(ierr); 168f6ced4a3SJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 1690298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,iflags,NULL,&nrecvs);CHKERRQ(ierr); 170f6ced4a3SJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 171f6ced4a3SJed Brown 172a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 173c3a59b84SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 174c3a59b84SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 175f6ced4a3SJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 176f6ced4a3SJed Brown tdata = (char*)todata; 177dcca6d9dSJed Brown ierr = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr); 17806926cafSJed Brown sendreqs = reqs + nrecvs; 179f6ced4a3SJed Brown for (i=0; i<nrecvs; i++) { 180f6ced4a3SJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 181f6ced4a3SJed Brown } 182f6ced4a3SJed Brown for (i=0; i<nto; i++) { 183f6ced4a3SJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 184f6ced4a3SJed Brown } 185f6ced4a3SJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 186785e854fSJed Brown ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr); 187a297a907SKarl Rupp for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 18806926cafSJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 189a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 190f6ced4a3SJed Brown 191f6ced4a3SJed Brown *nfrom = nrecvs; 192f6ced4a3SJed Brown *fromranks = franks; 193f6ced4a3SJed Brown *(void**)fromdata = fdata; 194f6ced4a3SJed Brown PetscFunctionReturn(0); 195f6ced4a3SJed Brown } 196f6ced4a3SJed Brown 1971654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 198cf4b5b4fSJed 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) 1991654bf6bSJed Brown { 2001654bf6bSJed Brown PetscErrorCode ierr; 201c3a59b84SJed Brown PetscMPIInt size,*iflags,nrecvs,tag,*franks,i; 202c3a59b84SJed Brown MPI_Aint lb,unitbytes; 2031654bf6bSJed Brown char *tdata,*fdata; 2041654bf6bSJed Brown MPI_Request *reqs,*sendreqs; 2051654bf6bSJed Brown MPI_Status *statuses; 2061654bf6bSJed Brown 2071654bf6bSJed Brown PetscFunctionBegin; 2081654bf6bSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 2091654bf6bSJed Brown ierr = PetscMalloc1(size,&iflags);CHKERRQ(ierr); 2101654bf6bSJed Brown ierr = PetscMemzero(iflags,size*sizeof(*iflags));CHKERRQ(ierr); 2111654bf6bSJed Brown for (i=0; i<nto; i++) iflags[toranks[i]] = 1; 2121654bf6bSJed Brown ierr = MPI_Reduce_scatter_block(iflags,&nrecvs,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 2131654bf6bSJed Brown ierr = PetscFree(iflags);CHKERRQ(ierr); 2141654bf6bSJed Brown 215a502f807SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 216c3a59b84SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 217c3a59b84SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 2181654bf6bSJed Brown ierr = PetscMalloc(nrecvs*count*unitbytes,&fdata);CHKERRQ(ierr); 2191654bf6bSJed Brown tdata = (char*)todata; 2201654bf6bSJed Brown ierr = PetscMalloc2(nto+nrecvs,&reqs,nto+nrecvs,&statuses);CHKERRQ(ierr); 2211654bf6bSJed Brown sendreqs = reqs + nrecvs; 2221654bf6bSJed Brown for (i=0; i<nrecvs; i++) { 2231654bf6bSJed Brown ierr = MPI_Irecv((void*)(fdata+count*unitbytes*i),count,dtype,MPI_ANY_SOURCE,tag,comm,reqs+i);CHKERRQ(ierr); 2241654bf6bSJed Brown } 2251654bf6bSJed Brown for (i=0; i<nto; i++) { 2261654bf6bSJed Brown ierr = MPI_Isend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 2271654bf6bSJed Brown } 2281654bf6bSJed Brown ierr = MPI_Waitall(nto+nrecvs,reqs,statuses);CHKERRQ(ierr); 2291654bf6bSJed Brown ierr = PetscMalloc1(nrecvs,&franks);CHKERRQ(ierr); 2301654bf6bSJed Brown for (i=0; i<nrecvs; i++) franks[i] = statuses[i].MPI_SOURCE; 2311654bf6bSJed Brown ierr = PetscFree2(reqs,statuses);CHKERRQ(ierr); 232a502f807SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 2331654bf6bSJed Brown 2341654bf6bSJed Brown *nfrom = nrecvs; 2351654bf6bSJed Brown *fromranks = franks; 2361654bf6bSJed Brown *(void**)fromdata = fdata; 2371654bf6bSJed Brown PetscFunctionReturn(0); 2381654bf6bSJed Brown } 2391654bf6bSJed Brown #endif 2401654bf6bSJed Brown 241f6ced4a3SJed Brown /*@C 242f6ced4a3SJed Brown PetscCommBuildTwoSided - discovers communicating ranks given one-sided information, moving constant-sized data in the process (often message lengths) 243f6ced4a3SJed Brown 244f6ced4a3SJed Brown Collective on MPI_Comm 245f6ced4a3SJed Brown 246f6ced4a3SJed Brown Input Arguments: 247f6ced4a3SJed Brown + comm - communicator 248f6ced4a3SJed Brown . count - number of entries to send/receive (must match on all ranks) 249f6ced4a3SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 250f6ced4a3SJed Brown . nto - number of ranks to send data to 251f6ced4a3SJed Brown . toranks - ranks to send to (array of length nto) 252f6ced4a3SJed Brown - todata - data to send to each rank (packed) 253f6ced4a3SJed Brown 254f6ced4a3SJed Brown Output Arguments: 255f6ced4a3SJed Brown + nfrom - number of ranks receiving messages from 256f6ced4a3SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 257f6ced4a3SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 258f6ced4a3SJed Brown 259f6ced4a3SJed Brown Level: developer 260f6ced4a3SJed Brown 2611654bf6bSJed Brown Options Database Keys: 2621654bf6bSJed Brown . -build_twosided <allreduce|ibarrier|redscatter> - algorithm to set up two-sided communication 2631654bf6bSJed Brown 264f6ced4a3SJed Brown Notes: 265f6ced4a3SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 266f6ced4a3SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other constant-size data. 267f6ced4a3SJed Brown 268f6ced4a3SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 269f6ced4a3SJed Brown 270f6ced4a3SJed Brown References: 27196a0c994SBarry Smith . 1. - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in 2724f02bc6aSBarry Smith Scalable communication protocols for dynamic sparse data exchange, 2010. 273f6ced4a3SJed Brown 274f6ced4a3SJed Brown .seealso: PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 275f6ced4a3SJed Brown @*/ 276cf4b5b4fSJed 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) 277f6ced4a3SJed Brown { 278f6ced4a3SJed Brown PetscErrorCode ierr; 2796145cd65SJed Brown PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET; 280f6ced4a3SJed Brown 281f6ced4a3SJed Brown PetscFunctionBegin; 2823b3561c8SJed Brown ierr = PetscSysInitializePackage();CHKERRQ(ierr); 2833b3561c8SJed Brown ierr = PetscLogEventBegin(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr); 2846145cd65SJed Brown ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr); 285f6ced4a3SJed Brown switch (buildtype) { 2866145cd65SJed Brown case PETSC_BUILDTWOSIDED_IBARRIER: 287b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 2886145cd65SJed Brown ierr = PetscCommBuildTwoSided_Ibarrier(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 2896145cd65SJed Brown #else 2906145cd65SJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)"); 291f6ced4a3SJed Brown #endif 2926145cd65SJed Brown break; 2936145cd65SJed Brown case PETSC_BUILDTWOSIDED_ALLREDUCE: 2946145cd65SJed Brown ierr = PetscCommBuildTwoSided_Allreduce(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 295f6ced4a3SJed Brown break; 2961654bf6bSJed Brown case PETSC_BUILDTWOSIDED_REDSCATTER: 2971654bf6bSJed Brown #if defined(PETSC_HAVE_MPI_REDUCE_SCATTER_BLOCK) 2981654bf6bSJed Brown ierr = PetscCommBuildTwoSided_RedScatter(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 2991654bf6bSJed Brown #else 3001654bf6bSJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Reduce_scatter_block (part of MPI-2.2)"); 3011654bf6bSJed Brown #endif 3021654bf6bSJed Brown break; 303f6ced4a3SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 304f6ced4a3SJed Brown } 3053b3561c8SJed Brown ierr = PetscLogEventEnd(PETSC_BuildTwoSided,0,0,0,0);CHKERRQ(ierr); 306f6ced4a3SJed Brown PetscFunctionReturn(0); 307f6ced4a3SJed Brown } 308d815da10SJed Brown 3096afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Reference(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata, 3106afbe7ddSJed Brown PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs, 311d815da10SJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 312d815da10SJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 313d815da10SJed Brown { 314d815da10SJed Brown PetscErrorCode ierr; 315d815da10SJed Brown PetscMPIInt i,*tag; 316d815da10SJed Brown MPI_Aint lb,unitbytes; 317d815da10SJed Brown MPI_Request *sendreq,*recvreq; 318d815da10SJed Brown 319d815da10SJed Brown PetscFunctionBegin; 320d815da10SJed Brown ierr = PetscMalloc1(ntags,&tag);CHKERRQ(ierr); 321d815da10SJed Brown if (ntags > 0) { 322d815da10SJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag[0]);CHKERRQ(ierr); 323d815da10SJed Brown } 324d815da10SJed Brown for (i=1; i<ntags; i++) { 325d815da10SJed Brown ierr = PetscCommGetNewTag(comm,&tag[i]);CHKERRQ(ierr); 326d815da10SJed Brown } 327d815da10SJed Brown 328d815da10SJed Brown /* Perform complete initial rendezvous */ 329d815da10SJed Brown ierr = PetscCommBuildTwoSided(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata);CHKERRQ(ierr); 330d815da10SJed Brown 3316afbe7ddSJed Brown ierr = PetscMalloc1(nto*ntags,&sendreq);CHKERRQ(ierr); 3326afbe7ddSJed Brown ierr = PetscMalloc1(*nfrom*ntags,&recvreq);CHKERRQ(ierr); 333d815da10SJed Brown 334d815da10SJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 335d815da10SJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 336d815da10SJed Brown for (i=0; i<nto; i++) { 33789157a57SJed Brown PetscMPIInt k; 33889157a57SJed Brown for (k=0; k<ntags; k++) sendreq[i*ntags+k] = MPI_REQUEST_NULL; 339d815da10SJed Brown ierr = (*send)(comm,tag,i,toranks[i],((char*)todata)+count*unitbytes*i,sendreq+i*ntags,ctx);CHKERRQ(ierr); 340d815da10SJed Brown } 341d815da10SJed Brown for (i=0; i<*nfrom; i++) { 342d815da10SJed Brown void *header = (*(char**)fromdata) + count*unitbytes*i; 34389157a57SJed Brown PetscMPIInt k; 34489157a57SJed Brown for (k=0; k<ntags; k++) recvreq[i*ntags+k] = MPI_REQUEST_NULL; 345d815da10SJed Brown ierr = (*recv)(comm,tag,(*fromranks)[i],header,recvreq+i*ntags,ctx);CHKERRQ(ierr); 346d815da10SJed Brown } 347d815da10SJed Brown ierr = PetscFree(tag);CHKERRQ(ierr); 348d815da10SJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 3496afbe7ddSJed Brown *toreqs = sendreq; 3506afbe7ddSJed Brown *fromreqs = recvreq; 351d815da10SJed Brown PetscFunctionReturn(0); 352d815da10SJed Brown } 353d815da10SJed Brown 354b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 355afb254cdSJed Brown 3566afbe7ddSJed Brown static PetscErrorCode PetscCommBuildTwoSidedFReq_Ibarrier(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata, 3576afbe7ddSJed Brown PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs, 358afb254cdSJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 359afb254cdSJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 360afb254cdSJed Brown { 361afb254cdSJed Brown PetscErrorCode ierr; 362afb254cdSJed Brown PetscMPIInt nrecvs,tag,*tags,done,i; 363afb254cdSJed Brown MPI_Aint lb,unitbytes; 364afb254cdSJed Brown char *tdata; 365afb254cdSJed Brown MPI_Request *sendreqs,*usendreqs,*req,barrier; 366afb254cdSJed Brown PetscSegBuffer segrank,segdata,segreq; 36746e50bb6SJed Brown PetscBool barrier_started; 368afb254cdSJed Brown 369afb254cdSJed Brown PetscFunctionBegin; 370afb254cdSJed Brown ierr = PetscCommDuplicate(comm,&comm,&tag);CHKERRQ(ierr); 371afb254cdSJed Brown ierr = PetscMalloc1(ntags,&tags);CHKERRQ(ierr); 372afb254cdSJed Brown for (i=0; i<ntags; i++) { 373afb254cdSJed Brown ierr = PetscCommGetNewTag(comm,&tags[i]);CHKERRQ(ierr); 374afb254cdSJed Brown } 375afb254cdSJed Brown ierr = MPI_Type_get_extent(dtype,&lb,&unitbytes);CHKERRQ(ierr); 376afb254cdSJed Brown if (lb != 0) SETERRQ1(comm,PETSC_ERR_SUP,"Datatype with nonzero lower bound %ld\n",(long)lb); 377afb254cdSJed Brown tdata = (char*)todata; 3786afbe7ddSJed Brown ierr = PetscMalloc1(nto,&sendreqs);CHKERRQ(ierr); 3796afbe7ddSJed Brown ierr = PetscMalloc1(nto*ntags,&usendreqs);CHKERRQ(ierr); 380afb254cdSJed Brown /* Post synchronous sends */ 381afb254cdSJed Brown for (i=0; i<nto; i++) { 382afb254cdSJed Brown ierr = MPI_Issend((void*)(tdata+count*unitbytes*i),count,dtype,toranks[i],tag,comm,sendreqs+i);CHKERRQ(ierr); 383afb254cdSJed Brown } 384afb254cdSJed Brown /* Post actual payloads. These are typically larger messages. Hopefully sending these later does not slow down the 385afb254cdSJed Brown * synchronous messages above. */ 386afb254cdSJed Brown for (i=0; i<nto; i++) { 38789157a57SJed Brown PetscMPIInt k; 38889157a57SJed Brown for (k=0; k<ntags; k++) usendreqs[i*ntags+k] = MPI_REQUEST_NULL; 389afb254cdSJed Brown ierr = (*send)(comm,tags,i,toranks[i],tdata+count*unitbytes*i,usendreqs+i*ntags,ctx);CHKERRQ(ierr); 390afb254cdSJed Brown } 391afb254cdSJed Brown 392afb254cdSJed Brown ierr = PetscSegBufferCreate(sizeof(PetscMPIInt),4,&segrank);CHKERRQ(ierr); 393afb254cdSJed Brown ierr = PetscSegBufferCreate(unitbytes,4*count,&segdata);CHKERRQ(ierr); 394afb254cdSJed Brown ierr = PetscSegBufferCreate(sizeof(MPI_Request),4,&segreq);CHKERRQ(ierr); 395afb254cdSJed Brown 396afb254cdSJed Brown nrecvs = 0; 397afb254cdSJed Brown barrier = MPI_REQUEST_NULL; 39846e50bb6SJed Brown /* MPICH-3.2 sometimes does not create a request in some "optimized" cases. This is arguably a standard violation, 39946e50bb6SJed Brown * but we need to work around it. */ 40046e50bb6SJed Brown barrier_started = PETSC_FALSE; 401afb254cdSJed Brown for (done=0; !done; ) { 402afb254cdSJed Brown PetscMPIInt flag; 403afb254cdSJed Brown MPI_Status status; 404afb254cdSJed Brown ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag,comm,&flag,&status);CHKERRQ(ierr); 405afb254cdSJed Brown if (flag) { /* incoming message */ 40689157a57SJed Brown PetscMPIInt *recvrank,k; 407afb254cdSJed Brown void *buf; 408afb254cdSJed Brown ierr = PetscSegBufferGet(segrank,1,&recvrank);CHKERRQ(ierr); 409afb254cdSJed Brown ierr = PetscSegBufferGet(segdata,count,&buf);CHKERRQ(ierr); 410afb254cdSJed Brown *recvrank = status.MPI_SOURCE; 411afb254cdSJed Brown ierr = MPI_Recv(buf,count,dtype,status.MPI_SOURCE,tag,comm,MPI_STATUS_IGNORE);CHKERRQ(ierr); 412afb254cdSJed Brown ierr = PetscSegBufferGet(segreq,ntags,&req);CHKERRQ(ierr); 41389157a57SJed Brown for (k=0; k<ntags; k++) req[k] = MPI_REQUEST_NULL; 414afb254cdSJed Brown ierr = (*recv)(comm,tags,status.MPI_SOURCE,buf,req,ctx);CHKERRQ(ierr); 415afb254cdSJed Brown nrecvs++; 416afb254cdSJed Brown } 41746e50bb6SJed Brown if (!barrier_started) { 418afb254cdSJed Brown PetscMPIInt sent,nsends; 419afb254cdSJed Brown ierr = PetscMPIIntCast(nto,&nsends);CHKERRQ(ierr); 420afb254cdSJed Brown ierr = MPI_Testall(nsends,sendreqs,&sent,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 421afb254cdSJed Brown if (sent) { 422b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) 423afb254cdSJed Brown ierr = MPI_Ibarrier(comm,&barrier);CHKERRQ(ierr); 424b5a01e9fSJed Brown #elif defined(PETSC_HAVE_MPIX_IBARRIER) 425b5a01e9fSJed Brown ierr = MPIX_Ibarrier(comm,&barrier);CHKERRQ(ierr); 426b5a01e9fSJed Brown #endif 42746e50bb6SJed Brown barrier_started = PETSC_TRUE; 428afb254cdSJed Brown } 429afb254cdSJed Brown } else { 430afb254cdSJed Brown ierr = MPI_Test(&barrier,&done,MPI_STATUS_IGNORE);CHKERRQ(ierr); 431afb254cdSJed Brown } 432afb254cdSJed Brown } 433afb254cdSJed Brown *nfrom = nrecvs; 434afb254cdSJed Brown ierr = PetscSegBufferExtractAlloc(segrank,fromranks);CHKERRQ(ierr); 435afb254cdSJed Brown ierr = PetscSegBufferDestroy(&segrank);CHKERRQ(ierr); 436afb254cdSJed Brown ierr = PetscSegBufferExtractAlloc(segdata,fromdata);CHKERRQ(ierr); 437afb254cdSJed Brown ierr = PetscSegBufferDestroy(&segdata);CHKERRQ(ierr); 4386afbe7ddSJed Brown *toreqs = usendreqs; 4396afbe7ddSJed Brown ierr = PetscSegBufferExtractAlloc(segreq,fromreqs);CHKERRQ(ierr); 440afb254cdSJed Brown ierr = PetscSegBufferDestroy(&segreq);CHKERRQ(ierr); 4416afbe7ddSJed Brown ierr = PetscFree(sendreqs);CHKERRQ(ierr); 442afb254cdSJed Brown ierr = PetscFree(tags);CHKERRQ(ierr); 443afb254cdSJed Brown ierr = PetscCommDestroy(&comm);CHKERRQ(ierr); 444afb254cdSJed Brown PetscFunctionReturn(0); 445afb254cdSJed Brown } 446afb254cdSJed Brown #endif 447afb254cdSJed Brown 448d815da10SJed Brown /*@C 449d815da10SJed Brown PetscCommBuildTwoSidedF - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous 450d815da10SJed Brown 451d815da10SJed Brown Collective on MPI_Comm 452d815da10SJed Brown 453d815da10SJed Brown Input Arguments: 454d815da10SJed Brown + comm - communicator 455d815da10SJed Brown . count - number of entries to send/receive in initial rendezvous (must match on all ranks) 456d815da10SJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 457d815da10SJed Brown . nto - number of ranks to send data to 458d815da10SJed Brown . toranks - ranks to send to (array of length nto) 459d815da10SJed Brown . todata - data to send to each rank (packed) 460d815da10SJed Brown . ntags - number of tags needed by send/recv callbacks 461d815da10SJed Brown . send - callback invoked on sending process when ready to send primary payload 462d815da10SJed Brown . recv - callback invoked on receiving process after delivery of rendezvous message 463d815da10SJed Brown - ctx - context for callbacks 464d815da10SJed Brown 465d815da10SJed Brown Output Arguments: 466d815da10SJed Brown + nfrom - number of ranks receiving messages from 467d815da10SJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 468d815da10SJed Brown - fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 469d815da10SJed Brown 470d815da10SJed Brown Level: developer 471d815da10SJed Brown 472d815da10SJed Brown Notes: 473d815da10SJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 474d815da10SJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data. 475d815da10SJed Brown 476d815da10SJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 477d815da10SJed Brown 478d815da10SJed Brown References: 47996a0c994SBarry Smith . 1. - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in 4804f02bc6aSBarry Smith Scalable communication protocols for dynamic sparse data exchange, 2010. 481d815da10SJed Brown 4826afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedFReq(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 483d815da10SJed Brown @*/ 484d815da10SJed 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, 485d815da10SJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 486d815da10SJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 487d815da10SJed Brown { 4886afbe7ddSJed Brown PetscErrorCode ierr; 4896afbe7ddSJed Brown MPI_Request *toreqs,*fromreqs; 4906afbe7ddSJed Brown 4916afbe7ddSJed Brown PetscFunctionBegin; 4926afbe7ddSJed Brown ierr = PetscCommBuildTwoSidedFReq(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,&toreqs,&fromreqs,send,recv,ctx);CHKERRQ(ierr); 4936afbe7ddSJed Brown ierr = MPI_Waitall(nto*ntags,toreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 4946afbe7ddSJed Brown ierr = MPI_Waitall(*nfrom*ntags,fromreqs,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 4956afbe7ddSJed Brown ierr = PetscFree(toreqs);CHKERRQ(ierr); 4966afbe7ddSJed Brown ierr = PetscFree(fromreqs);CHKERRQ(ierr); 4976afbe7ddSJed Brown PetscFunctionReturn(0); 4986afbe7ddSJed Brown } 4996afbe7ddSJed Brown 5006afbe7ddSJed Brown /*@C 5016afbe7ddSJed Brown PetscCommBuildTwoSidedFReq - discovers communicating ranks given one-sided information, calling user-defined functions during rendezvous, returns requests 5026afbe7ddSJed Brown 5036afbe7ddSJed Brown Collective on MPI_Comm 5046afbe7ddSJed Brown 5056afbe7ddSJed Brown Input Arguments: 5066afbe7ddSJed Brown + comm - communicator 5076afbe7ddSJed Brown . count - number of entries to send/receive in initial rendezvous (must match on all ranks) 5086afbe7ddSJed Brown . dtype - datatype to send/receive from each rank (must match on all ranks) 5096afbe7ddSJed Brown . nto - number of ranks to send data to 5106afbe7ddSJed Brown . toranks - ranks to send to (array of length nto) 5116afbe7ddSJed Brown . todata - data to send to each rank (packed) 5126afbe7ddSJed Brown . ntags - number of tags needed by send/recv callbacks 5136afbe7ddSJed Brown . send - callback invoked on sending process when ready to send primary payload 5146afbe7ddSJed Brown . recv - callback invoked on receiving process after delivery of rendezvous message 5156afbe7ddSJed Brown - ctx - context for callbacks 5166afbe7ddSJed Brown 5176afbe7ddSJed Brown Output Arguments: 5186afbe7ddSJed Brown + nfrom - number of ranks receiving messages from 5196afbe7ddSJed Brown . fromranks - ranks receiving messages from (length nfrom; caller should PetscFree()) 5206afbe7ddSJed Brown . fromdata - packed data from each rank, each with count entries of type dtype (length nfrom, caller responsible for PetscFree()) 5216afbe7ddSJed Brown . toreqs - array of nto*ntags sender requests (caller must wait on these, then PetscFree()) 5226afbe7ddSJed Brown - fromreqs - array of nfrom*ntags receiver requests (caller must wait on these, then PetscFree()) 5236afbe7ddSJed Brown 5246afbe7ddSJed Brown Level: developer 5256afbe7ddSJed Brown 5266afbe7ddSJed Brown Notes: 5276afbe7ddSJed Brown This memory-scalable interface is an alternative to calling PetscGatherNumberOfMessages() and 5286afbe7ddSJed Brown PetscGatherMessageLengths(), possibly with a subsequent round of communication to send other data. 5296afbe7ddSJed Brown 5306afbe7ddSJed Brown Basic data types as well as contiguous types are supported, but non-contiguous (e.g., strided) types are not. 5316afbe7ddSJed Brown 5326afbe7ddSJed Brown References: 53396a0c994SBarry Smith . 1. - Hoefler, Siebert and Lumsdaine, The MPI_Ibarrier implementation uses the algorithm in 5344f02bc6aSBarry Smith Scalable communication protocols for dynamic sparse data exchange, 2010. 5356afbe7ddSJed Brown 5366afbe7ddSJed Brown .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedF(), PetscGatherNumberOfMessages(), PetscGatherMessageLengths() 5376afbe7ddSJed Brown @*/ 5386afbe7ddSJed Brown PetscErrorCode PetscCommBuildTwoSidedFReq(MPI_Comm comm,PetscMPIInt count,MPI_Datatype dtype,PetscMPIInt nto,const PetscMPIInt *toranks,const void *todata, 5396afbe7ddSJed Brown PetscMPIInt *nfrom,PetscMPIInt **fromranks,void *fromdata,PetscMPIInt ntags,MPI_Request **toreqs,MPI_Request **fromreqs, 5406afbe7ddSJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 5416afbe7ddSJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx) 5426afbe7ddSJed Brown { 5436afbe7ddSJed Brown PetscErrorCode ierr,(*f)(MPI_Comm,PetscMPIInt,MPI_Datatype,PetscMPIInt,const PetscMPIInt[],const void*, 5446afbe7ddSJed Brown PetscMPIInt*,PetscMPIInt**,void*,PetscMPIInt,MPI_Request**,MPI_Request**, 545d815da10SJed Brown PetscErrorCode (*send)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,PetscMPIInt,void*,MPI_Request[],void*), 546d815da10SJed Brown PetscErrorCode (*recv)(MPI_Comm,const PetscMPIInt[],PetscMPIInt,void*,MPI_Request[],void*),void *ctx); 547d815da10SJed Brown PetscBuildTwoSidedType buildtype = PETSC_BUILDTWOSIDED_NOTSET; 5483461502dSJed Brown PetscMPIInt i,size; 549d815da10SJed Brown 550d815da10SJed Brown PetscFunctionBegin; 551d815da10SJed Brown ierr = PetscSysInitializePackage();CHKERRQ(ierr); 5523461502dSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 5533461502dSJed Brown for (i=0; i<nto; i++) { 5543461502dSJed Brown if (toranks[i] < 0 || size <= toranks[i]) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"toranks[%d] %d not in comm size %d",i,toranks[i],size); 5553461502dSJed Brown } 55627a32ea5SJed Brown ierr = PetscLogEventBegin(PETSC_BuildTwoSidedF,0,0,0,0);CHKERRQ(ierr); 557d815da10SJed Brown ierr = PetscCommBuildTwoSidedGetType(comm,&buildtype);CHKERRQ(ierr); 558d815da10SJed Brown switch (buildtype) { 559d815da10SJed Brown case PETSC_BUILDTWOSIDED_IBARRIER: 560b5a01e9fSJed Brown #if defined(PETSC_HAVE_MPI_IBARRIER) || defined(PETSC_HAVE_MPIX_IBARRIER) 5616afbe7ddSJed Brown f = PetscCommBuildTwoSidedFReq_Ibarrier; 562afb254cdSJed Brown #else 563afb254cdSJed Brown SETERRQ(comm,PETSC_ERR_PLIB,"MPI implementation does not provide MPI_Ibarrier (part of MPI-3)"); 564afb254cdSJed Brown #endif 565afb254cdSJed Brown break; 566d815da10SJed Brown case PETSC_BUILDTWOSIDED_ALLREDUCE: 567d815da10SJed Brown case PETSC_BUILDTWOSIDED_REDSCATTER: 5686afbe7ddSJed Brown f = PetscCommBuildTwoSidedFReq_Reference; 569d815da10SJed Brown break; 570d815da10SJed Brown default: SETERRQ(comm,PETSC_ERR_PLIB,"Unknown method for building two-sided communication"); 571d815da10SJed Brown } 5726afbe7ddSJed Brown ierr = (*f)(comm,count,dtype,nto,toranks,todata,nfrom,fromranks,fromdata,ntags,toreqs,fromreqs,send,recv,ctx);CHKERRQ(ierr); 57327a32ea5SJed Brown ierr = PetscLogEventEnd(PETSC_BuildTwoSidedF,0,0,0,0);CHKERRQ(ierr); 574d815da10SJed Brown PetscFunctionReturn(0); 575d815da10SJed Brown } 576