1*dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h> 2*dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/allgather/sfallgather.h> 3*dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/gatherv/sfgatherv.h> 4*dd5b3ca6SJunchao Zhang 5*dd5b3ca6SJunchao Zhang typedef PetscSFPack_Allgatherv PetscSFPack_Alltoall; 6*dd5b3ca6SJunchao Zhang #define PetscSFPackGet_Alltoall PetscSFPackGet_Allgatherv 7*dd5b3ca6SJunchao Zhang 8*dd5b3ca6SJunchao Zhang /* Reuse the type. The difference is some fields (i.e., displs, recvcounts) are not used, which is not a big deal */ 9*dd5b3ca6SJunchao Zhang typedef PetscSF_Allgatherv PetscSF_Alltoall; 10*dd5b3ca6SJunchao Zhang 11*dd5b3ca6SJunchao Zhang /*===================================================================================*/ 12*dd5b3ca6SJunchao Zhang /* Implementations of SF public APIs */ 13*dd5b3ca6SJunchao Zhang /*===================================================================================*/ 14*dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFGetGraph_Alltoall(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 15*dd5b3ca6SJunchao Zhang { 16*dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 17*dd5b3ca6SJunchao Zhang PetscInt i; 18*dd5b3ca6SJunchao Zhang 19*dd5b3ca6SJunchao Zhang PetscFunctionBegin; 20*dd5b3ca6SJunchao Zhang if (nroots) *nroots = sf->nroots; 21*dd5b3ca6SJunchao Zhang if (nleaves) *nleaves = sf->nleaves; 22*dd5b3ca6SJunchao Zhang if (ilocal) *ilocal = NULL; /* Contiguous local indices */ 23*dd5b3ca6SJunchao Zhang if (iremote) { 24*dd5b3ca6SJunchao Zhang if (!sf->remote) { 25*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(sf->nleaves,&sf->remote);CHKERRQ(ierr); 26*dd5b3ca6SJunchao Zhang sf->remote_alloc = sf->remote; 27*dd5b3ca6SJunchao Zhang for (i=0; i<sf->nleaves; i++) { 28*dd5b3ca6SJunchao Zhang sf->remote[i].rank = i; 29*dd5b3ca6SJunchao Zhang sf->remote[i].index = i; 30*dd5b3ca6SJunchao Zhang } 31*dd5b3ca6SJunchao Zhang } 32*dd5b3ca6SJunchao Zhang *iremote = sf->remote; 33*dd5b3ca6SJunchao Zhang } 34*dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 35*dd5b3ca6SJunchao Zhang } 36*dd5b3ca6SJunchao Zhang 37*dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFBcastAndOpBegin_Alltoall(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op) 38*dd5b3ca6SJunchao Zhang { 39*dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 40*dd5b3ca6SJunchao Zhang PetscSFPack_Alltoall link; 41*dd5b3ca6SJunchao Zhang MPI_Comm comm; 42*dd5b3ca6SJunchao Zhang void *recvbuf; 43*dd5b3ca6SJunchao Zhang 44*dd5b3ca6SJunchao Zhang PetscFunctionBegin; 45*dd5b3ca6SJunchao Zhang ierr = PetscSFPackGet_Alltoall(sf,unit,rootdata,&link);CHKERRQ(ierr); 46*dd5b3ca6SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 47*dd5b3ca6SJunchao Zhang 48*dd5b3ca6SJunchao Zhang if (op != MPIU_REPLACE) { 49*dd5b3ca6SJunchao Zhang if (!link->leaf) {ierr = PetscMalloc(sf->nleaves*link->unitbytes,&link->leaf);CHKERRQ(ierr);} 50*dd5b3ca6SJunchao Zhang recvbuf = link->leaf; 51*dd5b3ca6SJunchao Zhang } else { 52*dd5b3ca6SJunchao Zhang recvbuf = leafdata; 53*dd5b3ca6SJunchao Zhang } 54*dd5b3ca6SJunchao Zhang ierr = MPIU_Ialltoall(rootdata,1,unit,recvbuf,1,unit,comm,&link->request);CHKERRQ(ierr); 55*dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 56*dd5b3ca6SJunchao Zhang } 57*dd5b3ca6SJunchao Zhang 58*dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFReduceBegin_Alltoall(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op) 59*dd5b3ca6SJunchao Zhang { 60*dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 61*dd5b3ca6SJunchao Zhang PetscSFPack_Alltoall link; 62*dd5b3ca6SJunchao Zhang MPI_Comm comm; 63*dd5b3ca6SJunchao Zhang void *recvbuf; 64*dd5b3ca6SJunchao Zhang 65*dd5b3ca6SJunchao Zhang PetscFunctionBegin; 66*dd5b3ca6SJunchao Zhang ierr = PetscSFPackGet_Alltoall(sf,unit,leafdata,&link);CHKERRQ(ierr); 67*dd5b3ca6SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 68*dd5b3ca6SJunchao Zhang 69*dd5b3ca6SJunchao Zhang if (op != MPIU_REPLACE) { 70*dd5b3ca6SJunchao Zhang if (!link->root) {ierr = PetscMalloc(sf->nroots*link->unitbytes,&link->root);CHKERRQ(ierr);} 71*dd5b3ca6SJunchao Zhang recvbuf = link->root; 72*dd5b3ca6SJunchao Zhang } else { 73*dd5b3ca6SJunchao Zhang recvbuf = rootdata; 74*dd5b3ca6SJunchao Zhang } 75*dd5b3ca6SJunchao Zhang ierr = MPIU_Ialltoall(leafdata,1,unit,recvbuf,1,unit,comm,&link->request);CHKERRQ(ierr); 76*dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 77*dd5b3ca6SJunchao Zhang } 78*dd5b3ca6SJunchao Zhang 79*dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFCreateLocalSF_Alltoall(PetscSF sf,PetscSF *out) 80*dd5b3ca6SJunchao Zhang { 81*dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 82*dd5b3ca6SJunchao Zhang PetscInt nroots=1,nleaves=1,*ilocal; 83*dd5b3ca6SJunchao Zhang PetscSFNode *iremote = NULL; 84*dd5b3ca6SJunchao Zhang PetscSF lsf; 85*dd5b3ca6SJunchao Zhang PetscMPIInt rank; 86*dd5b3ca6SJunchao Zhang 87*dd5b3ca6SJunchao Zhang PetscFunctionBegin; 88*dd5b3ca6SJunchao Zhang nroots = 1; 89*dd5b3ca6SJunchao Zhang nleaves = 1; 90*dd5b3ca6SJunchao Zhang ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); 91*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 92*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 93*dd5b3ca6SJunchao Zhang ilocal[0] = rank; 94*dd5b3ca6SJunchao Zhang iremote[0].rank = 0; /* rank in PETSC_COMM_SELF */ 95*dd5b3ca6SJunchao Zhang iremote[0].index = rank; /* LocalSF is an embedded SF. Indices are not remapped */ 96*dd5b3ca6SJunchao Zhang 97*dd5b3ca6SJunchao Zhang ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr); 98*dd5b3ca6SJunchao Zhang ierr = PetscSFSetGraph(lsf,nroots,nleaves,NULL/*contiguous leaves*/,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 99*dd5b3ca6SJunchao Zhang ierr = PetscSFSetUp(lsf);CHKERRQ(ierr); 100*dd5b3ca6SJunchao Zhang *out = lsf; 101*dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 102*dd5b3ca6SJunchao Zhang } 103*dd5b3ca6SJunchao Zhang 104*dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFCreateEmbeddedSF_Alltoall(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf) 105*dd5b3ca6SJunchao Zhang { 106*dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 107*dd5b3ca6SJunchao Zhang PetscInt i,*tmproots,*ilocal,ndranks,ndiranks; 108*dd5b3ca6SJunchao Zhang PetscSFNode *iremote; 109*dd5b3ca6SJunchao Zhang PetscMPIInt nroots,*roots,nleaves,*leaves,rank; 110*dd5b3ca6SJunchao Zhang MPI_Comm comm; 111*dd5b3ca6SJunchao Zhang PetscSF_Basic *bas; 112*dd5b3ca6SJunchao Zhang PetscSF esf; 113*dd5b3ca6SJunchao Zhang 114*dd5b3ca6SJunchao Zhang PetscFunctionBegin; 115*dd5b3ca6SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 116*dd5b3ca6SJunchao Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 117*dd5b3ca6SJunchao Zhang 118*dd5b3ca6SJunchao Zhang /* Uniq selected[] and store the result in roots[] */ 119*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nselected,&tmproots);CHKERRQ(ierr); 120*dd5b3ca6SJunchao Zhang ierr = PetscMemcpy(tmproots,selected,nselected*sizeof(PetscInt));CHKERRQ(ierr); 121*dd5b3ca6SJunchao Zhang ierr = PetscSortRemoveDupsInt(&nselected,tmproots);CHKERRQ(ierr); /* nselected might be changed */ 122*dd5b3ca6SJunchao Zhang if (tmproots[0] < 0 || tmproots[nselected-1] >= sf->nroots) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max root indices %D/%D are not in [0,%D)",tmproots[0],tmproots[nselected-1],sf->nroots); 123*dd5b3ca6SJunchao Zhang nroots = nselected; /* For Alltoall, we know root indices will not overflow MPI_INT */ 124*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nselected,&roots);CHKERRQ(ierr); 125*dd5b3ca6SJunchao Zhang for (i=0; i<nselected; i++) roots[i] = tmproots[i]; 126*dd5b3ca6SJunchao Zhang ierr = PetscFree(tmproots);CHKERRQ(ierr); 127*dd5b3ca6SJunchao Zhang 128*dd5b3ca6SJunchao Zhang /* Find out which leaves are still connected to roots in the embedded sf. Expect PetscCommBuildTwoSided is more scalable than MPI_Alltoall */ 129*dd5b3ca6SJunchao Zhang ierr = PetscCommBuildTwoSided(comm,0/*empty msg*/,MPI_INT/*fake*/,nroots,roots,NULL/*todata*/,&nleaves,&leaves,NULL/*fromdata*/);CHKERRQ(ierr); 130*dd5b3ca6SJunchao Zhang 131*dd5b3ca6SJunchao Zhang /* Move myself ahead if rank is in leaves[], since I am a distinguished rank */ 132*dd5b3ca6SJunchao Zhang ndranks = 0; 133*dd5b3ca6SJunchao Zhang for (i=0; i<nleaves; i++) { 134*dd5b3ca6SJunchao Zhang if (leaves[i] == rank) {leaves[i] = -rank; ndranks = 1; break;} 135*dd5b3ca6SJunchao Zhang } 136*dd5b3ca6SJunchao Zhang ierr = PetscSortMPIInt(nleaves,leaves);CHKERRQ(ierr); 137*dd5b3ca6SJunchao Zhang if (nleaves && leaves[0] < 0) leaves[0] = rank; 138*dd5b3ca6SJunchao Zhang 139*dd5b3ca6SJunchao Zhang /* Build esf and fill its fields manually (without calling PetscSFSetUp) */ 140*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 141*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 142*dd5b3ca6SJunchao Zhang for (i=0; i<nleaves; i++) { /* 1:1 map from roots to leaves */ 143*dd5b3ca6SJunchao Zhang ilocal[i] = leaves[i]; 144*dd5b3ca6SJunchao Zhang iremote[i].rank = leaves[i]; 145*dd5b3ca6SJunchao Zhang iremote[i].index = leaves[i]; 146*dd5b3ca6SJunchao Zhang } 147*dd5b3ca6SJunchao Zhang ierr = PetscSFCreate(comm,&esf);CHKERRQ(ierr); 148*dd5b3ca6SJunchao Zhang ierr = PetscSFSetType(esf,PETSCSFBASIC);CHKERRQ(ierr); /* This optimized routine can only create a basic sf */ 149*dd5b3ca6SJunchao Zhang ierr = PetscSFSetGraph(esf,sf->nleaves,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 150*dd5b3ca6SJunchao Zhang 151*dd5b3ca6SJunchao Zhang /* As if we are calling PetscSFSetUpRanks(esf,self's group) */ 152*dd5b3ca6SJunchao Zhang ierr = PetscMalloc4(nleaves,&esf->ranks,nleaves+1,&esf->roffset,nleaves,&esf->rmine,nleaves,&esf->rremote);CHKERRQ(ierr); 153*dd5b3ca6SJunchao Zhang esf->nranks = nleaves; 154*dd5b3ca6SJunchao Zhang esf->ndranks = ndranks; 155*dd5b3ca6SJunchao Zhang esf->roffset[0] = 0; 156*dd5b3ca6SJunchao Zhang for (i=0; i<nleaves; i++) { 157*dd5b3ca6SJunchao Zhang esf->ranks[i] = leaves[i]; 158*dd5b3ca6SJunchao Zhang esf->roffset[i+1] = i+1; 159*dd5b3ca6SJunchao Zhang esf->rmine[i] = leaves[i]; 160*dd5b3ca6SJunchao Zhang esf->rremote[i] = leaves[i]; 161*dd5b3ca6SJunchao Zhang } 162*dd5b3ca6SJunchao Zhang 163*dd5b3ca6SJunchao Zhang /* Set up esf->data, the incoming communication (i.e., recv info), which is usually done by PetscSFSetUp_Basic */ 164*dd5b3ca6SJunchao Zhang bas = (PetscSF_Basic*)esf->data; 165*dd5b3ca6SJunchao Zhang ierr = PetscMalloc2(nroots,&bas->iranks,nroots+1,&bas->ioffset);CHKERRQ(ierr); 166*dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nroots,&bas->irootloc);CHKERRQ(ierr); 167*dd5b3ca6SJunchao Zhang /* Move myself ahead if rank is in roots[], since I am a distinguished irank */ 168*dd5b3ca6SJunchao Zhang ndiranks = 0; 169*dd5b3ca6SJunchao Zhang for (i=0; i<nroots; i++) { 170*dd5b3ca6SJunchao Zhang if (roots[i] == rank) {roots[i] = -rank; ndiranks = 1; break;} 171*dd5b3ca6SJunchao Zhang } 172*dd5b3ca6SJunchao Zhang ierr = PetscSortMPIInt(nroots,roots);CHKERRQ(ierr); 173*dd5b3ca6SJunchao Zhang if (nroots && roots[0] < 0) roots[0] = rank; 174*dd5b3ca6SJunchao Zhang 175*dd5b3ca6SJunchao Zhang bas->niranks = nroots; 176*dd5b3ca6SJunchao Zhang bas->ndiranks = ndiranks; 177*dd5b3ca6SJunchao Zhang bas->ioffset[0] = 0; 178*dd5b3ca6SJunchao Zhang bas->itotal = nroots; 179*dd5b3ca6SJunchao Zhang for (i=0; i<nroots; i++) { 180*dd5b3ca6SJunchao Zhang bas->iranks[i] = roots[i]; 181*dd5b3ca6SJunchao Zhang bas->ioffset[i+1] = i+1; 182*dd5b3ca6SJunchao Zhang bas->irootloc[i] = roots[i]; 183*dd5b3ca6SJunchao Zhang } 184*dd5b3ca6SJunchao Zhang 185*dd5b3ca6SJunchao Zhang *newsf = esf; 186*dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 187*dd5b3ca6SJunchao Zhang } 188*dd5b3ca6SJunchao Zhang 189*dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFCreate_Alltoall(PetscSF sf) 190*dd5b3ca6SJunchao Zhang { 191*dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 192*dd5b3ca6SJunchao Zhang PetscSF_Alltoall *dat = (PetscSF_Alltoall*)sf->data; 193*dd5b3ca6SJunchao Zhang 194*dd5b3ca6SJunchao Zhang PetscFunctionBegin; 195*dd5b3ca6SJunchao Zhang /* Inherit from Allgatherv. It is astonishing Alltoall can inherit so much from Allgather(v) */ 196*dd5b3ca6SJunchao Zhang sf->ops->Destroy = PetscSFDestroy_Allgatherv; 197*dd5b3ca6SJunchao Zhang sf->ops->Reset = PetscSFReset_Allgatherv; 198*dd5b3ca6SJunchao Zhang sf->ops->BcastAndOpEnd = PetscSFBcastAndOpEnd_Allgatherv; 199*dd5b3ca6SJunchao Zhang sf->ops->ReduceEnd = PetscSFReduceEnd_Allgatherv; 200*dd5b3ca6SJunchao Zhang sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Allgatherv; 201*dd5b3ca6SJunchao Zhang sf->ops->GetRootRanks = PetscSFGetRootRanks_Allgatherv; 202*dd5b3ca6SJunchao Zhang 203*dd5b3ca6SJunchao Zhang /* Inherit from Allgather. Every process gathers equal-sized data from others, which enables this inheritance. */ 204*dd5b3ca6SJunchao Zhang sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Allgatherv; 205*dd5b3ca6SJunchao Zhang 206*dd5b3ca6SJunchao Zhang /* Inherit from Gatherv. Each root has only one leaf connected, which enables this inheritance */ 207*dd5b3ca6SJunchao Zhang sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Gatherv; 208*dd5b3ca6SJunchao Zhang 209*dd5b3ca6SJunchao Zhang /* Alltoall stuff */ 210*dd5b3ca6SJunchao Zhang sf->ops->GetGraph = PetscSFGetGraph_Alltoall; 211*dd5b3ca6SJunchao Zhang sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Alltoall; 212*dd5b3ca6SJunchao Zhang sf->ops->ReduceBegin = PetscSFReduceBegin_Alltoall; 213*dd5b3ca6SJunchao Zhang sf->ops->CreateLocalSF = PetscSFCreateLocalSF_Alltoall; 214*dd5b3ca6SJunchao Zhang sf->ops->CreateEmbeddedSF = PetscSFCreateEmbeddedSF_Alltoall; 215*dd5b3ca6SJunchao Zhang 216*dd5b3ca6SJunchao Zhang ierr = PetscNewLog(sf,&dat);CHKERRQ(ierr); 217*dd5b3ca6SJunchao Zhang sf->data = (void*)dat; 218*dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 219*dd5b3ca6SJunchao Zhang } 220*dd5b3ca6SJunchao Zhang 221*dd5b3ca6SJunchao Zhang 222