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