xref: /petsc/src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.c (revision e91c04dfc8a52dee1965211bb1cc8e5bf775178f)
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        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 (PetscMPIInt 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(MPIU_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(MPIU_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));
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));
151     if (dat->bcast_pattern) {
152       PetscMPIInt nleavesi;
153 
154       PetscCall(PetscMPIIntCast(sf->nleaves, &nleavesi));
155 #if defined(PETSC_HAVE_OPENMPI) /* Workaround: cuda-aware Open MPI 4.1.3 does not support MPI_Ireduce() with device buffers */
156       *req = MPI_REQUEST_NULL;  /* Set NULL so that we can safely MPI_Wait(req) */
157       PetscCallMPI(MPIU_Reduce(leafbuf, rootbuf, nleavesi, unit, op, dat->bcast_root, comm));
158 #else
159       PetscCallMPI(MPIU_Ireduce(leafbuf, rootbuf, nleavesi, unit, op, dat->bcast_root, comm, req));
160 #endif
161     } else { /* Reduce leafdata, then scatter to rootdata */
162       PetscCallMPI(MPI_Comm_rank(comm, &rank));
163       PetscCall(PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE], &recvcount));
164       /* Allocate a separate leaf buffer on rank 0 */
165       if (rank == 0 && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) {
166         PetscCall(PetscSFMalloc(sf, link->leafmtype_mpi, sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes, (void **)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]));
167       }
168       /* 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 */
169       if (rank == 0 && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE;
170       PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
171       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 */
172       PetscCallMPI(MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], dat->recvcounts, dat->displs, unit, rootbuf, recvcount, unit, 0, comm, req));
173     }
174   }
175   PetscFunctionReturn(PETSC_SUCCESS);
176 }
177 
178 PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
179 {
180   PetscSFLink link;
181 
182   PetscFunctionBegin;
183   if (op == MPI_REPLACE) {
184     /* A rare case happens when op is MPI_REPLACE, using GPUs but no GPU aware MPI. In PetscSFReduceBegin_Allgather(v),
185       we did a device to device copy and in effect finished the communication. But in PetscSFLinkFinishCommunication()
186       of PetscSFReduceEnd_Basic(), it thinks since there is rootbuf, it calls PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI().
187       It does a host to device memory copy on rootbuf, wrongly overwriting the results. So we don't overload
188       PetscSFReduceEnd_Basic() in this case, and just reclaim the link.
189      */
190     PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
191     PetscCall(PetscSFLinkReclaim(sf, &link));
192   } else {
193     PetscCall(PetscSFReduceEnd_Basic(sf, unit, leafdata, rootdata, op));
194   }
195   PetscFunctionReturn(PETSC_SUCCESS);
196 }
197 
198 static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata)
199 {
200   PetscSFLink         link;
201   PetscMPIInt         rank;
202   PetscMPIInt         sendcount;
203   MPI_Comm            comm;
204   PetscSF_Allgatherv *dat     = (PetscSF_Allgatherv *)sf->data;
205   void               *rootbuf = NULL, *leafbuf = NULL; /* buffer seen by MPI */
206   MPI_Request        *req = NULL;
207 
208   PetscFunctionBegin;
209   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, MPI_REPLACE, PETSCSF_BCAST, &link));
210   PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
211   PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
212   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
213   PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
214   PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
215   PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
216   PetscCallMPI(MPIU_Igatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, 0 /*rank 0*/, comm, req));
217 
218   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
219   PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF));
220   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
221   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));
222   PetscCall(PetscSFLinkReclaim(sf, &link));
223   PetscFunctionReturn(PETSC_SUCCESS);
224 }
225 
226 /* 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).
227 
228    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
229    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
230    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
231    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
232    0,1,2 is 1,3,6 respectively. root is 10.
233 
234    We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
235              rank-0   rank-1    rank-2
236         Root     1
237         Leaf     2       3         4
238      Leafupdate  2       3         4
239 
240    Do MPI_Exscan on leafupdate,
241              rank-0   rank-1    rank-2
242         Root     1
243         Leaf     2       3         4
244      Leafupdate  2       2         5
245 
246    BcastAndOp from root to leafupdate,
247              rank-0   rank-1    rank-2
248         Root     1
249         Leaf     2       3         4
250      Leafupdate  3       3         6
251 
252    Copy root to leafupdate on rank-0
253              rank-0   rank-1    rank-2
254         Root     1
255         Leaf     2       3         4
256      Leafupdate  1       3         6
257 
258    Reduce from leaf to root,
259              rank-0   rank-1    rank-2
260         Root     10
261         Leaf     2       3         4
262      Leafupdate  1       3         6
263 */
264 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, void *leafupdate, MPI_Op op)
265 {
266   PetscSFLink link;
267   MPI_Comm    comm;
268   PetscMPIInt count;
269 
270   PetscFunctionBegin;
271   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
272   PetscCheck(!PetscMemTypeDevice(rootmtype) && !PetscMemTypeDevice(leafmtype), comm, PETSC_ERR_SUP, "Do FetchAndOp on device");
273   /* Copy leafdata to leafupdate */
274   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_FETCH, &link));
275   PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata)); /* Sync the device */
276   PetscCall((*link->Memcpy)(link, leafmtype, leafupdate, leafmtype, leafdata, sf->nleaves * link->unitbytes));
277   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
278 
279   /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
280   if (op == MPI_REPLACE) {
281     PetscMPIInt size, rank, prev, next;
282     PetscCallMPI(MPI_Comm_rank(comm, &rank));
283     PetscCallMPI(MPI_Comm_size(comm, &size));
284     prev = rank ? rank - 1 : MPI_PROC_NULL;
285     next = (rank < size - 1) ? rank + 1 : MPI_PROC_NULL;
286     PetscCall(PetscMPIIntCast(sf->nleaves, &count));
287     PetscCallMPI(MPI_Sendrecv_replace(leafupdate, count, unit, next, link->tag, prev, link->tag, comm, MPI_STATUSES_IGNORE));
288   } else {
289     PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
290     PetscCallMPI(MPI_Exscan(MPI_IN_PLACE, leafupdate, count, link->basicunit, op, comm));
291   }
292   PetscCall(PetscSFLinkReclaim(sf, &link));
293   PetscCall(PetscSFBcastBegin(sf, unit, rootdata, leafupdate, op));
294   PetscCall(PetscSFBcastEnd(sf, unit, rootdata, leafupdate, op));
295 
296   /* Bcast roots to rank 0's leafupdate */
297   PetscCall(PetscSFBcastToZero_Private(sf, unit, rootdata, leafupdate)); /* Using this line makes Allgather SFs able to inherit this routine */
298 
299   /* Reduce leafdata to rootdata */
300   PetscCall(PetscSFReduceBegin(sf, unit, leafdata, rootdata, op));
301   PetscFunctionReturn(PETSC_SUCCESS);
302 }
303 
304 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
305 {
306   PetscFunctionBegin;
307   PetscCall(PetscSFReduceEnd(sf, unit, leafdata, rootdata, op));
308   PetscFunctionReturn(PETSC_SUCCESS);
309 }
310 
311 /* Get root ranks accessing my leaves */
312 PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf, PetscMPIInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote)
313 {
314   PetscInt        j, k, size;
315   const PetscInt *range;
316 
317   PetscFunctionBegin;
318   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
319   if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
320     size = sf->nranks;
321     PetscCall(PetscLayoutGetRanges(sf->map, &range));
322     PetscCall(PetscMalloc4(size, &sf->ranks, size + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
323     for (PetscMPIInt i = 0; i < size; i++) sf->ranks[i] = i;
324     PetscCall(PetscArraycpy(sf->roffset, range, size + 1));
325     for (PetscInt i = 0; i < sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
326     for (PetscMPIInt i = 0; i < size; i++) {
327       for (j = range[i], k = 0; j < range[i + 1]; j++, k++) sf->rremote[j] = k;
328     }
329   }
330 
331   if (nranks) *nranks = sf->nranks;
332   if (ranks) *ranks = sf->ranks;
333   if (roffset) *roffset = sf->roffset;
334   if (rmine) *rmine = sf->rmine;
335   if (rremote) *rremote = sf->rremote;
336   PetscFunctionReturn(PETSC_SUCCESS);
337 }
338 
339 /* Get leaf ranks accessing my roots */
340 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf, PetscMPIInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc)
341 {
342   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
343   MPI_Comm            comm;
344   PetscMPIInt         size, rank;
345 
346   PetscFunctionBegin;
347   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
348   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
349   PetscCallMPI(MPI_Comm_size(comm, &size));
350   PetscCallMPI(MPI_Comm_rank(comm, &rank));
351   if (niranks) *niranks = size;
352 
353   /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and
354      sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why.
355    */
356   if (iranks) {
357     if (!dat->iranks) {
358       PetscCall(PetscMalloc1(size, &dat->iranks));
359       dat->iranks[0] = rank;
360       for (PetscMPIInt i = 0, j = 1; i < size; i++) {
361         if (i == rank) continue;
362         dat->iranks[j++] = i;
363       }
364     }
365     *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNew */
366   }
367 
368   if (ioffset) {
369     if (!dat->ioffset) {
370       PetscCall(PetscMalloc1(size + 1, &dat->ioffset));
371       for (PetscMPIInt i = 0; i <= size; i++) dat->ioffset[i] = i * sf->nroots;
372     }
373     *ioffset = dat->ioffset;
374   }
375 
376   if (irootloc) {
377     if (!dat->irootloc) {
378       PetscCall(PetscMalloc1(sf->nleaves, &dat->irootloc));
379       for (PetscMPIInt i = 0; i < size; i++) {
380         for (PetscInt j = 0; j < sf->nroots; j++) dat->irootloc[i * sf->nroots + j] = j;
381       }
382     }
383     *irootloc = dat->irootloc;
384   }
385   PetscFunctionReturn(PETSC_SUCCESS);
386 }
387 
388 PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf, PetscSF *out)
389 {
390   PetscInt     i, nroots, nleaves, rstart, *ilocal;
391   PetscSFNode *iremote;
392   PetscSF      lsf;
393 
394   PetscFunctionBegin;
395   nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */
396   nroots  = nleaves;
397   PetscCall(PetscMalloc1(nleaves, &ilocal));
398   PetscCall(PetscMalloc1(nleaves, &iremote));
399   PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
400 
401   for (i = 0; i < nleaves; i++) {
402     ilocal[i]        = rstart + i; /* lsf does not change leave indices */
403     iremote[i].rank  = 0;          /* rank in PETSC_COMM_SELF */
404     iremote[i].index = i;          /* root index */
405   }
406 
407   PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
408   PetscCall(PetscSFSetGraph(lsf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
409   PetscCall(PetscSFSetUp(lsf));
410   *out = lsf;
411   PetscFunctionReturn(PETSC_SUCCESS);
412 }
413 
414 PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf)
415 {
416   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
417 
418   PetscFunctionBegin;
419   sf->ops->BcastEnd  = PetscSFBcastEnd_Basic;
420   sf->ops->ReduceEnd = PetscSFReduceEnd_Allgatherv;
421 
422   sf->ops->SetUp           = PetscSFSetUp_Allgatherv;
423   sf->ops->Reset           = PetscSFReset_Allgatherv;
424   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
425   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
426   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
427   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
428   sf->ops->BcastBegin      = PetscSFBcastBegin_Allgatherv;
429   sf->ops->ReduceBegin     = PetscSFReduceBegin_Allgatherv;
430   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
431   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
432   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
433   sf->ops->BcastToZero     = PetscSFBcastToZero_Allgatherv;
434 
435   sf->collective = PETSC_TRUE;
436 
437   PetscCall(PetscNew(&dat));
438   dat->bcast_root = -1;
439   sf->data        = (void *)dat;
440   PetscFunctionReturn(PETSC_SUCCESS);
441 }
442