1dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h> 2dd5b3ca6SJunchao Zhang 3eb02082bSJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFBcastAndOpBegin_Gatherv(PetscSF,MPI_Datatype,PetscMemType,const void*,PetscMemType,void*,MPI_Op); 4dd5b3ca6SJunchao Zhang 5dd5b3ca6SJunchao Zhang /* PetscSFGetGraph is non-collective. An implementation should not have collective calls */ 6dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFGetGraph_Allgatherv(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote) 7dd5b3ca6SJunchao Zhang { 8dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 9dd5b3ca6SJunchao Zhang PetscInt i,j,k; 10dd5b3ca6SJunchao Zhang const PetscInt *range; 11dd5b3ca6SJunchao Zhang PetscMPIInt size; 12dd5b3ca6SJunchao Zhang 13dd5b3ca6SJunchao Zhang PetscFunctionBegin; 14ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRMPI(ierr); 15dd5b3ca6SJunchao Zhang if (nroots) *nroots = sf->nroots; 16dd5b3ca6SJunchao Zhang if (nleaves) *nleaves = sf->nleaves; 17dd5b3ca6SJunchao Zhang if (ilocal) *ilocal = NULL; /* Contiguous leaves */ 18dd5b3ca6SJunchao Zhang if (iremote) { 19dd5b3ca6SJunchao Zhang if (!sf->remote && sf->nleaves) { /* The && sf->nleaves makes sfgatherv able to inherit this routine */ 20dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr); 21dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(sf->nleaves,&sf->remote);CHKERRQ(ierr); 22dd5b3ca6SJunchao Zhang sf->remote_alloc = sf->remote; 23dd5b3ca6SJunchao Zhang for (i=0; i<size; i++) { 24dd5b3ca6SJunchao Zhang for (j=range[i],k=0; j<range[i+1]; j++,k++) { 25dd5b3ca6SJunchao Zhang sf->remote[j].rank = i; 26dd5b3ca6SJunchao Zhang sf->remote[j].index = k; 27dd5b3ca6SJunchao Zhang } 28dd5b3ca6SJunchao Zhang } 29dd5b3ca6SJunchao Zhang } 30dd5b3ca6SJunchao Zhang *iremote = sf->remote; 31dd5b3ca6SJunchao Zhang } 32dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 33dd5b3ca6SJunchao Zhang } 34dd5b3ca6SJunchao Zhang 35dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFSetUp_Allgatherv(PetscSF sf) 36dd5b3ca6SJunchao Zhang { 37dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 38dd5b3ca6SJunchao Zhang PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data; 39dd5b3ca6SJunchao Zhang PetscMPIInt size; 40dd5b3ca6SJunchao Zhang PetscInt i; 41dd5b3ca6SJunchao Zhang const PetscInt *range; 42dd5b3ca6SJunchao Zhang 43dd5b3ca6SJunchao Zhang PetscFunctionBegin; 44cd620004SJunchao Zhang ierr = PetscSFSetUp_Allgather(sf);CHKERRQ(ierr); 45ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRMPI(ierr); 46dd5b3ca6SJunchao Zhang if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */ 47dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(size,&dat->recvcounts);CHKERRQ(ierr); 48dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(size,&dat->displs);CHKERRQ(ierr); 49dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr); 50dd5b3ca6SJunchao Zhang 51dd5b3ca6SJunchao Zhang for (i=0; i<size; i++) { 52dd5b3ca6SJunchao Zhang ierr = PetscMPIIntCast(range[i],&dat->displs[i]);CHKERRQ(ierr); 53dd5b3ca6SJunchao Zhang ierr = PetscMPIIntCast(range[i+1]-range[i],&dat->recvcounts[i]);CHKERRQ(ierr); 54dd5b3ca6SJunchao Zhang } 55dd5b3ca6SJunchao Zhang } 56dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 57dd5b3ca6SJunchao Zhang } 58dd5b3ca6SJunchao Zhang 59dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf) 60dd5b3ca6SJunchao Zhang { 61dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 62eb02082bSJunchao Zhang PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data; 63dd5b3ca6SJunchao Zhang 64dd5b3ca6SJunchao Zhang PetscFunctionBegin; 65dd5b3ca6SJunchao Zhang ierr = PetscFree(dat->iranks);CHKERRQ(ierr); 66dd5b3ca6SJunchao Zhang ierr = PetscFree(dat->ioffset);CHKERRQ(ierr); 67dd5b3ca6SJunchao Zhang ierr = PetscFree(dat->irootloc);CHKERRQ(ierr); 68dd5b3ca6SJunchao Zhang ierr = PetscFree(dat->recvcounts);CHKERRQ(ierr); 69dd5b3ca6SJunchao Zhang ierr = PetscFree(dat->displs);CHKERRQ(ierr); 70dd5b3ca6SJunchao Zhang if (dat->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed"); 71cd620004SJunchao Zhang ierr = PetscSFLinkDestroy(sf,&dat->avail);CHKERRQ(ierr); 72dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 73dd5b3ca6SJunchao Zhang } 74dd5b3ca6SJunchao Zhang 75dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf) 76dd5b3ca6SJunchao Zhang { 77dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 78dd5b3ca6SJunchao Zhang 79dd5b3ca6SJunchao Zhang PetscFunctionBegin; 80dd5b3ca6SJunchao Zhang ierr = PetscSFReset_Allgatherv(sf);CHKERRQ(ierr); 81dd5b3ca6SJunchao Zhang ierr = PetscFree(sf->data);CHKERRQ(ierr); 82dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 83dd5b3ca6SJunchao Zhang } 84dd5b3ca6SJunchao Zhang 85eb02082bSJunchao Zhang static PetscErrorCode PetscSFBcastAndOpBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op) 86dd5b3ca6SJunchao Zhang { 87dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 88cd620004SJunchao Zhang PetscSFLink link; 89dd5b3ca6SJunchao Zhang PetscMPIInt sendcount; 90dd5b3ca6SJunchao Zhang MPI_Comm comm; 91cd620004SJunchao Zhang void *rootbuf = NULL,*leafbuf = NULL; 92cd620004SJunchao Zhang MPI_Request *req; 93dd5b3ca6SJunchao Zhang PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data; 94dd5b3ca6SJunchao Zhang 95dd5b3ca6SJunchao Zhang PetscFunctionBegin; 96cd620004SJunchao Zhang ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_BCAST,&link);CHKERRQ(ierr); 97cd620004SJunchao Zhang ierr = PetscSFLinkPackRootData(sf,link,PETSCSF_REMOTE,rootdata);CHKERRQ(ierr); 98dd5b3ca6SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 99dd5b3ca6SJunchao Zhang ierr = PetscMPIIntCast(sf->nroots,&sendcount);CHKERRQ(ierr); 100cd620004SJunchao Zhang ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_ROOT2LEAF,&rootbuf,&leafbuf,&req,NULL);CHKERRQ(ierr); 101cd620004SJunchao Zhang ierr = MPIU_Iallgatherv(rootbuf,sendcount,unit,leafbuf,dat->recvcounts,dat->displs,unit,comm,req);CHKERRQ(ierr); 102855db38dSJunchao Zhang PetscFunctionReturn(0); 103855db38dSJunchao Zhang } 104855db38dSJunchao Zhang 105eb02082bSJunchao Zhang static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op) 106eb02082bSJunchao Zhang { 107eb02082bSJunchao Zhang PetscErrorCode ierr; 108cd620004SJunchao Zhang PetscSFLink link; 109dd5b3ca6SJunchao Zhang PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data; 110dd5b3ca6SJunchao Zhang PetscInt rstart; 111cd620004SJunchao Zhang PetscMPIInt rank,count,recvcount; 112dd5b3ca6SJunchao Zhang MPI_Comm comm; 113cd620004SJunchao Zhang void *rootbuf = NULL,*leafbuf = NULL; 114cd620004SJunchao Zhang MPI_Request *req; 115dd5b3ca6SJunchao Zhang 116dd5b3ca6SJunchao Zhang PetscFunctionBegin; 117cd620004SJunchao Zhang ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_REDUCE,&link);CHKERRQ(ierr); 118dd5b3ca6SJunchao Zhang if (op == MPIU_REPLACE) { 119cd620004SJunchao Zhang /* REPLACE is only meaningful when all processes have the same leafdata to reduce. Therefore copying from local leafdata is fine */ 120dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetRange(sf->map,&rstart,NULL);CHKERRQ(ierr); 12120c24465SJunchao Zhang ierr = (*link->Memcpy)(link,rootmtype,rootdata,leafmtype,(const char*)leafdata+(size_t)rstart*link->unitbytes,(size_t)sf->nroots*link->unitbytes);CHKERRQ(ierr); 122*e7a27f33SJunchao Zhang #if defined(PETSC_HAVE_DEVICE) 123*e7a27f33SJunchao Zhang if (PetscMemTypeHost(rootmtype) && PetscMemTypeDevice(leafmtype)) {ierr = (*link->d_SyncStream)(link);CHKERRQ(ierr);} 124*e7a27f33SJunchao Zhang #endif 125dd5b3ca6SJunchao Zhang } else { 126cd620004SJunchao Zhang /* Reduce leafdata, then scatter to rootdata */ 127cd620004SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 128ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 129cd620004SJunchao Zhang ierr = PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);CHKERRQ(ierr); 130cd620004SJunchao Zhang ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_LEAF2ROOT,&rootbuf,&leafbuf,&req,NULL);CHKERRQ(ierr); 131cd620004SJunchao Zhang ierr = PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE],&recvcount);CHKERRQ(ierr); 132cd620004SJunchao Zhang /* Allocate a separate leaf buffer on rank 0 */ 133cd620004SJunchao Zhang if (!rank && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) { 13420c24465SJunchao Zhang ierr = PetscSFMalloc(sf,link->leafmtype_mpi,sf->leafbuflen[PETSCSF_REMOTE]*link->unitbytes,(void**)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]);CHKERRQ(ierr); 135dd5b3ca6SJunchao Zhang } 136cd620004SJunchao Zhang /* In case we already copied leafdata from device to host (i.e., no use_gpu_aware_mpi), we need to adjust leafbuf on rank 0 */ 137cd620004SJunchao Zhang if (!rank && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE; 138cd620004SJunchao Zhang ierr = PetscMPIIntCast(sf->nleaves*link->bs,&count);CHKERRQ(ierr); 139ffc4695bSBarry Smith ierr = MPI_Reduce(leafbuf,link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi],count,link->basicunit,op,0,comm);CHKERRMPI(ierr); 140cd620004SJunchao Zhang ierr = MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi],dat->recvcounts,dat->displs,unit,rootbuf,recvcount,unit,0,comm,req);CHKERRQ(ierr); 141dd5b3ca6SJunchao Zhang } 142eb02082bSJunchao Zhang PetscFunctionReturn(0); 143eb02082bSJunchao Zhang } 144eb02082bSJunchao Zhang 145eb02082bSJunchao Zhang static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata) 146eb02082bSJunchao Zhang { 147eb02082bSJunchao Zhang PetscErrorCode ierr; 148cd620004SJunchao Zhang PetscSFLink link; 149855db38dSJunchao Zhang PetscMPIInt rank; 150eb02082bSJunchao Zhang 151eb02082bSJunchao Zhang PetscFunctionBegin; 152eb02082bSJunchao Zhang ierr = PetscSFBcastAndOpBegin_Gatherv(sf,unit,rootmtype,rootdata,leafmtype,leafdata,MPIU_REPLACE);CHKERRQ(ierr); 153cd620004SJunchao Zhang ierr = PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr); 154cd620004SJunchao Zhang ierr = PetscSFLinkMPIWaitall(sf,link,PETSCSF_ROOT2LEAF);CHKERRQ(ierr); 155ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr); 156c2a741eeSJunchao Zhang if (!rank && leafmtype == PETSC_MEMTYPE_DEVICE && !sf->use_gpu_aware_mpi) { 15720c24465SJunchao Zhang ierr = (*link->Memcpy)(link,PETSC_MEMTYPE_DEVICE,leafdata,PETSC_MEMTYPE_HOST,link->leafbuf[PETSC_MEMTYPE_HOST],sf->leafbuflen[PETSCSF_REMOTE]*link->unitbytes);CHKERRQ(ierr); 158855db38dSJunchao Zhang } 159cd620004SJunchao Zhang ierr = PetscSFLinkReclaim(sf,&link);CHKERRQ(ierr); 160dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 161dd5b3ca6SJunchao Zhang } 162dd5b3ca6SJunchao Zhang 163dd5b3ca6SJunchao Zhang /* This routine is very tricky (I believe it is rarely used with this kind of graph so just provide a simple but not-optimal implementation). 164dd5b3ca6SJunchao Zhang 165dd5b3ca6SJunchao Zhang Suppose we have three ranks. Rank 0 has a root with value 1. Rank 0,1,2 has a leaf with value 2,3,4 respectively. The leaves are connected 166dd5b3ca6SJunchao Zhang to the root on rank 0. Suppose op=MPI_SUM and rank 0,1,2 gets root state in their rank order. By definition of this routine, rank 0 sees 1 167dd5b3ca6SJunchao Zhang in root, fetches it into its leafupate, then updates root to 1 + 2 = 3; rank 1 sees 3 in root, fetches it into its leafupate, then updates 168dd5b3ca6SJunchao Zhang root to 3 + 3 = 6; rank 2 sees 6 in root, fetches it into its leafupdate, then updates root to 6 + 4 = 10. At the end, leafupdate on rank 169dd5b3ca6SJunchao Zhang 0,1,2 is 1,3,6 respectively. root is 10. 170dd5b3ca6SJunchao Zhang 171dd5b3ca6SJunchao Zhang We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate 172dd5b3ca6SJunchao Zhang rank-0 rank-1 rank-2 173dd5b3ca6SJunchao Zhang Root 1 174dd5b3ca6SJunchao Zhang Leaf 2 3 4 175dd5b3ca6SJunchao Zhang Leafupdate 2 3 4 176dd5b3ca6SJunchao Zhang 177dd5b3ca6SJunchao Zhang Do MPI_Exscan on leafupdate, 178dd5b3ca6SJunchao Zhang rank-0 rank-1 rank-2 179dd5b3ca6SJunchao Zhang Root 1 180dd5b3ca6SJunchao Zhang Leaf 2 3 4 181dd5b3ca6SJunchao Zhang Leafupdate 2 2 5 182dd5b3ca6SJunchao Zhang 183dd5b3ca6SJunchao Zhang BcastAndOp from root to leafupdate, 184dd5b3ca6SJunchao Zhang rank-0 rank-1 rank-2 185dd5b3ca6SJunchao Zhang Root 1 186dd5b3ca6SJunchao Zhang Leaf 2 3 4 187dd5b3ca6SJunchao Zhang Leafupdate 3 3 6 188dd5b3ca6SJunchao Zhang 189dd5b3ca6SJunchao Zhang Copy root to leafupdate on rank-0 190dd5b3ca6SJunchao Zhang rank-0 rank-1 rank-2 191dd5b3ca6SJunchao Zhang Root 1 192dd5b3ca6SJunchao Zhang Leaf 2 3 4 193dd5b3ca6SJunchao Zhang Leafupdate 1 3 6 194dd5b3ca6SJunchao Zhang 195dd5b3ca6SJunchao Zhang Reduce from leaf to root, 196dd5b3ca6SJunchao Zhang rank-0 rank-1 rank-2 197dd5b3ca6SJunchao Zhang Root 10 198dd5b3ca6SJunchao Zhang Leaf 2 3 4 199dd5b3ca6SJunchao Zhang Leafupdate 1 3 6 200dd5b3ca6SJunchao Zhang */ 201eb02082bSJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,void *rootdata,PetscMemType leafmtype,const void *leafdata,void *leafupdate,MPI_Op op) 202dd5b3ca6SJunchao Zhang { 203dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 204cd620004SJunchao Zhang PetscSFLink link; 205dd5b3ca6SJunchao Zhang MPI_Comm comm; 206dd5b3ca6SJunchao Zhang PetscMPIInt count; 207dd5b3ca6SJunchao Zhang 208dd5b3ca6SJunchao Zhang PetscFunctionBegin; 209855db38dSJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 210c2a741eeSJunchao Zhang if (!sf->use_gpu_aware_mpi && (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE)) SETERRQ(comm,PETSC_ERR_SUP,"Do FetchAndOp on device but not use gpu-aware MPI"); 211dd5b3ca6SJunchao Zhang /* Copy leafdata to leafupdate */ 212cd620004SJunchao Zhang ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_FETCH,&link);CHKERRQ(ierr); 213cd620004SJunchao Zhang ierr = PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);CHKERRQ(ierr); /* Sync the device */ 21420c24465SJunchao Zhang ierr = (*link->Memcpy)(link,leafmtype,leafupdate,leafmtype,leafdata,sf->nleaves*link->unitbytes);CHKERRQ(ierr); 215cd620004SJunchao Zhang ierr = PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);CHKERRQ(ierr); 216dd5b3ca6SJunchao Zhang 217dd5b3ca6SJunchao Zhang /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */ 218dd5b3ca6SJunchao Zhang if (op == MPIU_REPLACE) { 219dd5b3ca6SJunchao Zhang PetscMPIInt size,rank,prev,next; 220ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 221ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 222dd5b3ca6SJunchao Zhang prev = rank ? rank-1 : MPI_PROC_NULL; 223dd5b3ca6SJunchao Zhang next = (rank < size-1) ? rank+1 : MPI_PROC_NULL; 224cd620004SJunchao Zhang ierr = PetscMPIIntCast(sf->nleaves,&count);CHKERRQ(ierr); 225ffc4695bSBarry Smith ierr = MPI_Sendrecv_replace(leafupdate,count,unit,next,link->tag,prev,link->tag,comm,MPI_STATUSES_IGNORE);CHKERRMPI(ierr); 226cd620004SJunchao Zhang } else { 227cd620004SJunchao Zhang ierr = PetscMPIIntCast(sf->nleaves*link->bs,&count);CHKERRQ(ierr); 228ffc4695bSBarry Smith ierr = MPI_Exscan(MPI_IN_PLACE,leafupdate,count,link->basicunit,op,comm);CHKERRMPI(ierr); 229cd620004SJunchao Zhang } 230cd620004SJunchao Zhang ierr = PetscSFLinkReclaim(sf,&link);CHKERRQ(ierr); 231dd5b3ca6SJunchao Zhang ierr = PetscSFBcastAndOpBegin(sf,unit,rootdata,leafupdate,op);CHKERRQ(ierr); 232dd5b3ca6SJunchao Zhang ierr = PetscSFBcastAndOpEnd(sf,unit,rootdata,leafupdate,op);CHKERRQ(ierr); 233dd5b3ca6SJunchao Zhang 234dd5b3ca6SJunchao Zhang /* Bcast roots to rank 0's leafupdate */ 235dd5b3ca6SJunchao Zhang ierr = PetscSFBcastToZero_Private(sf,unit,rootdata,leafupdate);CHKERRQ(ierr); /* Using this line makes Allgather SFs able to inherit this routine */ 236dd5b3ca6SJunchao Zhang 237dd5b3ca6SJunchao Zhang /* Reduce leafdata to rootdata */ 238dd5b3ca6SJunchao Zhang ierr = PetscSFReduceBegin(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 239dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 240dd5b3ca6SJunchao Zhang } 241dd5b3ca6SJunchao Zhang 24200816365SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op) 243dd5b3ca6SJunchao Zhang { 244dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 245dd5b3ca6SJunchao Zhang 246dd5b3ca6SJunchao Zhang PetscFunctionBegin; 247dd5b3ca6SJunchao Zhang ierr = PetscSFReduceEnd(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr); 248dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 249dd5b3ca6SJunchao Zhang } 250dd5b3ca6SJunchao Zhang 251dd5b3ca6SJunchao Zhang /* Get root ranks accessing my leaves */ 252dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote) 253dd5b3ca6SJunchao Zhang { 254dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 255dd5b3ca6SJunchao Zhang PetscInt i,j,k,size; 256dd5b3ca6SJunchao Zhang const PetscInt *range; 257dd5b3ca6SJunchao Zhang 258dd5b3ca6SJunchao Zhang PetscFunctionBegin; 259dd5b3ca6SJunchao Zhang /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */ 260dd5b3ca6SJunchao Zhang if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */ 261dd5b3ca6SJunchao Zhang size = sf->nranks; 262dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr); 263dd5b3ca6SJunchao Zhang ierr = PetscMalloc4(size,&sf->ranks,size+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr); 264dd5b3ca6SJunchao Zhang for (i=0; i<size; i++) sf->ranks[i] = i; 265da2e4c71SJunchao Zhang ierr = PetscArraycpy(sf->roffset,range,size+1);CHKERRQ(ierr); 266dd5b3ca6SJunchao Zhang for (i=0; i<sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */ 267dd5b3ca6SJunchao Zhang for (i=0; i<size; i++) { 268dd5b3ca6SJunchao Zhang for (j=range[i],k=0; j<range[i+1]; j++,k++) sf->rremote[j] = k; 269dd5b3ca6SJunchao Zhang } 270dd5b3ca6SJunchao Zhang } 271dd5b3ca6SJunchao Zhang 272dd5b3ca6SJunchao Zhang if (nranks) *nranks = sf->nranks; 273dd5b3ca6SJunchao Zhang if (ranks) *ranks = sf->ranks; 274dd5b3ca6SJunchao Zhang if (roffset) *roffset = sf->roffset; 275dd5b3ca6SJunchao Zhang if (rmine) *rmine = sf->rmine; 276dd5b3ca6SJunchao Zhang if (rremote) *rremote = sf->rremote; 277dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 278dd5b3ca6SJunchao Zhang } 279dd5b3ca6SJunchao Zhang 280dd5b3ca6SJunchao Zhang /* Get leaf ranks accessing my roots */ 281dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc) 282dd5b3ca6SJunchao Zhang { 283dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 284dd5b3ca6SJunchao Zhang PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data; 285dd5b3ca6SJunchao Zhang MPI_Comm comm; 286dd5b3ca6SJunchao Zhang PetscMPIInt size,rank; 287dd5b3ca6SJunchao Zhang PetscInt i,j; 288dd5b3ca6SJunchao Zhang 289dd5b3ca6SJunchao Zhang PetscFunctionBegin; 290dd5b3ca6SJunchao Zhang /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */ 291dd5b3ca6SJunchao Zhang ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 292ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 293ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 294dd5b3ca6SJunchao Zhang if (niranks) *niranks = size; 295dd5b3ca6SJunchao Zhang 296dd5b3ca6SJunchao Zhang /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and 297dd5b3ca6SJunchao Zhang sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why. 298dd5b3ca6SJunchao Zhang */ 299dd5b3ca6SJunchao Zhang if (iranks) { 300dd5b3ca6SJunchao Zhang if (!dat->iranks) { 301dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(size,&dat->iranks);CHKERRQ(ierr); 302dd5b3ca6SJunchao Zhang dat->iranks[0] = rank; 303dd5b3ca6SJunchao Zhang for (i=0,j=1; i<size; i++) {if (i == rank) continue; dat->iranks[j++] = i;} 304dd5b3ca6SJunchao Zhang } 305dd5b3ca6SJunchao Zhang *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNewLog */ 306dd5b3ca6SJunchao Zhang } 307dd5b3ca6SJunchao Zhang 308dd5b3ca6SJunchao Zhang if (ioffset) { 309dd5b3ca6SJunchao Zhang if (!dat->ioffset) { 310dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(size+1,&dat->ioffset);CHKERRQ(ierr); 311dd5b3ca6SJunchao Zhang for (i=0; i<=size; i++) dat->ioffset[i] = i*sf->nroots; 312dd5b3ca6SJunchao Zhang } 313dd5b3ca6SJunchao Zhang *ioffset = dat->ioffset; 314dd5b3ca6SJunchao Zhang } 315dd5b3ca6SJunchao Zhang 316dd5b3ca6SJunchao Zhang if (irootloc) { 317dd5b3ca6SJunchao Zhang if (!dat->irootloc) { 318dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(sf->nleaves,&dat->irootloc);CHKERRQ(ierr); 319dd5b3ca6SJunchao Zhang for (i=0; i<size; i++) { 320dd5b3ca6SJunchao Zhang for (j=0; j<sf->nroots; j++) dat->irootloc[i*sf->nroots+j] = j; 321dd5b3ca6SJunchao Zhang } 322dd5b3ca6SJunchao Zhang } 323dd5b3ca6SJunchao Zhang *irootloc = dat->irootloc; 324dd5b3ca6SJunchao Zhang } 325dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 326dd5b3ca6SJunchao Zhang } 327dd5b3ca6SJunchao Zhang 328dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf,PetscSF *out) 329dd5b3ca6SJunchao Zhang { 330dd5b3ca6SJunchao Zhang PetscInt i,nroots,nleaves,rstart,*ilocal; 331dd5b3ca6SJunchao Zhang PetscSFNode *iremote; 332dd5b3ca6SJunchao Zhang PetscSF lsf; 333dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 334dd5b3ca6SJunchao Zhang 335dd5b3ca6SJunchao Zhang PetscFunctionBegin; 336eb02082bSJunchao Zhang nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */ 337eb02082bSJunchao Zhang nroots = nleaves; 338dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 339dd5b3ca6SJunchao Zhang ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 340dd5b3ca6SJunchao Zhang ierr = PetscLayoutGetRange(sf->map,&rstart,NULL);CHKERRQ(ierr); 341dd5b3ca6SJunchao Zhang 342dd5b3ca6SJunchao Zhang for (i=0; i<nleaves; i++) { 343dd5b3ca6SJunchao Zhang ilocal[i] = rstart + i; /* lsf does not change leave indices */ 344dd5b3ca6SJunchao Zhang iremote[i].rank = 0; /* rank in PETSC_COMM_SELF */ 345dd5b3ca6SJunchao Zhang iremote[i].index = i; /* root index */ 346dd5b3ca6SJunchao Zhang } 347dd5b3ca6SJunchao Zhang 348dd5b3ca6SJunchao Zhang ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr); 349dd5b3ca6SJunchao Zhang ierr = PetscSFSetGraph(lsf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 350dd5b3ca6SJunchao Zhang ierr = PetscSFSetUp(lsf);CHKERRQ(ierr); 351dd5b3ca6SJunchao Zhang *out = lsf; 352dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 353dd5b3ca6SJunchao Zhang } 354dd5b3ca6SJunchao Zhang 355dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf) 356dd5b3ca6SJunchao Zhang { 357dd5b3ca6SJunchao Zhang PetscErrorCode ierr; 358dd5b3ca6SJunchao Zhang PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data; 359dd5b3ca6SJunchao Zhang 360dd5b3ca6SJunchao Zhang PetscFunctionBegin; 361cd620004SJunchao Zhang sf->ops->BcastAndOpEnd = PetscSFBcastAndOpEnd_Basic; 362cd620004SJunchao Zhang sf->ops->ReduceEnd = PetscSFReduceEnd_Basic; 363cd620004SJunchao Zhang 364dd5b3ca6SJunchao Zhang sf->ops->SetUp = PetscSFSetUp_Allgatherv; 365dd5b3ca6SJunchao Zhang sf->ops->Reset = PetscSFReset_Allgatherv; 366dd5b3ca6SJunchao Zhang sf->ops->Destroy = PetscSFDestroy_Allgatherv; 367dd5b3ca6SJunchao Zhang sf->ops->GetRootRanks = PetscSFGetRootRanks_Allgatherv; 368dd5b3ca6SJunchao Zhang sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Allgatherv; 369dd5b3ca6SJunchao Zhang sf->ops->GetGraph = PetscSFGetGraph_Allgatherv; 370dd5b3ca6SJunchao Zhang sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Allgatherv; 371dd5b3ca6SJunchao Zhang sf->ops->ReduceBegin = PetscSFReduceBegin_Allgatherv; 372dd5b3ca6SJunchao Zhang sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv; 373dd5b3ca6SJunchao Zhang sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Allgatherv; 374dd5b3ca6SJunchao Zhang sf->ops->CreateLocalSF = PetscSFCreateLocalSF_Allgatherv; 375dd5b3ca6SJunchao Zhang sf->ops->BcastToZero = PetscSFBcastToZero_Allgatherv; 376dd5b3ca6SJunchao Zhang 377dd5b3ca6SJunchao Zhang ierr = PetscNewLog(sf,&dat);CHKERRQ(ierr); 378dd5b3ca6SJunchao Zhang sf->data = (void*)dat; 379dd5b3ca6SJunchao Zhang PetscFunctionReturn(0); 380dd5b3ca6SJunchao Zhang } 381