xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision cab5436489e834181b0ebe29f8c17afa2aa8e653)
1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I  "petscsf.h"   I*/
2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h>
3b0c7db22SLisandro Dalcin 
4eb58282aSVaclav Hapla /*@
5*cab54364SBarry Smith    PetscSFSetGraphLayout - Set a parallel star forest via global indices and a `PetscLayout`
6b0c7db22SLisandro Dalcin 
7b0c7db22SLisandro Dalcin    Collective
8b0c7db22SLisandro Dalcin 
94165533cSJose E. Roman    Input Parameters:
10b0c7db22SLisandro Dalcin +  sf - star forest
11*cab54364SBarry Smith .  layout - `PetscLayout` defining the global space for roots
12b0c7db22SLisandro Dalcin .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
13b0c7db22SLisandro Dalcin .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
14b0c7db22SLisandro Dalcin .  localmode - copy mode for ilocal
15eb58282aSVaclav Hapla -  gremote - root vertices in global numbering corresponding to leaves in ilocal
16b0c7db22SLisandro Dalcin 
17b0c7db22SLisandro Dalcin    Level: intermediate
18b0c7db22SLisandro Dalcin 
19*cab54364SBarry Smith    Note:
2086c56f52SVaclav Hapla    Global indices must lie in [0, N) where N is the global size of layout.
21*cab54364SBarry Smith    Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`.
2286c56f52SVaclav Hapla 
23*cab54364SBarry Smith    Developer Note:
2486c56f52SVaclav Hapla    Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
25*cab54364SBarry Smith    encode contiguous storage. In such case, if localmode is `PETSC_OWN_POINTER`, the memory is deallocated as it is not
26b0c7db22SLisandro Dalcin    needed
27b0c7db22SLisandro Dalcin 
28*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
29b0c7db22SLisandro Dalcin @*/
30d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, const PetscInt *gremote)
31d71ae5a4SJacob Faibussowitsch {
3238a25198SStefano Zampini   const PetscInt *range;
3338a25198SStefano Zampini   PetscInt        i, nroots, ls = -1, ln = -1;
3438a25198SStefano Zampini   PetscMPIInt     lr = -1;
35b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
36b0c7db22SLisandro Dalcin 
37b0c7db22SLisandro Dalcin   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &nroots));
399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRanges(layout, &range));
409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &remote));
41eb58282aSVaclav Hapla   if (nleaves) ls = gremote[0] + 1;
42b0c7db22SLisandro Dalcin   for (i = 0; i < nleaves; i++) {
43eb58282aSVaclav Hapla     const PetscInt idx = gremote[i] - ls;
4438a25198SStefano Zampini     if (idx < 0 || idx >= ln) { /* short-circuit the search */
45eb58282aSVaclav Hapla       PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index));
4638a25198SStefano Zampini       remote[i].rank = lr;
4738a25198SStefano Zampini       ls             = range[lr];
4838a25198SStefano Zampini       ln             = range[lr + 1] - ls;
4938a25198SStefano Zampini     } else {
5038a25198SStefano Zampini       remote[i].rank  = lr;
5138a25198SStefano Zampini       remote[i].index = idx;
5238a25198SStefano Zampini     }
53b0c7db22SLisandro Dalcin   }
549566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER));
55b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
56b0c7db22SLisandro Dalcin }
57b0c7db22SLisandro Dalcin 
58eb58282aSVaclav Hapla /*@C
59*cab54364SBarry Smith    PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe this star forest
60eb58282aSVaclav Hapla 
61eb58282aSVaclav Hapla    Collective
62eb58282aSVaclav Hapla 
63eb58282aSVaclav Hapla    Input Parameter:
64eb58282aSVaclav Hapla .  sf - star forest
65eb58282aSVaclav Hapla 
66eb58282aSVaclav Hapla    Output Parameters:
67*cab54364SBarry Smith +  layout - `PetscLayout` defining the global space for roots
68eb58282aSVaclav Hapla .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
69eb58282aSVaclav Hapla .  ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage
70eb58282aSVaclav Hapla -  gremote - root vertices in global numbering corresponding to leaves in ilocal
71eb58282aSVaclav Hapla 
72eb58282aSVaclav Hapla    Level: intermediate
73eb58282aSVaclav Hapla 
74eb58282aSVaclav Hapla    Notes:
75eb58282aSVaclav Hapla    The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest.
76eb58282aSVaclav Hapla    The outputs layout and gremote are freshly created each time this function is called,
77eb58282aSVaclav Hapla    so they need to be freed by user and cannot be qualified as const.
78eb58282aSVaclav Hapla 
79*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
80eb58282aSVaclav Hapla @*/
81d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[])
82d71ae5a4SJacob Faibussowitsch {
83eb58282aSVaclav Hapla   PetscInt           nr, nl;
84eb58282aSVaclav Hapla   const PetscSFNode *ir;
85eb58282aSVaclav Hapla   PetscLayout        lt;
86eb58282aSVaclav Hapla 
87eb58282aSVaclav Hapla   PetscFunctionBegin;
88eb58282aSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir));
89eb58282aSVaclav Hapla   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, &lt));
90eb58282aSVaclav Hapla   if (gremote) {
91eb58282aSVaclav Hapla     PetscInt        i;
92eb58282aSVaclav Hapla     const PetscInt *range;
93eb58282aSVaclav Hapla     PetscInt       *gr;
94eb58282aSVaclav Hapla 
95eb58282aSVaclav Hapla     PetscCall(PetscLayoutGetRanges(lt, &range));
96eb58282aSVaclav Hapla     PetscCall(PetscMalloc1(nl, &gr));
97eb58282aSVaclav Hapla     for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index;
98eb58282aSVaclav Hapla     *gremote = gr;
99eb58282aSVaclav Hapla   }
100eb58282aSVaclav Hapla   if (nleaves) *nleaves = nl;
101eb58282aSVaclav Hapla   if (layout) *layout = lt;
102eb58282aSVaclav Hapla   else PetscCall(PetscLayoutDestroy(&lt));
103eb58282aSVaclav Hapla   PetscFunctionReturn(0);
104eb58282aSVaclav Hapla }
105eb58282aSVaclav Hapla 
106b0c7db22SLisandro Dalcin /*@
107*cab54364SBarry Smith   PetscSFSetGraphSection - Sets the `PetscSF` graph encoding the parallel dof overlap based upon the `PetscSection` describing the data layout.
108b0c7db22SLisandro Dalcin 
109b0c7db22SLisandro Dalcin   Input Parameters:
110*cab54364SBarry Smith + sf - The `PetscSF`
111*cab54364SBarry Smith . localSection - `PetscSection` describing the local data layout
112*cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout
113b0c7db22SLisandro Dalcin 
114b0c7db22SLisandro Dalcin   Level: developer
115b0c7db22SLisandro Dalcin 
116*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()`
117b0c7db22SLisandro Dalcin @*/
118d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
119d71ae5a4SJacob Faibussowitsch {
120b0c7db22SLisandro Dalcin   MPI_Comm        comm;
121b0c7db22SLisandro Dalcin   PetscLayout     layout;
122b0c7db22SLisandro Dalcin   const PetscInt *ranges;
123b0c7db22SLisandro Dalcin   PetscInt       *local;
124b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
125b0c7db22SLisandro Dalcin   PetscInt        pStart, pEnd, p, nroots, nleaves = 0, l;
126b0c7db22SLisandro Dalcin   PetscMPIInt     size, rank;
127b0c7db22SLisandro Dalcin 
128b0c7db22SLisandro Dalcin   PetscFunctionBegin;
129b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
130b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
131b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
132b0c7db22SLisandro Dalcin 
1339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
1349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1369566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
1379566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
1389566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &layout));
1399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(layout, 1));
1409566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
1419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
1429566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRanges(layout, &ranges));
143b0c7db22SLisandro Dalcin   for (p = pStart; p < pEnd; ++p) {
144b0c7db22SLisandro Dalcin     PetscInt gdof, gcdof;
145b0c7db22SLisandro Dalcin 
1469566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
1479566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
148c9cc58a2SBarry Smith     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));
149b0c7db22SLisandro Dalcin     nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
150b0c7db22SLisandro Dalcin   }
1519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &local));
1529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &remote));
153b0c7db22SLisandro Dalcin   for (p = pStart, l = 0; p < pEnd; ++p) {
154b0c7db22SLisandro Dalcin     const PetscInt *cind;
155b0c7db22SLisandro Dalcin     PetscInt        dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
156b0c7db22SLisandro Dalcin 
1579566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(localSection, p, &dof));
1589566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(localSection, p, &off));
1599566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
1609566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind));
1619566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
1629566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
1639566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
164b0c7db22SLisandro Dalcin     if (!gdof) continue; /* Censored point */
165b0c7db22SLisandro Dalcin     gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
166b0c7db22SLisandro Dalcin     if (gsize != dof - cdof) {
16708401ef6SPierre Jolivet       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);
168b0c7db22SLisandro Dalcin       cdof = 0; /* Ignore constraints */
169b0c7db22SLisandro Dalcin     }
170b0c7db22SLisandro Dalcin     for (d = 0, c = 0; d < dof; ++d) {
1719371c9d4SSatish Balay       if ((c < cdof) && (cind[c] == d)) {
1729371c9d4SSatish Balay         ++c;
1739371c9d4SSatish Balay         continue;
1749371c9d4SSatish Balay       }
175b0c7db22SLisandro Dalcin       local[l + d - c] = off + d;
176b0c7db22SLisandro Dalcin     }
17708401ef6SPierre Jolivet     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);
178b0c7db22SLisandro Dalcin     if (gdof < 0) {
179b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
180b0c7db22SLisandro Dalcin         PetscInt offset = -(goff + 1) + d, r;
181b0c7db22SLisandro Dalcin 
1829566063dSJacob Faibussowitsch         PetscCall(PetscFindInt(offset, size + 1, ranges, &r));
183b0c7db22SLisandro Dalcin         if (r < 0) r = -(r + 2);
18408401ef6SPierre Jolivet         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);
185b0c7db22SLisandro Dalcin         remote[l].rank  = r;
186b0c7db22SLisandro Dalcin         remote[l].index = offset - ranges[r];
187b0c7db22SLisandro Dalcin       }
188b0c7db22SLisandro Dalcin     } else {
189b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
190b0c7db22SLisandro Dalcin         remote[l].rank  = rank;
191b0c7db22SLisandro Dalcin         remote[l].index = goff + d - ranges[rank];
192b0c7db22SLisandro Dalcin       }
193b0c7db22SLisandro Dalcin     }
194b0c7db22SLisandro Dalcin   }
19508401ef6SPierre Jolivet   PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
1969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
1979566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
198b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
199b0c7db22SLisandro Dalcin }
200b0c7db22SLisandro Dalcin 
201b0c7db22SLisandro Dalcin /*@C
202*cab54364SBarry Smith   PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF`
203b0c7db22SLisandro Dalcin 
204b0c7db22SLisandro Dalcin   Collective on sf
205b0c7db22SLisandro Dalcin 
206b0c7db22SLisandro Dalcin   Input Parameters:
207*cab54364SBarry Smith + sf - The `PetscSF`
208b0c7db22SLisandro Dalcin - rootSection - Section defined on root space
209b0c7db22SLisandro Dalcin 
210b0c7db22SLisandro Dalcin   Output Parameters:
211b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL
212b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space
213b0c7db22SLisandro Dalcin 
214b0c7db22SLisandro Dalcin   Level: advanced
215b0c7db22SLisandro Dalcin 
216*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
217b0c7db22SLisandro Dalcin @*/
218d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
219d71ae5a4SJacob Faibussowitsch {
220b0c7db22SLisandro Dalcin   PetscSF         embedSF;
221b0c7db22SLisandro Dalcin   const PetscInt *indices;
222b0c7db22SLisandro Dalcin   IS              selected;
223b0c7db22SLisandro Dalcin   PetscInt        numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c;
224b0c7db22SLisandro Dalcin   PetscBool      *sub, hasc;
225b0c7db22SLisandro Dalcin 
226b0c7db22SLisandro Dalcin   PetscFunctionBegin;
2279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0));
2289566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(rootSection, &numFields));
229029e8818Sksagiyam   if (numFields) {
230029e8818Sksagiyam     IS perm;
231029e8818Sksagiyam 
232029e8818Sksagiyam     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
233029e8818Sksagiyam        leafSection->perm. To keep this permutation set by the user, we grab
234029e8818Sksagiyam        the reference before calling PetscSectionSetNumFields() and set it
235029e8818Sksagiyam        back after. */
2369566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetPermutation(leafSection, &perm));
2379566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)perm));
2389566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetNumFields(leafSection, numFields));
2399566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetPermutation(leafSection, perm));
2409566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&perm));
241029e8818Sksagiyam   }
2429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numFields + 2, &sub));
243b0c7db22SLisandro Dalcin   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
244b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
2452ee82556SMatthew G. Knepley     PetscSectionSym sym, dsym = NULL;
246b0c7db22SLisandro Dalcin     const char     *name    = NULL;
247b0c7db22SLisandro Dalcin     PetscInt        numComp = 0;
248b0c7db22SLisandro Dalcin 
249b0c7db22SLisandro Dalcin     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
2509566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp));
2519566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(rootSection, f, &name));
2529566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym));
2539566063dSJacob Faibussowitsch     if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym));
2549566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp));
2559566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(leafSection, f, name));
2569566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym));
2579566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymDestroy(&dsym));
258b0c7db22SLisandro Dalcin     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
2599566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name));
2609566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(leafSection, f, c, name));
261b0c7db22SLisandro Dalcin     }
262b0c7db22SLisandro Dalcin   }
2639566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
2649566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
265b0c7db22SLisandro Dalcin   rpEnd = PetscMin(rpEnd, nroots);
266b0c7db22SLisandro Dalcin   rpEnd = PetscMax(rpStart, rpEnd);
267b0c7db22SLisandro Dalcin   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
268b0c7db22SLisandro Dalcin   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
2691c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf)));
270b0c7db22SLisandro Dalcin   if (sub[0]) {
2719566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
2729566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(selected, &indices));
2739566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
2749566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(selected, &indices));
2759566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&selected));
276b0c7db22SLisandro Dalcin   } else {
2779566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)sf));
278b0c7db22SLisandro Dalcin     embedSF = sf;
279b0c7db22SLisandro Dalcin   }
2809566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd));
281b0c7db22SLisandro Dalcin   lpEnd++;
282b0c7db22SLisandro Dalcin 
2839566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd));
284b0c7db22SLisandro Dalcin 
285b0c7db22SLisandro Dalcin   /* Constrained dof section */
286b0c7db22SLisandro Dalcin   hasc = sub[1];
287b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]);
288b0c7db22SLisandro Dalcin 
289b0c7db22SLisandro Dalcin   /* Could fuse these at the cost of copies and extra allocation */
2909566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE));
2919566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE));
292b0c7db22SLisandro Dalcin   if (sub[1]) {
2937c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(rootSection));
2947c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(leafSection));
2959566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
2969566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
297b0c7db22SLisandro Dalcin   }
298b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
2999566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE));
3009566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE));
301b0c7db22SLisandro Dalcin     if (sub[2 + f]) {
3027c0883d5SVaclav Hapla       PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f]));
3037c0883d5SVaclav Hapla       PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f]));
3049566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
3059566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
306b0c7db22SLisandro Dalcin     }
307b0c7db22SLisandro Dalcin   }
308b0c7db22SLisandro Dalcin   if (remoteOffsets) {
3099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets));
3109566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
3119566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
312b0c7db22SLisandro Dalcin   }
31369c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection));
3149566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSection));
315b0c7db22SLisandro Dalcin   if (hasc) { /* need to communicate bcIndices */
316b0c7db22SLisandro Dalcin     PetscSF   bcSF;
317b0c7db22SLisandro Dalcin     PetscInt *rOffBc;
318b0c7db22SLisandro Dalcin 
3199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc));
320b0c7db22SLisandro Dalcin     if (sub[1]) {
3219566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3229566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3239566063dSJacob Faibussowitsch       PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF));
3249566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
3259566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
3269566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&bcSF));
327b0c7db22SLisandro Dalcin     }
328b0c7db22SLisandro Dalcin     for (f = 0; f < numFields; ++f) {
329b0c7db22SLisandro Dalcin       if (sub[2 + f]) {
3309566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3319566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3329566063dSJacob Faibussowitsch         PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF));
3339566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
3349566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
3359566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&bcSF));
336b0c7db22SLisandro Dalcin       }
337b0c7db22SLisandro Dalcin     }
3389566063dSJacob Faibussowitsch     PetscCall(PetscFree(rOffBc));
339b0c7db22SLisandro Dalcin   }
3409566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&embedSF));
3419566063dSJacob Faibussowitsch   PetscCall(PetscFree(sub));
3429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0));
343b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
344b0c7db22SLisandro Dalcin }
345b0c7db22SLisandro Dalcin 
346b0c7db22SLisandro Dalcin /*@C
347b0c7db22SLisandro Dalcin   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
348b0c7db22SLisandro Dalcin 
349b0c7db22SLisandro Dalcin   Collective on sf
350b0c7db22SLisandro Dalcin 
351b0c7db22SLisandro Dalcin   Input Parameters:
352*cab54364SBarry Smith + sf          - The `PetscSF`
353b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
354b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
355b0c7db22SLisandro Dalcin 
356b0c7db22SLisandro Dalcin   Output Parameter:
357b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
358b0c7db22SLisandro Dalcin 
359b0c7db22SLisandro Dalcin   Level: developer
360b0c7db22SLisandro Dalcin 
361*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
362b0c7db22SLisandro Dalcin @*/
363d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
364d71ae5a4SJacob Faibussowitsch {
365b0c7db22SLisandro Dalcin   PetscSF         embedSF;
366b0c7db22SLisandro Dalcin   const PetscInt *indices;
367b0c7db22SLisandro Dalcin   IS              selected;
368b0c7db22SLisandro Dalcin   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
369b0c7db22SLisandro Dalcin 
370b0c7db22SLisandro Dalcin   PetscFunctionBegin;
371b0c7db22SLisandro Dalcin   *remoteOffsets = NULL;
3729566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL));
373b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
3749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0));
3759566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
3769566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
3779566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
3789566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected, &indices));
3799566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
3809566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(selected, &indices));
3819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selected));
3829566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets));
3839566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
3849566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
3859566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&embedSF));
3869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0));
387b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
388b0c7db22SLisandro Dalcin }
389b0c7db22SLisandro Dalcin 
390b0c7db22SLisandro Dalcin /*@C
391*cab54364SBarry Smith   PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points
392b0c7db22SLisandro Dalcin 
393b0c7db22SLisandro Dalcin   Collective on sf
394b0c7db22SLisandro Dalcin 
395b0c7db22SLisandro Dalcin   Input Parameters:
396*cab54364SBarry Smith + sf - The `PetscSF`
397b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
398b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
399b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is the distributed section)
400b0c7db22SLisandro Dalcin 
401b0c7db22SLisandro Dalcin   Output Parameters:
402*cab54364SBarry Smith - sectionSF - The new `PetscSF`
403b0c7db22SLisandro Dalcin 
404b0c7db22SLisandro Dalcin   Level: advanced
405b0c7db22SLisandro Dalcin 
406*cab54364SBarry Smith   Note:
407*cab54364SBarry Smith   Either rootSection or remoteOffsets can be specified
408*cab54364SBarry Smith 
409*cab54364SBarry Smith .seealso:  `PetscSF`, `PetscSFCreate()`
410b0c7db22SLisandro Dalcin @*/
411d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
412d71ae5a4SJacob Faibussowitsch {
413b0c7db22SLisandro Dalcin   MPI_Comm           comm;
414b0c7db22SLisandro Dalcin   const PetscInt    *localPoints;
415b0c7db22SLisandro Dalcin   const PetscSFNode *remotePoints;
416b0c7db22SLisandro Dalcin   PetscInt           lpStart, lpEnd;
417b0c7db22SLisandro Dalcin   PetscInt           numRoots, numSectionRoots, numPoints, numIndices = 0;
418b0c7db22SLisandro Dalcin   PetscInt          *localIndices;
419b0c7db22SLisandro Dalcin   PetscSFNode       *remoteIndices;
420b0c7db22SLisandro Dalcin   PetscInt           i, ind;
421b0c7db22SLisandro Dalcin 
422b0c7db22SLisandro Dalcin   PetscFunctionBegin;
423b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
424b0c7db22SLisandro Dalcin   PetscValidPointer(rootSection, 2);
425b0c7db22SLisandro Dalcin   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
426b0c7db22SLisandro Dalcin   PetscValidPointer(leafSection, 4);
427b0c7db22SLisandro Dalcin   PetscValidPointer(sectionSF, 5);
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
4299566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, sectionSF));
4309566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
4319566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots));
4329566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints));
433b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
4349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0));
435b0c7db22SLisandro Dalcin   for (i = 0; i < numPoints; ++i) {
436b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
437b0c7db22SLisandro Dalcin     PetscInt dof;
438b0c7db22SLisandro Dalcin 
439b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
4409566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
44122a18c9fSMatthew G. Knepley       numIndices += dof < 0 ? 0 : dof;
442b0c7db22SLisandro Dalcin     }
443b0c7db22SLisandro Dalcin   }
4449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numIndices, &localIndices));
4459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numIndices, &remoteIndices));
446b0c7db22SLisandro Dalcin   /* Create new index graph */
447b0c7db22SLisandro Dalcin   for (i = 0, ind = 0; i < numPoints; ++i) {
448b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
449b0c7db22SLisandro Dalcin     PetscInt rank       = remotePoints[i].rank;
450b0c7db22SLisandro Dalcin 
451b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
452b0c7db22SLisandro Dalcin       PetscInt remoteOffset = remoteOffsets[localPoint - lpStart];
453b0c7db22SLisandro Dalcin       PetscInt loff, dof, d;
454b0c7db22SLisandro Dalcin 
4559566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff));
4569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
457b0c7db22SLisandro Dalcin       for (d = 0; d < dof; ++d, ++ind) {
458b0c7db22SLisandro Dalcin         localIndices[ind]        = loff + d;
459b0c7db22SLisandro Dalcin         remoteIndices[ind].rank  = rank;
460b0c7db22SLisandro Dalcin         remoteIndices[ind].index = remoteOffset + d;
461b0c7db22SLisandro Dalcin       }
462b0c7db22SLisandro Dalcin     }
463b0c7db22SLisandro Dalcin   }
46408401ef6SPierre Jolivet   PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
4659566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER));
4669566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(*sectionSF));
4679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0));
468b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
469b0c7db22SLisandro Dalcin }
4708fa9e22eSVaclav Hapla 
4718fa9e22eSVaclav Hapla /*@C
472*cab54364SBarry Smith    PetscSFCreateFromLayouts - Creates a parallel star forest mapping two `PetscLayout` objects
4738fa9e22eSVaclav Hapla 
4748fa9e22eSVaclav Hapla    Collective
4758fa9e22eSVaclav Hapla 
4764165533cSJose E. Roman    Input Parameters:
477*cab54364SBarry Smith +  rmap - `PetscLayout` defining the global root space
478*cab54364SBarry Smith -  lmap - `PetscLayout` defining the global leaf space
4798fa9e22eSVaclav Hapla 
4804165533cSJose E. Roman    Output Parameter:
4818fa9e22eSVaclav Hapla .  sf - The parallel star forest
4828fa9e22eSVaclav Hapla 
4838fa9e22eSVaclav Hapla    Level: intermediate
4848fa9e22eSVaclav Hapla 
485*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()`
4868fa9e22eSVaclav Hapla @*/
487d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf)
488d71ae5a4SJacob Faibussowitsch {
4898fa9e22eSVaclav Hapla   PetscInt     i, nroots, nleaves = 0;
4908fa9e22eSVaclav Hapla   PetscInt     rN, lst, len;
4918fa9e22eSVaclav Hapla   PetscMPIInt  owner = -1;
4928fa9e22eSVaclav Hapla   PetscSFNode *remote;
4938fa9e22eSVaclav Hapla   MPI_Comm     rcomm = rmap->comm;
4948fa9e22eSVaclav Hapla   MPI_Comm     lcomm = lmap->comm;
4958fa9e22eSVaclav Hapla   PetscMPIInt  flg;
4968fa9e22eSVaclav Hapla 
4978fa9e22eSVaclav Hapla   PetscFunctionBegin;
4988fa9e22eSVaclav Hapla   PetscValidPointer(sf, 3);
49928b400f6SJacob Faibussowitsch   PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup");
50028b400f6SJacob Faibussowitsch   PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup");
5019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg));
502c9cc58a2SBarry Smith   PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators");
5039566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(rcomm, sf));
5049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(rmap, &nroots));
5059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(rmap, &rN));
5069566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(lmap, &lst, &len));
5079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len - lst, &remote));
5088fa9e22eSVaclav Hapla   for (i = lst; i < len && i < rN; i++) {
50948a46eb9SPierre Jolivet     if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner));
5108fa9e22eSVaclav Hapla     remote[nleaves].rank  = owner;
5118fa9e22eSVaclav Hapla     remote[nleaves].index = i - rmap->range[owner];
5128fa9e22eSVaclav Hapla     nleaves++;
5138fa9e22eSVaclav Hapla   }
5149566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES));
5159566063dSJacob Faibussowitsch   PetscCall(PetscFree(remote));
5168fa9e22eSVaclav Hapla   PetscFunctionReturn(0);
5178fa9e22eSVaclav Hapla }
5188fa9e22eSVaclav Hapla 
5198fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
520d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs)
521d71ae5a4SJacob Faibussowitsch {
5228fa9e22eSVaclav Hapla   PetscInt    *owners = map->range;
5238fa9e22eSVaclav Hapla   PetscInt     n      = map->n;
5248fa9e22eSVaclav Hapla   PetscSF      sf;
5258fa9e22eSVaclav Hapla   PetscInt    *lidxs, *work = NULL;
5268fa9e22eSVaclav Hapla   PetscSFNode *ridxs;
5278fa9e22eSVaclav Hapla   PetscMPIInt  rank, p = 0;
5288fa9e22eSVaclav Hapla   PetscInt     r, len  = 0;
5298fa9e22eSVaclav Hapla 
5308fa9e22eSVaclav Hapla   PetscFunctionBegin;
5318fa9e22eSVaclav Hapla   if (on) *on = 0; /* squelch -Wmaybe-uninitialized */
5328fa9e22eSVaclav Hapla   /* Create SF where leaves are input idxs and roots are owned idxs */
5339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
5349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &lidxs));
5358fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r) lidxs[r] = -1;
5369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N, &ridxs));
5378fa9e22eSVaclav Hapla   for (r = 0; r < N; ++r) {
5388fa9e22eSVaclav Hapla     const PetscInt idx = idxs[r];
539c9cc58a2SBarry Smith     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);
5408fa9e22eSVaclav Hapla     if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */
5419566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(map, idx, &p));
5428fa9e22eSVaclav Hapla     }
5438fa9e22eSVaclav Hapla     ridxs[r].rank  = p;
5448fa9e22eSVaclav Hapla     ridxs[r].index = idxs[r] - owners[p];
5458fa9e22eSVaclav Hapla   }
5469566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(map->comm, &sf));
5479566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER));
5489566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
5499566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
5508fa9e22eSVaclav Hapla   if (ogidxs) { /* communicate global idxs */
5518fa9e22eSVaclav Hapla     PetscInt cum = 0, start, *work2;
5528fa9e22eSVaclav Hapla 
5539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &work));
5549566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(N, &work2));
5559371c9d4SSatish Balay     for (r = 0; r < N; ++r)
5569371c9d4SSatish Balay       if (idxs[r] >= 0) cum++;
5579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm));
5588fa9e22eSVaclav Hapla     start -= cum;
5598fa9e22eSVaclav Hapla     cum = 0;
5609371c9d4SSatish Balay     for (r = 0; r < N; ++r)
5619371c9d4SSatish Balay       if (idxs[r] >= 0) work2[r] = start + cum++;
5629566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE));
5639566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE));
5649566063dSJacob Faibussowitsch     PetscCall(PetscFree(work2));
5658fa9e22eSVaclav Hapla   }
5669566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5678fa9e22eSVaclav Hapla   /* Compress and put in indices */
5688fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r)
5698fa9e22eSVaclav Hapla     if (lidxs[r] >= 0) {
5708fa9e22eSVaclav Hapla       if (work) work[len] = work[r];
5718fa9e22eSVaclav Hapla       lidxs[len++] = r;
5728fa9e22eSVaclav Hapla     }
5738fa9e22eSVaclav Hapla   if (on) *on = len;
5748fa9e22eSVaclav Hapla   if (oidxs) *oidxs = lidxs;
5758fa9e22eSVaclav Hapla   if (ogidxs) *ogidxs = work;
5768fa9e22eSVaclav Hapla   PetscFunctionReturn(0);
5778fa9e22eSVaclav Hapla }
578deffd5ebSksagiyam 
579deffd5ebSksagiyam /*@
580*cab54364SBarry Smith   PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices
581deffd5ebSksagiyam 
582deffd5ebSksagiyam   Collective
583deffd5ebSksagiyam 
584deffd5ebSksagiyam   Input Parameters:
585*cab54364SBarry Smith + layout - `PetscLayout` defining the global index space and the rank that brokers each index
586deffd5ebSksagiyam . numRootIndices - size of rootIndices
587*cab54364SBarry Smith . rootIndices - `PetscInt` array of global indices of which this process requests ownership
588deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation)
589deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices
590deffd5ebSksagiyam . numLeafIndices - size of leafIndices
591*cab54364SBarry Smith . leafIndices - `PetscInt` array of global indices with which this process requires data associated
592deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation)
593deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices
594deffd5ebSksagiyam 
595d8d19677SJose E. Roman   Output Parameters:
596deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed)
597deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space
598deffd5ebSksagiyam 
599*cab54364SBarry Smith   Level: advanced
600*cab54364SBarry Smith 
601deffd5ebSksagiyam   Example 1:
602*cab54364SBarry Smith .vb
603*cab54364SBarry Smith   rank             : 0            1            2
604*cab54364SBarry Smith   rootIndices      : [1 0 2]      [3]          [3]
605*cab54364SBarry Smith   rootLocalOffset  : 100          200          300
606*cab54364SBarry Smith   layout           : [0 1]        [2]          [3]
607*cab54364SBarry Smith   leafIndices      : [0]          [2]          [0 3]
608*cab54364SBarry Smith   leafLocalOffset  : 400          500          600
609*cab54364SBarry Smith 
610*cab54364SBarry Smith would build the following PetscSF
611*cab54364SBarry Smith 
612*cab54364SBarry Smith   [0] 400 <- (0,101)
613*cab54364SBarry Smith   [1] 500 <- (0,102)
614*cab54364SBarry Smith   [2] 600 <- (0,101)
615*cab54364SBarry Smith   [2] 601 <- (2,300)
616*cab54364SBarry Smith .ve
617*cab54364SBarry Smith 
618deffd5ebSksagiyam   Example 2:
619*cab54364SBarry Smith .vb
620*cab54364SBarry Smith   rank             : 0               1               2
621*cab54364SBarry Smith   rootIndices      : [1 0 2]         [3]             [3]
622*cab54364SBarry Smith   rootLocalOffset  : 100             200             300
623*cab54364SBarry Smith   layout           : [0 1]           [2]             [3]
624*cab54364SBarry Smith   leafIndices      : rootIndices     rootIndices     rootIndices
625*cab54364SBarry Smith   leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset
626*cab54364SBarry Smith 
627*cab54364SBarry Smith would build the following PetscSF
628*cab54364SBarry Smith 
629*cab54364SBarry Smith   [1] 200 <- (2,300)
630*cab54364SBarry Smith .ve
631*cab54364SBarry Smith 
632deffd5ebSksagiyam   Example 3:
633*cab54364SBarry Smith .vb
634*cab54364SBarry Smith   No process requests ownership of global index 1, but no process needs it.
635*cab54364SBarry Smith 
636*cab54364SBarry Smith   rank             : 0            1            2
637*cab54364SBarry Smith   numRootIndices   : 2            1            1
638*cab54364SBarry Smith   rootIndices      : [0 2]        [3]          [3]
639*cab54364SBarry Smith   rootLocalOffset  : 100          200          300
640*cab54364SBarry Smith   layout           : [0 1]        [2]          [3]
641*cab54364SBarry Smith   numLeafIndices   : 1            1            2
642*cab54364SBarry Smith   leafIndices      : [0]          [2]          [0 3]
643*cab54364SBarry Smith   leafLocalOffset  : 400          500          600
644*cab54364SBarry Smith 
645*cab54364SBarry Smith would build the following PetscSF
646*cab54364SBarry Smith 
647*cab54364SBarry Smith   [0] 400 <- (0,100)
648*cab54364SBarry Smith   [1] 500 <- (0,101)
649*cab54364SBarry Smith   [2] 600 <- (0,100)
650*cab54364SBarry Smith   [2] 601 <- (2,300)
651*cab54364SBarry Smith .ve
652deffd5ebSksagiyam 
653deffd5ebSksagiyam   Notes:
654deffd5ebSksagiyam   The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its
655*cab54364SBarry Smith   local size can be set to `PETSC_DECIDE`.
656*cab54364SBarry Smith 
657deffd5ebSksagiyam   If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests
658deffd5ebSksagiyam   ownership of x and sends its own rank and the local index of x to process i.
659deffd5ebSksagiyam   If multiple processes request ownership of x, the one with the highest rank is to own x.
660deffd5ebSksagiyam   Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the
661deffd5ebSksagiyam   ownership information of x.
662deffd5ebSksagiyam   The output sf is constructed by associating each leaf point to a root point in this way.
663deffd5ebSksagiyam 
664deffd5ebSksagiyam   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
665*cab54364SBarry Smith   The optional output `PetscSF` object sfA can be used to push such data to leaf points.
666deffd5ebSksagiyam 
667deffd5ebSksagiyam   All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices
668deffd5ebSksagiyam   must cover that of leafIndices, but need not cover the entire layout.
669deffd5ebSksagiyam 
670deffd5ebSksagiyam   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
671deffd5ebSksagiyam   star forest is almost identity, so will only include non-trivial part of the map.
672deffd5ebSksagiyam 
673*cab54364SBarry Smith   Developer Note:
674deffd5ebSksagiyam   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
675deffd5ebSksagiyam   hash(rank, root_local_index) as the bid for the ownership determination.
676deffd5ebSksagiyam 
677*cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
678deffd5ebSksagiyam @*/
679d71ae5a4SJacob Faibussowitsch 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)
680d71ae5a4SJacob Faibussowitsch {
681deffd5ebSksagiyam   MPI_Comm     comm = layout->comm;
682deffd5ebSksagiyam   PetscMPIInt  size, rank;
683deffd5ebSksagiyam   PetscSF      sf1;
684deffd5ebSksagiyam   PetscSFNode *owners, *buffer, *iremote;
685deffd5ebSksagiyam   PetscInt    *ilocal, nleaves, N, n, i;
686deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
687deffd5ebSksagiyam   PetscInt N1;
688deffd5ebSksagiyam #endif
689deffd5ebSksagiyam   PetscBool flag;
690deffd5ebSksagiyam 
691deffd5ebSksagiyam   PetscFunctionBegin;
692deffd5ebSksagiyam   if (rootIndices) PetscValidIntPointer(rootIndices, 3);
693deffd5ebSksagiyam   if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4);
694deffd5ebSksagiyam   if (leafIndices) PetscValidIntPointer(leafIndices, 7);
695deffd5ebSksagiyam   if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8);
696deffd5ebSksagiyam   if (sfA) PetscValidPointer(sfA, 10);
697deffd5ebSksagiyam   PetscValidPointer(sf, 11);
69808401ef6SPierre Jolivet   PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
69908401ef6SPierre Jolivet   PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
7009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
7039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(layout, &N));
7049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &n));
705deffd5ebSksagiyam   flag = (PetscBool)(leafIndices == rootIndices);
7061c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm));
707c9cc58a2SBarry Smith   PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
708deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
709deffd5ebSksagiyam   N1 = PETSC_MIN_INT;
7109371c9d4SSatish Balay   for (i = 0; i < numRootIndices; i++)
7119371c9d4SSatish Balay     if (rootIndices[i] > N1) N1 = rootIndices[i];
7129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
71308401ef6SPierre Jolivet   PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
714deffd5ebSksagiyam   if (!flag) {
715deffd5ebSksagiyam     N1 = PETSC_MIN_INT;
7169371c9d4SSatish Balay     for (i = 0; i < numLeafIndices; i++)
7179371c9d4SSatish Balay       if (leafIndices[i] > N1) N1 = leafIndices[i];
7189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
71908401ef6SPierre Jolivet     PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
720deffd5ebSksagiyam   }
721deffd5ebSksagiyam #endif
722deffd5ebSksagiyam   /* Reduce: owners -> buffer */
7239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &buffer));
7249566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, &sf1));
7259566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf1));
7269566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices));
7279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootIndices, &owners));
728deffd5ebSksagiyam   for (i = 0; i < numRootIndices; ++i) {
729deffd5ebSksagiyam     owners[i].rank  = rank;
730deffd5ebSksagiyam     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
731deffd5ebSksagiyam   }
732deffd5ebSksagiyam   for (i = 0; i < n; ++i) {
733deffd5ebSksagiyam     buffer[i].index = -1;
734deffd5ebSksagiyam     buffer[i].rank  = -1;
735deffd5ebSksagiyam   }
7369566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
7379566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
738deffd5ebSksagiyam   /* Bcast: buffer -> owners */
739deffd5ebSksagiyam   if (!flag) {
740deffd5ebSksagiyam     /* leafIndices is different from rootIndices */
7419566063dSJacob Faibussowitsch     PetscCall(PetscFree(owners));
7429566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices));
7439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeafIndices, &owners));
744deffd5ebSksagiyam   }
7459566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
7469566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
74708401ef6SPierre Jolivet   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]);
7489566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
7499371c9d4SSatish Balay   if (sfA) {
7509371c9d4SSatish Balay     *sfA = sf1;
7519371c9d4SSatish Balay   } else PetscCall(PetscSFDestroy(&sf1));
752deffd5ebSksagiyam   /* Create sf */
753deffd5ebSksagiyam   if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
754deffd5ebSksagiyam     /* leaf space == root space */
7559371c9d4SSatish Balay     for (i = 0, nleaves = 0; i < numLeafIndices; ++i)
7569371c9d4SSatish Balay       if (owners[i].rank != rank) ++nleaves;
7579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &ilocal));
7589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &iremote));
759deffd5ebSksagiyam     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
760deffd5ebSksagiyam       if (owners[i].rank != rank) {
761deffd5ebSksagiyam         ilocal[nleaves]        = leafLocalOffset + i;
762deffd5ebSksagiyam         iremote[nleaves].rank  = owners[i].rank;
763deffd5ebSksagiyam         iremote[nleaves].index = owners[i].index;
764deffd5ebSksagiyam         ++nleaves;
765deffd5ebSksagiyam       }
766deffd5ebSksagiyam     }
7679566063dSJacob Faibussowitsch     PetscCall(PetscFree(owners));
768deffd5ebSksagiyam   } else {
769deffd5ebSksagiyam     nleaves = numLeafIndices;
7709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &ilocal));
771ad540459SPierre Jolivet     for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);
772deffd5ebSksagiyam     iremote = owners;
773deffd5ebSksagiyam   }
7749566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, sf));
7759566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sf));
7769566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
777deffd5ebSksagiyam   PetscFunctionReturn(0);
778deffd5ebSksagiyam }
779