xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision 750b007cd8d816cecd9de99077bb0a703b4cf61a) !
1 #include <petsc/private/sfimpl.h> /*I  "petscsf.h"   I*/
2 #include <petsc/private/sectionimpl.h>
3 
4 /*@
5    PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout
6 
7    Collective
8 
9    Input Parameters:
10 +  sf - star forest
11 .  layout - PetscLayout defining the global space for roots
12 .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
13 .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
14 .  localmode - copy mode for ilocal
15 -  gremote - root vertices in global numbering corresponding to leaves in ilocal
16 
17    Level: intermediate
18 
19    Notes:
20    Global indices must lie in [0, N) where N is the global size of layout.
21    Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is PETSC_OWN_POINTER.
22 
23    Developer Notes:
24    Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
25    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
26    needed
27 
28 .seealso: `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
29 @*/
30 PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *gremote) {
31   const PetscInt *range;
32   PetscInt        i, nroots, ls = -1, ln = -1;
33   PetscMPIInt     lr = -1;
34   PetscSFNode    *remote;
35 
36   PetscFunctionBegin;
37   PetscCall(PetscLayoutGetLocalSize(layout, &nroots));
38   PetscCall(PetscLayoutGetRanges(layout, &range));
39   PetscCall(PetscMalloc1(nleaves, &remote));
40   if (nleaves) ls = gremote[0] + 1;
41   for (i = 0; i < nleaves; i++) {
42     const PetscInt idx = gremote[i] - ls;
43     if (idx < 0 || idx >= ln) { /* short-circuit the search */
44       PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index));
45       remote[i].rank = lr;
46       ls             = range[lr];
47       ln             = range[lr + 1] - ls;
48     } else {
49       remote[i].rank  = lr;
50       remote[i].index = idx;
51     }
52   }
53   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER));
54   PetscFunctionReturn(0);
55 }
56 
57 /*@C
58    PetscSFGetGraphLayout - Get the global indices and PetscLayout that describe this star forest
59 
60    Collective
61 
62    Input Parameter:
63 .  sf - star forest
64 
65    Output Parameters:
66 +  layout - PetscLayout defining the global space for roots
67 .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
68 .  ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage
69 -  gremote - root vertices in global numbering corresponding to leaves in ilocal
70 
71    Level: intermediate
72 
73    Notes:
74    The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest.
75    The outputs layout and gremote are freshly created each time this function is called,
76    so they need to be freed by user and cannot be qualified as const.
77 
78 
79 .seealso: `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
80 @*/
81 PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[]) {
82   PetscInt           nr, nl;
83   const PetscSFNode *ir;
84   PetscLayout        lt;
85 
86   PetscFunctionBegin;
87   PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir));
88   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, &lt));
89   if (gremote) {
90     PetscInt        i;
91     const PetscInt *range;
92     PetscInt       *gr;
93 
94     PetscCall(PetscLayoutGetRanges(lt, &range));
95     PetscCall(PetscMalloc1(nl, &gr));
96     for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index;
97     *gremote = gr;
98   }
99   if (nleaves) *nleaves = nl;
100   if (layout) *layout = lt;
101   else PetscCall(PetscLayoutDestroy(&lt));
102   PetscFunctionReturn(0);
103 }
104 
105 /*@
106   PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections
107   describing the data layout.
108 
109   Input Parameters:
110 + sf - The SF
111 . localSection - PetscSection describing the local data layout
112 - globalSection - PetscSection describing the global data layout
113 
114   Level: developer
115 
116 .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphLayout()`
117 @*/
118 PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection) {
119   MPI_Comm        comm;
120   PetscLayout     layout;
121   const PetscInt *ranges;
122   PetscInt       *local;
123   PetscSFNode    *remote;
124   PetscInt        pStart, pEnd, p, nroots, nleaves = 0, l;
125   PetscMPIInt     size, rank;
126 
127   PetscFunctionBegin;
128   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
129   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
130   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
131 
132   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
133   PetscCallMPI(MPI_Comm_size(comm, &size));
134   PetscCallMPI(MPI_Comm_rank(comm, &rank));
135   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
136   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
137   PetscCall(PetscLayoutCreate(comm, &layout));
138   PetscCall(PetscLayoutSetBlockSize(layout, 1));
139   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
140   PetscCall(PetscLayoutSetUp(layout));
141   PetscCall(PetscLayoutGetRanges(layout, &ranges));
142   for (p = pStart; p < pEnd; ++p) {
143     PetscInt gdof, gcdof;
144 
145     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
146     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
147     PetscCheck(gcdof <= (gdof < 0 ? -(gdof + 1) : gdof), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " has %" PetscInt_FMT " constraints > %" PetscInt_FMT " dof", p, gcdof, (gdof < 0 ? -(gdof + 1) : gdof));
148     nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
149   }
150   PetscCall(PetscMalloc1(nleaves, &local));
151   PetscCall(PetscMalloc1(nleaves, &remote));
152   for (p = pStart, l = 0; p < pEnd; ++p) {
153     const PetscInt *cind;
154     PetscInt        dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
155 
156     PetscCall(PetscSectionGetDof(localSection, p, &dof));
157     PetscCall(PetscSectionGetOffset(localSection, p, &off));
158     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
159     PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind));
160     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
161     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
162     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
163     if (!gdof) continue; /* Censored point */
164     gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
165     if (gsize != dof - cdof) {
166       PetscCheck(gsize == dof, comm, PETSC_ERR_ARG_WRONG, "Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " is neither the constrained size %" PetscInt_FMT ", nor the unconstrained %" PetscInt_FMT, gsize, p, dof - cdof, dof);
167       cdof = 0; /* Ignore constraints */
168     }
169     for (d = 0, c = 0; d < dof; ++d) {
170       if ((c < cdof) && (cind[c] == d)) {
171         ++c;
172         continue;
173       }
174       local[l + d - c] = off + d;
175     }
176     PetscCheck(d - c == gsize, comm, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT ": Global dof %" PetscInt_FMT " != %" PetscInt_FMT " size - number of constraints", p, gsize, d - c);
177     if (gdof < 0) {
178       for (d = 0; d < gsize; ++d, ++l) {
179         PetscInt offset = -(goff + 1) + d, r;
180 
181         PetscCall(PetscFindInt(offset, size + 1, ranges, &r));
182         if (r < 0) r = -(r + 2);
183         PetscCheck(!(r < 0) && !(r >= size), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff);
184         remote[l].rank  = r;
185         remote[l].index = offset - ranges[r];
186       }
187     } else {
188       for (d = 0; d < gsize; ++d, ++l) {
189         remote[l].rank  = rank;
190         remote[l].index = goff + d - ranges[rank];
191       }
192     }
193   }
194   PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
195   PetscCall(PetscLayoutDestroy(&layout));
196   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
197   PetscFunctionReturn(0);
198 }
199 
200 /*@C
201   PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF
202 
203   Collective on sf
204 
205   Input Parameters:
206 + sf - The SF
207 - rootSection - Section defined on root space
208 
209   Output Parameters:
210 + remoteOffsets - root offsets in leaf storage, or NULL
211 - leafSection - Section defined on the leaf space
212 
213   Level: advanced
214 
215 .seealso: `PetscSFCreate()`
216 @*/
217 PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection) {
218   PetscSF         embedSF;
219   const PetscInt *indices;
220   IS              selected;
221   PetscInt        numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c;
222   PetscBool      *sub, hasc;
223 
224   PetscFunctionBegin;
225   PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0));
226   PetscCall(PetscSectionGetNumFields(rootSection, &numFields));
227   if (numFields) {
228     IS perm;
229 
230     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
231        leafSection->perm. To keep this permutation set by the user, we grab
232        the reference before calling PetscSectionSetNumFields() and set it
233        back after. */
234     PetscCall(PetscSectionGetPermutation(leafSection, &perm));
235     PetscCall(PetscObjectReference((PetscObject)perm));
236     PetscCall(PetscSectionSetNumFields(leafSection, numFields));
237     PetscCall(PetscSectionSetPermutation(leafSection, perm));
238     PetscCall(ISDestroy(&perm));
239   }
240   PetscCall(PetscMalloc1(numFields + 2, &sub));
241   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
242   for (f = 0; f < numFields; ++f) {
243     PetscSectionSym sym, dsym = NULL;
244     const char     *name    = NULL;
245     PetscInt        numComp = 0;
246 
247     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
248     PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp));
249     PetscCall(PetscSectionGetFieldName(rootSection, f, &name));
250     PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym));
251     if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym));
252     PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp));
253     PetscCall(PetscSectionSetFieldName(leafSection, f, name));
254     PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym));
255     PetscCall(PetscSectionSymDestroy(&dsym));
256     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
257       PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name));
258       PetscCall(PetscSectionSetComponentName(leafSection, f, c, name));
259     }
260   }
261   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
262   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
263   rpEnd  = PetscMin(rpEnd, nroots);
264   rpEnd  = PetscMax(rpStart, rpEnd);
265   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
266   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
267   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf)));
268   if (sub[0]) {
269     PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
270     PetscCall(ISGetIndices(selected, &indices));
271     PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
272     PetscCall(ISRestoreIndices(selected, &indices));
273     PetscCall(ISDestroy(&selected));
274   } else {
275     PetscCall(PetscObjectReference((PetscObject)sf));
276     embedSF = sf;
277   }
278   PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd));
279   lpEnd++;
280 
281   PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd));
282 
283   /* Constrained dof section */
284   hasc = sub[1];
285   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]);
286 
287   /* Could fuse these at the cost of copies and extra allocation */
288   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE));
289   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE));
290   if (sub[1]) {
291     PetscCall(PetscSectionCheckConstraints_Private(rootSection));
292     PetscCall(PetscSectionCheckConstraints_Private(leafSection));
293     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
294     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
295   }
296   for (f = 0; f < numFields; ++f) {
297     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE));
298     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE));
299     if (sub[2 + f]) {
300       PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f]));
301       PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f]));
302       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
303       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
304     }
305   }
306   if (remoteOffsets) {
307     PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets));
308     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
309     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
310   }
311   PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection));
312   PetscCall(PetscSectionSetUp(leafSection));
313   if (hasc) { /* need to communicate bcIndices */
314     PetscSF   bcSF;
315     PetscInt *rOffBc;
316 
317     PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc));
318     if (sub[1]) {
319       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
320       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
321       PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF));
322       PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
323       PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
324       PetscCall(PetscSFDestroy(&bcSF));
325     }
326     for (f = 0; f < numFields; ++f) {
327       if (sub[2 + f]) {
328         PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
329         PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
330         PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF));
331         PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
332         PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
333         PetscCall(PetscSFDestroy(&bcSF));
334       }
335     }
336     PetscCall(PetscFree(rOffBc));
337   }
338   PetscCall(PetscSFDestroy(&embedSF));
339   PetscCall(PetscFree(sub));
340   PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0));
341   PetscFunctionReturn(0);
342 }
343 
344 /*@C
345   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
346 
347   Collective on sf
348 
349   Input Parameters:
350 + sf          - The SF
351 . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
352 - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
353 
354   Output Parameter:
355 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
356 
357   Level: developer
358 
359 .seealso: `PetscSFCreate()`
360 @*/
361 PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets) {
362   PetscSF         embedSF;
363   const PetscInt *indices;
364   IS              selected;
365   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
366 
367   PetscFunctionBegin;
368   *remoteOffsets = NULL;
369   PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL));
370   if (numRoots < 0) PetscFunctionReturn(0);
371   PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0));
372   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
373   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
374   PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
375   PetscCall(ISGetIndices(selected, &indices));
376   PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
377   PetscCall(ISRestoreIndices(selected, &indices));
378   PetscCall(ISDestroy(&selected));
379   PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets));
380   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
381   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
382   PetscCall(PetscSFDestroy(&embedSF));
383   PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0));
384   PetscFunctionReturn(0);
385 }
386 
387 /*@C
388   PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points
389 
390   Collective on sf
391 
392   Input Parameters:
393 + sf - The SF
394 . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
395 . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
396 - leafSection - Data layout of local points for incoming data  (this is the distributed section)
397 
398   Output Parameters:
399 - sectionSF - The new SF
400 
401   Note: Either rootSection or remoteOffsets can be specified
402 
403   Level: advanced
404 
405 .seealso: `PetscSFCreate()`
406 @*/
407 PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF) {
408   MPI_Comm           comm;
409   const PetscInt    *localPoints;
410   const PetscSFNode *remotePoints;
411   PetscInt           lpStart, lpEnd;
412   PetscInt           numRoots, numSectionRoots, numPoints, numIndices = 0;
413   PetscInt          *localIndices;
414   PetscSFNode       *remoteIndices;
415   PetscInt           i, ind;
416 
417   PetscFunctionBegin;
418   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
419   PetscValidPointer(rootSection, 2);
420   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
421   PetscValidPointer(leafSection, 4);
422   PetscValidPointer(sectionSF, 5);
423   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
424   PetscCall(PetscSFCreate(comm, sectionSF));
425   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
426   PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots));
427   PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints));
428   if (numRoots < 0) PetscFunctionReturn(0);
429   PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0));
430   for (i = 0; i < numPoints; ++i) {
431     PetscInt localPoint = localPoints ? localPoints[i] : i;
432     PetscInt dof;
433 
434     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
435       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
436       numIndices += dof < 0 ? 0 : dof;
437     }
438   }
439   PetscCall(PetscMalloc1(numIndices, &localIndices));
440   PetscCall(PetscMalloc1(numIndices, &remoteIndices));
441   /* Create new index graph */
442   for (i = 0, ind = 0; i < numPoints; ++i) {
443     PetscInt localPoint = localPoints ? localPoints[i] : i;
444     PetscInt rank       = remotePoints[i].rank;
445 
446     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
447       PetscInt remoteOffset = remoteOffsets[localPoint - lpStart];
448       PetscInt loff, dof, d;
449 
450       PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff));
451       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
452       for (d = 0; d < dof; ++d, ++ind) {
453         localIndices[ind]        = loff + d;
454         remoteIndices[ind].rank  = rank;
455         remoteIndices[ind].index = remoteOffset + d;
456       }
457     }
458   }
459   PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
460   PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER));
461   PetscCall(PetscSFSetUp(*sectionSF));
462   PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0));
463   PetscFunctionReturn(0);
464 }
465 
466 /*@C
467    PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects
468 
469    Collective
470 
471    Input Parameters:
472 +  rmap - PetscLayout defining the global root space
473 -  lmap - PetscLayout defining the global leaf space
474 
475    Output Parameter:
476 .  sf - The parallel star forest
477 
478    Level: intermediate
479 
480 .seealso: `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()`
481 @*/
482 PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf) {
483   PetscInt     i, nroots, nleaves = 0;
484   PetscInt     rN, lst, len;
485   PetscMPIInt  owner = -1;
486   PetscSFNode *remote;
487   MPI_Comm     rcomm = rmap->comm;
488   MPI_Comm     lcomm = lmap->comm;
489   PetscMPIInt  flg;
490 
491   PetscFunctionBegin;
492   PetscValidPointer(sf, 3);
493   PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup");
494   PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup");
495   PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg));
496   PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators");
497   PetscCall(PetscSFCreate(rcomm, sf));
498   PetscCall(PetscLayoutGetLocalSize(rmap, &nroots));
499   PetscCall(PetscLayoutGetSize(rmap, &rN));
500   PetscCall(PetscLayoutGetRange(lmap, &lst, &len));
501   PetscCall(PetscMalloc1(len - lst, &remote));
502   for (i = lst; i < len && i < rN; i++) {
503     if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner));
504     remote[nleaves].rank  = owner;
505     remote[nleaves].index = i - rmap->range[owner];
506     nleaves++;
507   }
508   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES));
509   PetscCall(PetscFree(remote));
510   PetscFunctionReturn(0);
511 }
512 
513 /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
514 PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs) {
515   PetscInt    *owners = map->range;
516   PetscInt     n      = map->n;
517   PetscSF      sf;
518   PetscInt    *lidxs, *work = NULL;
519   PetscSFNode *ridxs;
520   PetscMPIInt  rank, p = 0;
521   PetscInt     r, len  = 0;
522 
523   PetscFunctionBegin;
524   if (on) *on = 0; /* squelch -Wmaybe-uninitialized */
525   /* Create SF where leaves are input idxs and roots are owned idxs */
526   PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
527   PetscCall(PetscMalloc1(n, &lidxs));
528   for (r = 0; r < n; ++r) lidxs[r] = -1;
529   PetscCall(PetscMalloc1(N, &ridxs));
530   for (r = 0; r < N; ++r) {
531     const PetscInt idx = idxs[r];
532     PetscCheck(idx >= 0 && idx < map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")", idx, map->N);
533     if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */
534       PetscCall(PetscLayoutFindOwner(map, idx, &p));
535     }
536     ridxs[r].rank  = p;
537     ridxs[r].index = idxs[r] - owners[p];
538   }
539   PetscCall(PetscSFCreate(map->comm, &sf));
540   PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER));
541   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
542   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
543   if (ogidxs) { /* communicate global idxs */
544     PetscInt cum = 0, start, *work2;
545 
546     PetscCall(PetscMalloc1(n, &work));
547     PetscCall(PetscCalloc1(N, &work2));
548     for (r = 0; r < N; ++r)
549       if (idxs[r] >= 0) cum++;
550     PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm));
551     start -= cum;
552     cum = 0;
553     for (r = 0; r < N; ++r)
554       if (idxs[r] >= 0) work2[r] = start + cum++;
555     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE));
556     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE));
557     PetscCall(PetscFree(work2));
558   }
559   PetscCall(PetscSFDestroy(&sf));
560   /* Compress and put in indices */
561   for (r = 0; r < n; ++r)
562     if (lidxs[r] >= 0) {
563       if (work) work[len] = work[r];
564       lidxs[len++] = r;
565     }
566   if (on) *on = len;
567   if (oidxs) *oidxs = lidxs;
568   if (ogidxs) *ogidxs = work;
569   PetscFunctionReturn(0);
570 }
571 
572 /*@
573   PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices
574 
575   Collective
576 
577   Input Parameters:
578 + layout - PetscLayout defining the global index space and the rank that brokers each index
579 . numRootIndices - size of rootIndices
580 . rootIndices - PetscInt array of global indices of which this process requests ownership
581 . rootLocalIndices - root local index permutation (NULL if no permutation)
582 . rootLocalOffset - offset to be added to root local indices
583 . numLeafIndices - size of leafIndices
584 . leafIndices - PetscInt array of global indices with which this process requires data associated
585 . leafLocalIndices - leaf local index permutation (NULL if no permutation)
586 - leafLocalOffset - offset to be added to leaf local indices
587 
588   Output Parameters:
589 + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed)
590 - sf - star forest representing the communication pattern from the root space to the leaf space
591 
592   Example 1:
593 $
594 $  rank             : 0            1            2
595 $  rootIndices      : [1 0 2]      [3]          [3]
596 $  rootLocalOffset  : 100          200          300
597 $  layout           : [0 1]        [2]          [3]
598 $  leafIndices      : [0]          [2]          [0 3]
599 $  leafLocalOffset  : 400          500          600
600 $
601 would build the following SF
602 $
603 $  [0] 400 <- (0,101)
604 $  [1] 500 <- (0,102)
605 $  [2] 600 <- (0,101)
606 $  [2] 601 <- (2,300)
607 $
608   Example 2:
609 $
610 $  rank             : 0               1               2
611 $  rootIndices      : [1 0 2]         [3]             [3]
612 $  rootLocalOffset  : 100             200             300
613 $  layout           : [0 1]           [2]             [3]
614 $  leafIndices      : rootIndices     rootIndices     rootIndices
615 $  leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset
616 $
617 would build the following SF
618 $
619 $  [1] 200 <- (2,300)
620 $
621   Example 3:
622 $
623 $  No process requests ownership of global index 1, but no process needs it.
624 $
625 $  rank             : 0            1            2
626 $  numRootIndices   : 2            1            1
627 $  rootIndices      : [0 2]        [3]          [3]
628 $  rootLocalOffset  : 100          200          300
629 $  layout           : [0 1]        [2]          [3]
630 $  numLeafIndices   : 1            1            2
631 $  leafIndices      : [0]          [2]          [0 3]
632 $  leafLocalOffset  : 400          500          600
633 $
634 would build the following SF
635 $
636 $  [0] 400 <- (0,100)
637 $  [1] 500 <- (0,101)
638 $  [2] 600 <- (0,100)
639 $  [2] 601 <- (2,300)
640 $
641 
642   Notes:
643   The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its
644   local size can be set to PETSC_DECIDE.
645   If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests
646   ownership of x and sends its own rank and the local index of x to process i.
647   If multiple processes request ownership of x, the one with the highest rank is to own x.
648   Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the
649   ownership information of x.
650   The output sf is constructed by associating each leaf point to a root point in this way.
651 
652   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
653   The optional output PetscSF object sfA can be used to push such data to leaf points.
654 
655   All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices
656   must cover that of leafIndices, but need not cover the entire layout.
657 
658   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
659   star forest is almost identity, so will only include non-trivial part of the map.
660 
661   Developer Notes:
662   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
663   hash(rank, root_local_index) as the bid for the ownership determination.
664 
665   Level: advanced
666 
667 .seealso: `PetscSFCreate()`
668 @*/
669 PetscErrorCode PetscSFCreateByMatchingIndices(PetscLayout layout, PetscInt numRootIndices, const PetscInt *rootIndices, const PetscInt *rootLocalIndices, PetscInt rootLocalOffset, PetscInt numLeafIndices, const PetscInt *leafIndices, const PetscInt *leafLocalIndices, PetscInt leafLocalOffset, PetscSF *sfA, PetscSF *sf) {
670   MPI_Comm     comm = layout->comm;
671   PetscMPIInt  size, rank;
672   PetscSF      sf1;
673   PetscSFNode *owners, *buffer, *iremote;
674   PetscInt    *ilocal, nleaves, N, n, i;
675 #if defined(PETSC_USE_DEBUG)
676   PetscInt N1;
677 #endif
678   PetscBool flag;
679 
680   PetscFunctionBegin;
681   if (rootIndices) PetscValidIntPointer(rootIndices, 3);
682   if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4);
683   if (leafIndices) PetscValidIntPointer(leafIndices, 7);
684   if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8);
685   if (sfA) PetscValidPointer(sfA, 10);
686   PetscValidPointer(sf, 11);
687   PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
688   PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
689   PetscCallMPI(MPI_Comm_size(comm, &size));
690   PetscCallMPI(MPI_Comm_rank(comm, &rank));
691   PetscCall(PetscLayoutSetUp(layout));
692   PetscCall(PetscLayoutGetSize(layout, &N));
693   PetscCall(PetscLayoutGetLocalSize(layout, &n));
694   flag = (PetscBool)(leafIndices == rootIndices);
695   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm));
696   PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
697 #if defined(PETSC_USE_DEBUG)
698   N1 = PETSC_MIN_INT;
699   for (i = 0; i < numRootIndices; i++)
700     if (rootIndices[i] > N1) N1 = rootIndices[i];
701   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
702   PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
703   if (!flag) {
704     N1 = PETSC_MIN_INT;
705     for (i = 0; i < numLeafIndices; i++)
706       if (leafIndices[i] > N1) N1 = leafIndices[i];
707     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
708     PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
709   }
710 #endif
711   /* Reduce: owners -> buffer */
712   PetscCall(PetscMalloc1(n, &buffer));
713   PetscCall(PetscSFCreate(comm, &sf1));
714   PetscCall(PetscSFSetFromOptions(sf1));
715   PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices));
716   PetscCall(PetscMalloc1(numRootIndices, &owners));
717   for (i = 0; i < numRootIndices; ++i) {
718     owners[i].rank  = rank;
719     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
720   }
721   for (i = 0; i < n; ++i) {
722     buffer[i].index = -1;
723     buffer[i].rank  = -1;
724   }
725   PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
726   PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
727   /* Bcast: buffer -> owners */
728   if (!flag) {
729     /* leafIndices is different from rootIndices */
730     PetscCall(PetscFree(owners));
731     PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices));
732     PetscCall(PetscMalloc1(numLeafIndices, &owners));
733   }
734   PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
735   PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
736   for (i = 0; i < numLeafIndices; ++i) PetscCheck(owners[i].rank >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global point %" PetscInt_FMT " was unclaimed", leafIndices[i]);
737   PetscCall(PetscFree(buffer));
738   if (sfA) {
739     *sfA = sf1;
740   } else PetscCall(PetscSFDestroy(&sf1));
741   /* Create sf */
742   if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
743     /* leaf space == root space */
744     for (i = 0, nleaves = 0; i < numLeafIndices; ++i)
745       if (owners[i].rank != rank) ++nleaves;
746     PetscCall(PetscMalloc1(nleaves, &ilocal));
747     PetscCall(PetscMalloc1(nleaves, &iremote));
748     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
749       if (owners[i].rank != rank) {
750         ilocal[nleaves]        = leafLocalOffset + i;
751         iremote[nleaves].rank  = owners[i].rank;
752         iremote[nleaves].index = owners[i].index;
753         ++nleaves;
754       }
755     }
756     PetscCall(PetscFree(owners));
757   } else {
758     nleaves = numLeafIndices;
759     PetscCall(PetscMalloc1(nleaves, &ilocal));
760     for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);
761     iremote = owners;
762   }
763   PetscCall(PetscSFCreate(comm, sf));
764   PetscCall(PetscSFSetFromOptions(*sf));
765   PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
766   PetscFunctionReturn(0);
767 }
768