xref: /petsc/src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h>
2dd5b3ca6SJunchao Zhang 
3ad227feaSJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFBcastBegin_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 */
69371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFGetGraph_Allgatherv(PetscSF sf, PetscInt *nroots, PetscInt *nleaves, const PetscInt **ilocal, const PetscSFNode **iremote) {
7dd5b3ca6SJunchao Zhang   PetscInt        i, j, k;
8dd5b3ca6SJunchao Zhang   const PetscInt *range;
9dd5b3ca6SJunchao Zhang   PetscMPIInt     size;
10dd5b3ca6SJunchao Zhang 
11dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
13dd5b3ca6SJunchao Zhang   if (nroots) *nroots = sf->nroots;
14dd5b3ca6SJunchao Zhang   if (nleaves) *nleaves = sf->nleaves;
15dd5b3ca6SJunchao Zhang   if (ilocal) *ilocal = NULL; /* Contiguous leaves */
16dd5b3ca6SJunchao Zhang   if (iremote) {
17dd5b3ca6SJunchao Zhang     if (!sf->remote && sf->nleaves) { /* The && sf->nleaves makes sfgatherv able to inherit this routine */
189566063dSJacob Faibussowitsch       PetscCall(PetscLayoutGetRanges(sf->map, &range));
199566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sf->nleaves, &sf->remote));
20dd5b3ca6SJunchao Zhang       sf->remote_alloc = sf->remote;
21dd5b3ca6SJunchao Zhang       for (i = 0; i < size; i++) {
22dd5b3ca6SJunchao Zhang         for (j = range[i], k = 0; j < range[i + 1]; j++, k++) {
23dd5b3ca6SJunchao Zhang           sf->remote[j].rank  = i;
24dd5b3ca6SJunchao Zhang           sf->remote[j].index = k;
25dd5b3ca6SJunchao Zhang         }
26dd5b3ca6SJunchao Zhang       }
27dd5b3ca6SJunchao Zhang     }
28dd5b3ca6SJunchao Zhang     *iremote = sf->remote;
29dd5b3ca6SJunchao Zhang   }
30dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
31dd5b3ca6SJunchao Zhang }
32dd5b3ca6SJunchao Zhang 
339371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFSetUp_Allgatherv(PetscSF sf) {
34dd5b3ca6SJunchao Zhang   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
35dd5b3ca6SJunchao Zhang   PetscMPIInt         size;
36dd5b3ca6SJunchao Zhang   PetscInt            i;
37dd5b3ca6SJunchao Zhang   const PetscInt     *range;
38dd5b3ca6SJunchao Zhang 
39dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp_Allgather(sf));
419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
42dd5b3ca6SJunchao Zhang   if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */
439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &dat->recvcounts));
449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &dat->displs));
459566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(sf->map, &range));
46dd5b3ca6SJunchao Zhang 
47dd5b3ca6SJunchao Zhang     for (i = 0; i < size; i++) {
489566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(range[i], &dat->displs[i]));
499566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(range[i + 1] - range[i], &dat->recvcounts[i]));
50dd5b3ca6SJunchao Zhang     }
51dd5b3ca6SJunchao Zhang   }
52dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
53dd5b3ca6SJunchao Zhang }
54dd5b3ca6SJunchao Zhang 
559371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf) {
56eb02082bSJunchao Zhang   PetscSF_Allgatherv *dat  = (PetscSF_Allgatherv *)sf->data;
5771438e86SJunchao Zhang   PetscSFLink         link = dat->avail, next;
58dd5b3ca6SJunchao Zhang 
59dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(PetscFree(dat->iranks));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(dat->ioffset));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree(dat->irootloc));
639566063dSJacob Faibussowitsch   PetscCall(PetscFree(dat->recvcounts));
649566063dSJacob Faibussowitsch   PetscCall(PetscFree(dat->displs));
6528b400f6SJacob Faibussowitsch   PetscCheck(!dat->inuse, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Outstanding operation has not been completed");
669371c9d4SSatish Balay   for (; link; link = next) {
679371c9d4SSatish Balay     next = link->next;
689371c9d4SSatish Balay     PetscCall(PetscSFLinkDestroy(sf, link));
699371c9d4SSatish Balay   }
7071438e86SJunchao Zhang   dat->avail = NULL;
71dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
72dd5b3ca6SJunchao Zhang }
73dd5b3ca6SJunchao Zhang 
749371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf) {
75dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   PetscCall(PetscSFReset_Allgatherv(sf));
779566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->data));
78dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
79dd5b3ca6SJunchao Zhang }
80dd5b3ca6SJunchao Zhang 
819371c9d4SSatish Balay static PetscErrorCode PetscSFBcastBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op) {
82cd620004SJunchao Zhang   PetscSFLink         link;
83dd5b3ca6SJunchao Zhang   PetscMPIInt         sendcount;
84dd5b3ca6SJunchao Zhang   MPI_Comm            comm;
85cd620004SJunchao Zhang   void               *rootbuf = NULL, *leafbuf = NULL;
86cd620004SJunchao Zhang   MPI_Request        *req;
87dd5b3ca6SJunchao Zhang   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
88dd5b3ca6SJunchao Zhang 
89dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_BCAST, &link));
919566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
929566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
949566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
959566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link, PETSCSF_ROOT2LEAF));
979566063dSJacob Faibussowitsch   PetscCallMPI(MPIU_Iallgatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, comm, req));
98855db38dSJunchao Zhang   PetscFunctionReturn(0);
99855db38dSJunchao Zhang }
100855db38dSJunchao Zhang 
1019371c9d4SSatish Balay static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op) {
102cd620004SJunchao Zhang   PetscSFLink         link;
103dd5b3ca6SJunchao Zhang   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
104dd5b3ca6SJunchao Zhang   PetscInt            rstart;
105cd620004SJunchao Zhang   PetscMPIInt         rank, count, recvcount;
106dd5b3ca6SJunchao Zhang   MPI_Comm            comm;
107cd620004SJunchao Zhang   void               *rootbuf = NULL, *leafbuf = NULL;
108cd620004SJunchao Zhang   MPI_Request        *req;
109dd5b3ca6SJunchao Zhang 
110dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_REDUCE, &link));
11283df288dSJunchao Zhang   if (op == MPI_REPLACE) {
113cd620004SJunchao Zhang     /* REPLACE is only meaningful when all processes have the same leafdata to reduce. Therefore copying from local leafdata is fine */
1149566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
1159566063dSJacob Faibussowitsch     PetscCall((*link->Memcpy)(link, rootmtype, rootdata, leafmtype, (const char *)leafdata + (size_t)rstart * link->unitbytes, (size_t)sf->nroots * link->unitbytes));
1169566063dSJacob Faibussowitsch     if (PetscMemTypeDevice(leafmtype) && PetscMemTypeHost(rootmtype)) PetscCall((*link->SyncStream)(link));
117dd5b3ca6SJunchao Zhang   } else {
118cd620004SJunchao Zhang     /* Reduce leafdata, then scatter to rootdata */
1199566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
1209566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1219566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata));
1229566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
1239566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_LEAF2ROOT, &rootbuf, &leafbuf, &req, NULL));
1249566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE], &recvcount));
125cd620004SJunchao Zhang     /* Allocate a separate leaf buffer on rank 0 */
126dd400576SPatrick Sanan     if (rank == 0 && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) {
1279566063dSJacob Faibussowitsch       PetscCall(PetscSFMalloc(sf, link->leafmtype_mpi, sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes, (void **)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]));
128dd5b3ca6SJunchao Zhang     }
129cd620004SJunchao 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 */
130dd400576SPatrick Sanan     if (rank == 0 && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE;
1319566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
1329566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link, PETSCSF_LEAF2ROOT));
1339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(leafbuf, link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], count, link->basicunit, op, 0, comm)); /* Must do reduce with MPI builltin datatype basicunit */
1349566063dSJacob Faibussowitsch     PetscCallMPI(MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], dat->recvcounts, dat->displs, unit, rootbuf, recvcount, unit, 0, comm, req));
135dd5b3ca6SJunchao Zhang   }
136eb02082bSJunchao Zhang   PetscFunctionReturn(0);
137eb02082bSJunchao Zhang }
138eb02082bSJunchao Zhang 
1399371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op) {
1409319200aSJunchao Zhang   PetscSFLink link;
1419319200aSJunchao Zhang 
1429319200aSJunchao Zhang   PetscFunctionBegin;
1439319200aSJunchao Zhang   if (op == MPI_REPLACE) {
1449319200aSJunchao Zhang     /* A rare case happens when op is MPI_REPLACE, using GPUs but no GPU aware MPI. In PetscSFReduceBegin_Allgather(v),
1459319200aSJunchao Zhang       we did a device to device copy and in effect finished the communication. But in PetscSFLinkFinishCommunication()
1469319200aSJunchao Zhang       of PetscSFReduceEnd_Basic(), it thinks since there is rootbuf, it calls PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI().
1476aad120cSJose E. Roman       It does a host to device memory copy on rootbuf, wrongly overwriting the results. So we don't overload
1489319200aSJunchao Zhang       PetscSFReduceEnd_Basic() in this case, and just reclaim the link.
1499319200aSJunchao Zhang      */
1509566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
1519566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkReclaim(sf, &link));
1529319200aSJunchao Zhang   } else {
1539566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd_Basic(sf, unit, leafdata, rootdata, op));
1549319200aSJunchao Zhang   }
1559319200aSJunchao Zhang   PetscFunctionReturn(0);
1569319200aSJunchao Zhang }
1579319200aSJunchao Zhang 
1589371c9d4SSatish Balay static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata) {
159cd620004SJunchao Zhang   PetscSFLink link;
160855db38dSJunchao Zhang   PetscMPIInt rank;
161eb02082bSJunchao Zhang 
162eb02082bSJunchao Zhang   PetscFunctionBegin;
1639566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin_Gatherv(sf, unit, rootmtype, rootdata, leafmtype, leafdata, MPI_REPLACE));
1649566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
1659566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF));
1669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
167*48a46eb9SPierre Jolivet   if (rank == 0 && PetscMemTypeDevice(leafmtype) && !sf->use_gpu_aware_mpi) PetscCall((*link->Memcpy)(link, PETSC_MEMTYPE_DEVICE, leafdata, PETSC_MEMTYPE_HOST, link->leafbuf[PETSC_MEMTYPE_HOST], sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes));
1689566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkReclaim(sf, &link));
169dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
170dd5b3ca6SJunchao Zhang }
171dd5b3ca6SJunchao Zhang 
172dd5b3ca6SJunchao 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).
173dd5b3ca6SJunchao Zhang 
174dd5b3ca6SJunchao 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
175dd5b3ca6SJunchao 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
176dd5b3ca6SJunchao 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
177dd5b3ca6SJunchao 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
178dd5b3ca6SJunchao Zhang    0,1,2 is 1,3,6 respectively. root is 10.
179dd5b3ca6SJunchao Zhang 
180dd5b3ca6SJunchao Zhang    We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
181dd5b3ca6SJunchao Zhang              rank-0   rank-1    rank-2
182dd5b3ca6SJunchao Zhang         Root     1
183dd5b3ca6SJunchao Zhang         Leaf     2       3         4
184dd5b3ca6SJunchao Zhang      Leafupdate  2       3         4
185dd5b3ca6SJunchao Zhang 
186dd5b3ca6SJunchao Zhang    Do MPI_Exscan on leafupdate,
187dd5b3ca6SJunchao Zhang              rank-0   rank-1    rank-2
188dd5b3ca6SJunchao Zhang         Root     1
189dd5b3ca6SJunchao Zhang         Leaf     2       3         4
190dd5b3ca6SJunchao Zhang      Leafupdate  2       2         5
191dd5b3ca6SJunchao Zhang 
192dd5b3ca6SJunchao Zhang    BcastAndOp from root to leafupdate,
193dd5b3ca6SJunchao Zhang              rank-0   rank-1    rank-2
194dd5b3ca6SJunchao Zhang         Root     1
195dd5b3ca6SJunchao Zhang         Leaf     2       3         4
196dd5b3ca6SJunchao Zhang      Leafupdate  3       3         6
197dd5b3ca6SJunchao Zhang 
198dd5b3ca6SJunchao Zhang    Copy root to leafupdate on rank-0
199dd5b3ca6SJunchao Zhang              rank-0   rank-1    rank-2
200dd5b3ca6SJunchao Zhang         Root     1
201dd5b3ca6SJunchao Zhang         Leaf     2       3         4
202dd5b3ca6SJunchao Zhang      Leafupdate  1       3         6
203dd5b3ca6SJunchao Zhang 
204dd5b3ca6SJunchao Zhang    Reduce from leaf to root,
205dd5b3ca6SJunchao Zhang              rank-0   rank-1    rank-2
206dd5b3ca6SJunchao Zhang         Root     10
207dd5b3ca6SJunchao Zhang         Leaf     2       3         4
208dd5b3ca6SJunchao Zhang      Leafupdate  1       3         6
209dd5b3ca6SJunchao Zhang */
2109371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, void *leafupdate, MPI_Op op) {
211cd620004SJunchao Zhang   PetscSFLink link;
212dd5b3ca6SJunchao Zhang   MPI_Comm    comm;
213dd5b3ca6SJunchao Zhang   PetscMPIInt count;
214dd5b3ca6SJunchao Zhang 
215dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2169566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
21708401ef6SPierre Jolivet   PetscCheck(!PetscMemTypeDevice(rootmtype) && !PetscMemTypeDevice(leafmtype), comm, PETSC_ERR_SUP, "Do FetchAndOp on device");
218dd5b3ca6SJunchao Zhang   /* Copy leafdata to leafupdate */
2199566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_FETCH, &link));
2209566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata)); /* Sync the device */
2219566063dSJacob Faibussowitsch   PetscCall((*link->Memcpy)(link, leafmtype, leafupdate, leafmtype, leafdata, sf->nleaves * link->unitbytes));
2229566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
223dd5b3ca6SJunchao Zhang 
224dd5b3ca6SJunchao Zhang   /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
22583df288dSJunchao Zhang   if (op == MPI_REPLACE) {
226dd5b3ca6SJunchao Zhang     PetscMPIInt size, rank, prev, next;
2279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
2289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm, &size));
229dd5b3ca6SJunchao Zhang     prev = rank ? rank - 1 : MPI_PROC_NULL;
230dd5b3ca6SJunchao Zhang     next = (rank < size - 1) ? rank + 1 : MPI_PROC_NULL;
2319566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(sf->nleaves, &count));
2329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Sendrecv_replace(leafupdate, count, unit, next, link->tag, prev, link->tag, comm, MPI_STATUSES_IGNORE));
233cd620004SJunchao Zhang   } else {
2349566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
2359566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Exscan(MPI_IN_PLACE, leafupdate, count, link->basicunit, op, comm));
236cd620004SJunchao Zhang   }
2379566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkReclaim(sf, &link));
2389566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, unit, rootdata, leafupdate, op));
2399566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, unit, rootdata, leafupdate, op));
240dd5b3ca6SJunchao Zhang 
241dd5b3ca6SJunchao Zhang   /* Bcast roots to rank 0's leafupdate */
2429566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastToZero_Private(sf, unit, rootdata, leafupdate)); /* Using this line makes Allgather SFs able to inherit this routine */
243dd5b3ca6SJunchao Zhang 
244dd5b3ca6SJunchao Zhang   /* Reduce leafdata to rootdata */
2459566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, unit, leafdata, rootdata, op));
246dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
247dd5b3ca6SJunchao Zhang }
248dd5b3ca6SJunchao Zhang 
2499371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op) {
250dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2519566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, unit, leafdata, rootdata, op));
252dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
253dd5b3ca6SJunchao Zhang }
254dd5b3ca6SJunchao Zhang 
255dd5b3ca6SJunchao Zhang /* Get root ranks accessing my leaves */
2569371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf, PetscInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote) {
257dd5b3ca6SJunchao Zhang   PetscInt        i, j, k, size;
258dd5b3ca6SJunchao Zhang   const PetscInt *range;
259dd5b3ca6SJunchao Zhang 
260dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
261dd5b3ca6SJunchao Zhang   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
262dd5b3ca6SJunchao Zhang   if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
263dd5b3ca6SJunchao Zhang     size = sf->nranks;
2649566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(sf->map, &range));
2659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(size, &sf->ranks, size + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
266dd5b3ca6SJunchao Zhang     for (i = 0; i < size; i++) sf->ranks[i] = i;
2679566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(sf->roffset, range, size + 1));
268dd5b3ca6SJunchao Zhang     for (i = 0; i < sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
269dd5b3ca6SJunchao Zhang     for (i = 0; i < size; i++) {
270dd5b3ca6SJunchao Zhang       for (j = range[i], k = 0; j < range[i + 1]; j++, k++) sf->rremote[j] = k;
271dd5b3ca6SJunchao Zhang     }
272dd5b3ca6SJunchao Zhang   }
273dd5b3ca6SJunchao Zhang 
274dd5b3ca6SJunchao Zhang   if (nranks) *nranks = sf->nranks;
275dd5b3ca6SJunchao Zhang   if (ranks) *ranks = sf->ranks;
276dd5b3ca6SJunchao Zhang   if (roffset) *roffset = sf->roffset;
277dd5b3ca6SJunchao Zhang   if (rmine) *rmine = sf->rmine;
278dd5b3ca6SJunchao Zhang   if (rremote) *rremote = sf->rremote;
279dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
280dd5b3ca6SJunchao Zhang }
281dd5b3ca6SJunchao Zhang 
282dd5b3ca6SJunchao Zhang /* Get leaf ranks accessing my roots */
2839371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf, PetscInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc) {
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 */
2919566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
2929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
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) {
3019566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size, &dat->iranks));
302dd5b3ca6SJunchao Zhang       dat->iranks[0] = rank;
3039371c9d4SSatish Balay       for (i = 0, j = 1; i < size; i++) {
3049371c9d4SSatish Balay         if (i == rank) continue;
3059371c9d4SSatish Balay         dat->iranks[j++] = i;
3069371c9d4SSatish Balay       }
307dd5b3ca6SJunchao Zhang     }
308dd5b3ca6SJunchao Zhang     *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNewLog */
309dd5b3ca6SJunchao Zhang   }
310dd5b3ca6SJunchao Zhang 
311dd5b3ca6SJunchao Zhang   if (ioffset) {
312dd5b3ca6SJunchao Zhang     if (!dat->ioffset) {
3139566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size + 1, &dat->ioffset));
314dd5b3ca6SJunchao Zhang       for (i = 0; i <= size; i++) dat->ioffset[i] = i * sf->nroots;
315dd5b3ca6SJunchao Zhang     }
316dd5b3ca6SJunchao Zhang     *ioffset = dat->ioffset;
317dd5b3ca6SJunchao Zhang   }
318dd5b3ca6SJunchao Zhang 
319dd5b3ca6SJunchao Zhang   if (irootloc) {
320dd5b3ca6SJunchao Zhang     if (!dat->irootloc) {
3219566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sf->nleaves, &dat->irootloc));
322dd5b3ca6SJunchao Zhang       for (i = 0; i < size; i++) {
323dd5b3ca6SJunchao Zhang         for (j = 0; j < sf->nroots; j++) dat->irootloc[i * sf->nroots + j] = j;
324dd5b3ca6SJunchao Zhang       }
325dd5b3ca6SJunchao Zhang     }
326dd5b3ca6SJunchao Zhang     *irootloc = dat->irootloc;
327dd5b3ca6SJunchao Zhang   }
328dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
329dd5b3ca6SJunchao Zhang }
330dd5b3ca6SJunchao Zhang 
3319371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf, PetscSF *out) {
332dd5b3ca6SJunchao Zhang   PetscInt     i, nroots, nleaves, rstart, *ilocal;
333dd5b3ca6SJunchao Zhang   PetscSFNode *iremote;
334dd5b3ca6SJunchao Zhang   PetscSF      lsf;
335dd5b3ca6SJunchao Zhang 
336dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
337eb02082bSJunchao Zhang   nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */
338eb02082bSJunchao Zhang   nroots  = nleaves;
3399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &ilocal));
3409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &iremote));
3419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
342dd5b3ca6SJunchao Zhang 
343dd5b3ca6SJunchao Zhang   for (i = 0; i < nleaves; i++) {
344dd5b3ca6SJunchao Zhang     ilocal[i]        = rstart + i; /* lsf does not change leave indices */
345dd5b3ca6SJunchao Zhang     iremote[i].rank  = 0;          /* rank in PETSC_COMM_SELF */
346dd5b3ca6SJunchao Zhang     iremote[i].index = i;          /* root index */
347dd5b3ca6SJunchao Zhang   }
348dd5b3ca6SJunchao Zhang 
3499566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
3509566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(lsf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
3519566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(lsf));
352dd5b3ca6SJunchao Zhang   *out = lsf;
353dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
354dd5b3ca6SJunchao Zhang }
355dd5b3ca6SJunchao Zhang 
3569371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf) {
357dd5b3ca6SJunchao Zhang   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
358dd5b3ca6SJunchao Zhang 
359dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
360ad227feaSJunchao Zhang   sf->ops->BcastEnd  = PetscSFBcastEnd_Basic;
3619319200aSJunchao Zhang   sf->ops->ReduceEnd = PetscSFReduceEnd_Allgatherv;
362cd620004SJunchao Zhang 
363dd5b3ca6SJunchao Zhang   sf->ops->SetUp           = PetscSFSetUp_Allgatherv;
364dd5b3ca6SJunchao Zhang   sf->ops->Reset           = PetscSFReset_Allgatherv;
365dd5b3ca6SJunchao Zhang   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
366dd5b3ca6SJunchao Zhang   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
367dd5b3ca6SJunchao Zhang   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
368dd5b3ca6SJunchao Zhang   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
369ad227feaSJunchao Zhang   sf->ops->BcastBegin      = PetscSFBcastBegin_Allgatherv;
370dd5b3ca6SJunchao Zhang   sf->ops->ReduceBegin     = PetscSFReduceBegin_Allgatherv;
371dd5b3ca6SJunchao Zhang   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
372dd5b3ca6SJunchao Zhang   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
373dd5b3ca6SJunchao Zhang   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
374dd5b3ca6SJunchao Zhang   sf->ops->BcastToZero     = PetscSFBcastToZero_Allgatherv;
375dd5b3ca6SJunchao Zhang 
3769566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(sf, &dat));
377dd5b3ca6SJunchao Zhang   sf->data = (void *)dat;
378dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
379dd5b3ca6SJunchao Zhang }
380