xref: /petsc/src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.c (revision 4ad8454beace47809662cdae21ee081016eaa39a)
1 #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h>
2 
3 /* PetscSFGetGraph is non-collective. An implementation should not have collective calls */
4 PETSC_INTERN PetscErrorCode PetscSFGetGraph_Allgatherv(PetscSF sf, PetscInt *nroots, PetscInt *nleaves, const PetscInt **ilocal, const PetscSFNode **iremote)
5 {
6   PetscInt        i, j, k;
7   const PetscInt *range;
8   PetscMPIInt     size;
9 
10   PetscFunctionBegin;
11   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
12   if (nroots) *nroots = sf->nroots;
13   if (nleaves) *nleaves = sf->nleaves;
14   if (ilocal) *ilocal = NULL; /* Contiguous leaves */
15   if (iremote) {
16     if (!sf->remote && sf->nleaves) { /* The && sf->nleaves makes sfgatherv able to inherit this routine */
17       PetscCall(PetscLayoutGetRanges(sf->map, &range));
18       PetscCall(PetscMalloc1(sf->nleaves, &sf->remote));
19       sf->remote_alloc = sf->remote;
20       for (i = 0; i < size; i++) {
21         for (j = range[i], k = 0; j < range[i + 1]; j++, k++) {
22           sf->remote[j].rank  = i;
23           sf->remote[j].index = k;
24         }
25       }
26     }
27     *iremote = sf->remote;
28   }
29   PetscFunctionReturn(PETSC_SUCCESS);
30 }
31 
32 PETSC_INTERN PetscErrorCode PetscSFSetUp_Allgatherv(PetscSF sf)
33 {
34   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
35   PetscMPIInt         size;
36   PetscInt            i;
37   const PetscInt     *range;
38   MPI_Comm            comm;
39 
40   PetscFunctionBegin;
41   PetscCall(PetscSFSetUp_Allgather(sf));
42   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
43   PetscCallMPI(MPI_Comm_size(comm, &size));
44   if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */
45     PetscBool isallgatherv = PETSC_FALSE;
46 
47     PetscCall(PetscMalloc1(size, &dat->recvcounts));
48     PetscCall(PetscMalloc1(size, &dat->displs));
49     PetscCall(PetscLayoutGetRanges(sf->map, &range));
50 
51     for (i = 0; i < size; i++) {
52       PetscCall(PetscMPIIntCast(range[i], &dat->displs[i]));
53       PetscCall(PetscMPIIntCast(range[i + 1] - range[i], &dat->recvcounts[i]));
54     }
55 
56     /* check if we actually have a one-to-all pattern */
57     PetscCall(PetscObjectTypeCompare((PetscObject)sf, PETSCSFALLGATHERV, &isallgatherv));
58     if (isallgatherv) {
59       PetscMPIInt rank, nRanksWithZeroRoots;
60 
61       nRanksWithZeroRoots = (sf->nroots == 0) ? 1 : 0; /* I have no roots */
62       PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &nRanksWithZeroRoots, 1, MPI_INT, MPI_SUM, comm));
63       if (nRanksWithZeroRoots == size - 1) { /* Only one rank has roots, which indicates a bcast pattern */
64         dat->bcast_pattern = PETSC_TRUE;
65         PetscCallMPI(MPI_Comm_rank(comm, &rank));
66         dat->bcast_root = sf->nroots > 0 ? rank : -1;
67         PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &dat->bcast_root, 1, MPI_INT, MPI_MAX, comm));
68       }
69     }
70   }
71   PetscFunctionReturn(PETSC_SUCCESS);
72 }
73 
74 PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf)
75 {
76   PetscSF_Allgatherv *dat  = (PetscSF_Allgatherv *)sf->data;
77   PetscSFLink         link = dat->avail, next;
78 
79   PetscFunctionBegin;
80   PetscCall(PetscFree(dat->iranks));
81   PetscCall(PetscFree(dat->ioffset));
82   PetscCall(PetscFree(dat->irootloc));
83   PetscCall(PetscFree(dat->recvcounts));
84   PetscCall(PetscFree(dat->displs));
85   PetscCheck(!dat->inuse, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Outstanding operation has not been completed");
86   for (; link; link = next) {
87     next = link->next;
88     PetscCall(PetscSFLinkDestroy(sf, link));
89   }
90   dat->avail = NULL;
91   PetscFunctionReturn(PETSC_SUCCESS);
92 }
93 
94 PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf)
95 {
96   PetscFunctionBegin;
97   PetscCall(PetscSFReset_Allgatherv(sf));
98   PetscCall(PetscFree(sf->data));
99   PetscFunctionReturn(PETSC_SUCCESS);
100 }
101 
102 static PetscErrorCode PetscSFBcastBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op)
103 {
104   PetscSFLink         link;
105   PetscMPIInt         sendcount, rank;
106   MPI_Comm            comm;
107   void               *rootbuf = NULL, *leafbuf = NULL;
108   MPI_Request        *req = NULL;
109   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
110 
111   PetscFunctionBegin;
112   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_BCAST, &link));
113   PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
114   PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
115   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
116   PetscCallMPI(MPI_Comm_rank(comm, &rank));
117   PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
118   PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
119 
120   if (dat->bcast_pattern && rank == dat->bcast_root) PetscCall((*link->Memcpy)(link, link->leafmtype_mpi, leafbuf, link->rootmtype_mpi, rootbuf, (size_t)sendcount * link->unitbytes));
121   /* Ready the buffers for MPI */
122   PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link, PETSCSF_ROOT2LEAF));
123   if (dat->bcast_pattern) PetscCallMPI(MPIU_Ibcast(leafbuf, sf->nleaves, unit, dat->bcast_root, comm, req));
124   else PetscCallMPI(MPIU_Iallgatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, comm, req));
125   PetscFunctionReturn(PETSC_SUCCESS);
126 }
127 
128 static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op)
129 {
130   PetscSFLink         link;
131   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
132   PetscInt            rstart;
133   PetscMPIInt         rank, count, recvcount;
134   MPI_Comm            comm;
135   void               *rootbuf = NULL, *leafbuf = NULL;
136   MPI_Request        *req = NULL;
137 
138   PetscFunctionBegin;
139   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_REDUCE, &link));
140   if (op == MPI_REPLACE) {
141     /* REPLACE is only meaningful when all processes have the same leafdata to reduce. Therefore copying from local leafdata is fine */
142     PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
143     PetscCall((*link->Memcpy)(link, rootmtype, rootdata, leafmtype, (const char *)leafdata + (size_t)rstart * link->unitbytes, (size_t)sf->nroots * link->unitbytes));
144     if (PetscMemTypeDevice(leafmtype) && PetscMemTypeHost(rootmtype)) PetscCall((*link->SyncStream)(link));
145   } else {
146     PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
147     PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata));
148     PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
149     PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_LEAF2ROOT, &rootbuf, &leafbuf, &req, NULL));
150     PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link, PETSCSF_LEAF2ROOT));
151     if (dat->bcast_pattern) {
152 #if defined(PETSC_HAVE_OMPI_MAJOR_VERSION) /* Workaround: cuda-aware Open MPI 4.1.3 does not support MPI_Ireduce() with device buffers */
153       *req = MPI_REQUEST_NULL;             /* Set NULL so that we can safely MPI_Wait(req) */
154       PetscCallMPI(MPI_Reduce(leafbuf, rootbuf, sf->nleaves, unit, op, dat->bcast_root, comm));
155 #else
156       PetscCallMPI(MPIU_Ireduce(leafbuf, rootbuf, sf->nleaves, unit, op, dat->bcast_root, comm, req));
157 #endif
158     } else { /* Reduce leafdata, then scatter to rootdata */
159       PetscCallMPI(MPI_Comm_rank(comm, &rank));
160       PetscCall(PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE], &recvcount));
161       /* Allocate a separate leaf buffer on rank 0 */
162       if (rank == 0 && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) {
163         PetscCall(PetscSFMalloc(sf, link->leafmtype_mpi, sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes, (void **)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]));
164       }
165       /* 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 */
166       if (rank == 0 && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE;
167       PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
168       PetscCallMPI(MPI_Reduce(leafbuf, link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], count, link->basicunit, op, 0, comm)); /* Must do reduce with MPI builtin datatype basicunit */
169       PetscCallMPI(MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], dat->recvcounts, dat->displs, unit, rootbuf, recvcount, unit, 0, comm, req));
170     }
171   }
172   PetscFunctionReturn(PETSC_SUCCESS);
173 }
174 
175 PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
176 {
177   PetscSFLink link;
178 
179   PetscFunctionBegin;
180   if (op == MPI_REPLACE) {
181     /* A rare case happens when op is MPI_REPLACE, using GPUs but no GPU aware MPI. In PetscSFReduceBegin_Allgather(v),
182       we did a device to device copy and in effect finished the communication. But in PetscSFLinkFinishCommunication()
183       of PetscSFReduceEnd_Basic(), it thinks since there is rootbuf, it calls PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI().
184       It does a host to device memory copy on rootbuf, wrongly overwriting the results. So we don't overload
185       PetscSFReduceEnd_Basic() in this case, and just reclaim the link.
186      */
187     PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
188     PetscCall(PetscSFLinkReclaim(sf, &link));
189   } else {
190     PetscCall(PetscSFReduceEnd_Basic(sf, unit, leafdata, rootdata, op));
191   }
192   PetscFunctionReturn(PETSC_SUCCESS);
193 }
194 
195 static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata)
196 {
197   PetscSFLink         link;
198   PetscMPIInt         rank;
199   PetscMPIInt         sendcount;
200   MPI_Comm            comm;
201   PetscSF_Allgatherv *dat     = (PetscSF_Allgatherv *)sf->data;
202   void               *rootbuf = NULL, *leafbuf = NULL; /* buffer seen by MPI */
203   MPI_Request        *req = NULL;
204 
205   PetscFunctionBegin;
206   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, MPI_REPLACE, PETSCSF_BCAST, &link));
207   PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
208   PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
209   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
210   PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
211   PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
212   PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link, PETSCSF_ROOT2LEAF));
213   PetscCallMPI(MPIU_Igatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, 0 /*rank 0*/, comm, req));
214 
215   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
216   PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF));
217   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
218   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));
219   PetscCall(PetscSFLinkReclaim(sf, &link));
220   PetscFunctionReturn(PETSC_SUCCESS);
221 }
222 
223 /* 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).
224 
225    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
226    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
227    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
228    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
229    0,1,2 is 1,3,6 respectively. root is 10.
230 
231    We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
232              rank-0   rank-1    rank-2
233         Root     1
234         Leaf     2       3         4
235      Leafupdate  2       3         4
236 
237    Do MPI_Exscan on leafupdate,
238              rank-0   rank-1    rank-2
239         Root     1
240         Leaf     2       3         4
241      Leafupdate  2       2         5
242 
243    BcastAndOp from root to leafupdate,
244              rank-0   rank-1    rank-2
245         Root     1
246         Leaf     2       3         4
247      Leafupdate  3       3         6
248 
249    Copy root to leafupdate on rank-0
250              rank-0   rank-1    rank-2
251         Root     1
252         Leaf     2       3         4
253      Leafupdate  1       3         6
254 
255    Reduce from leaf to root,
256              rank-0   rank-1    rank-2
257         Root     10
258         Leaf     2       3         4
259      Leafupdate  1       3         6
260 */
261 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, void *leafupdate, MPI_Op op)
262 {
263   PetscSFLink link;
264   MPI_Comm    comm;
265   PetscMPIInt count;
266 
267   PetscFunctionBegin;
268   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
269   PetscCheck(!PetscMemTypeDevice(rootmtype) && !PetscMemTypeDevice(leafmtype), comm, PETSC_ERR_SUP, "Do FetchAndOp on device");
270   /* Copy leafdata to leafupdate */
271   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_FETCH, &link));
272   PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata)); /* Sync the device */
273   PetscCall((*link->Memcpy)(link, leafmtype, leafupdate, leafmtype, leafdata, sf->nleaves * link->unitbytes));
274   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
275 
276   /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
277   if (op == MPI_REPLACE) {
278     PetscMPIInt size, rank, prev, next;
279     PetscCallMPI(MPI_Comm_rank(comm, &rank));
280     PetscCallMPI(MPI_Comm_size(comm, &size));
281     prev = rank ? rank - 1 : MPI_PROC_NULL;
282     next = (rank < size - 1) ? rank + 1 : MPI_PROC_NULL;
283     PetscCall(PetscMPIIntCast(sf->nleaves, &count));
284     PetscCallMPI(MPI_Sendrecv_replace(leafupdate, count, unit, next, link->tag, prev, link->tag, comm, MPI_STATUSES_IGNORE));
285   } else {
286     PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
287     PetscCallMPI(MPI_Exscan(MPI_IN_PLACE, leafupdate, count, link->basicunit, op, comm));
288   }
289   PetscCall(PetscSFLinkReclaim(sf, &link));
290   PetscCall(PetscSFBcastBegin(sf, unit, rootdata, leafupdate, op));
291   PetscCall(PetscSFBcastEnd(sf, unit, rootdata, leafupdate, op));
292 
293   /* Bcast roots to rank 0's leafupdate */
294   PetscCall(PetscSFBcastToZero_Private(sf, unit, rootdata, leafupdate)); /* Using this line makes Allgather SFs able to inherit this routine */
295 
296   /* Reduce leafdata to rootdata */
297   PetscCall(PetscSFReduceBegin(sf, unit, leafdata, rootdata, op));
298   PetscFunctionReturn(PETSC_SUCCESS);
299 }
300 
301 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
302 {
303   PetscFunctionBegin;
304   PetscCall(PetscSFReduceEnd(sf, unit, leafdata, rootdata, op));
305   PetscFunctionReturn(PETSC_SUCCESS);
306 }
307 
308 /* Get root ranks accessing my leaves */
309 PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf, PetscInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote)
310 {
311   PetscInt        i, j, k, size;
312   const PetscInt *range;
313 
314   PetscFunctionBegin;
315   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
316   if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
317     size = sf->nranks;
318     PetscCall(PetscLayoutGetRanges(sf->map, &range));
319     PetscCall(PetscMalloc4(size, &sf->ranks, size + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
320     for (i = 0; i < size; i++) sf->ranks[i] = i;
321     PetscCall(PetscArraycpy(sf->roffset, range, size + 1));
322     for (i = 0; i < sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
323     for (i = 0; i < size; i++) {
324       for (j = range[i], k = 0; j < range[i + 1]; j++, k++) sf->rremote[j] = k;
325     }
326   }
327 
328   if (nranks) *nranks = sf->nranks;
329   if (ranks) *ranks = sf->ranks;
330   if (roffset) *roffset = sf->roffset;
331   if (rmine) *rmine = sf->rmine;
332   if (rremote) *rremote = sf->rremote;
333   PetscFunctionReturn(PETSC_SUCCESS);
334 }
335 
336 /* Get leaf ranks accessing my roots */
337 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf, PetscInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc)
338 {
339   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
340   MPI_Comm            comm;
341   PetscMPIInt         size, rank;
342   PetscInt            i, j;
343 
344   PetscFunctionBegin;
345   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
346   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
347   PetscCallMPI(MPI_Comm_size(comm, &size));
348   PetscCallMPI(MPI_Comm_rank(comm, &rank));
349   if (niranks) *niranks = size;
350 
351   /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and
352      sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why.
353    */
354   if (iranks) {
355     if (!dat->iranks) {
356       PetscCall(PetscMalloc1(size, &dat->iranks));
357       dat->iranks[0] = rank;
358       for (i = 0, j = 1; i < size; i++) {
359         if (i == rank) continue;
360         dat->iranks[j++] = i;
361       }
362     }
363     *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNew */
364   }
365 
366   if (ioffset) {
367     if (!dat->ioffset) {
368       PetscCall(PetscMalloc1(size + 1, &dat->ioffset));
369       for (i = 0; i <= size; i++) dat->ioffset[i] = i * sf->nroots;
370     }
371     *ioffset = dat->ioffset;
372   }
373 
374   if (irootloc) {
375     if (!dat->irootloc) {
376       PetscCall(PetscMalloc1(sf->nleaves, &dat->irootloc));
377       for (i = 0; i < size; i++) {
378         for (j = 0; j < sf->nroots; j++) dat->irootloc[i * sf->nroots + j] = j;
379       }
380     }
381     *irootloc = dat->irootloc;
382   }
383   PetscFunctionReturn(PETSC_SUCCESS);
384 }
385 
386 PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf, PetscSF *out)
387 {
388   PetscInt     i, nroots, nleaves, rstart, *ilocal;
389   PetscSFNode *iremote;
390   PetscSF      lsf;
391 
392   PetscFunctionBegin;
393   nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */
394   nroots  = nleaves;
395   PetscCall(PetscMalloc1(nleaves, &ilocal));
396   PetscCall(PetscMalloc1(nleaves, &iremote));
397   PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
398 
399   for (i = 0; i < nleaves; i++) {
400     ilocal[i]        = rstart + i; /* lsf does not change leave indices */
401     iremote[i].rank  = 0;          /* rank in PETSC_COMM_SELF */
402     iremote[i].index = i;          /* root index */
403   }
404 
405   PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
406   PetscCall(PetscSFSetGraph(lsf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
407   PetscCall(PetscSFSetUp(lsf));
408   *out = lsf;
409   PetscFunctionReturn(PETSC_SUCCESS);
410 }
411 
412 PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf)
413 {
414   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
415 
416   PetscFunctionBegin;
417   sf->ops->BcastEnd  = PetscSFBcastEnd_Basic;
418   sf->ops->ReduceEnd = PetscSFReduceEnd_Allgatherv;
419 
420   sf->ops->SetUp           = PetscSFSetUp_Allgatherv;
421   sf->ops->Reset           = PetscSFReset_Allgatherv;
422   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
423   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
424   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
425   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
426   sf->ops->BcastBegin      = PetscSFBcastBegin_Allgatherv;
427   sf->ops->ReduceBegin     = PetscSFReduceBegin_Allgatherv;
428   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
429   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
430   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
431   sf->ops->BcastToZero     = PetscSFBcastToZero_Allgatherv;
432 
433   sf->collective = PETSC_TRUE;
434 
435   PetscCall(PetscNew(&dat));
436   dat->bcast_root = -1;
437   sf->data        = (void *)dat;
438   PetscFunctionReturn(PETSC_SUCCESS);
439 }
440