xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision b114984a0359a26ae515feffb05b49f8a558770e)
1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I  "petscsf.h"   I*/
2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h>
3b0c7db22SLisandro Dalcin 
4eb58282aSVaclav Hapla /*@
5cab54364SBarry 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
11cab54364SBarry 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 
19cab54364SBarry Smith    Note:
2086c56f52SVaclav Hapla    Global indices must lie in [0, N) where N is the global size of layout.
21cab54364SBarry Smith    Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`.
2286c56f52SVaclav Hapla 
23cab54364SBarry Smith    Developer Note:
2486c56f52SVaclav Hapla    Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
25cab54364SBarry 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 
28cab54364SBarry 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;
38*b114984aSVaclav Hapla   PetscCall(PetscLayoutSetUp(layout));
399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &nroots));
409566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRanges(layout, &range));
419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &remote));
42eb58282aSVaclav Hapla   if (nleaves) ls = gremote[0] + 1;
43b0c7db22SLisandro Dalcin   for (i = 0; i < nleaves; i++) {
44eb58282aSVaclav Hapla     const PetscInt idx = gremote[i] - ls;
4538a25198SStefano Zampini     if (idx < 0 || idx >= ln) { /* short-circuit the search */
46eb58282aSVaclav Hapla       PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index));
4738a25198SStefano Zampini       remote[i].rank = lr;
4838a25198SStefano Zampini       ls             = range[lr];
4938a25198SStefano Zampini       ln             = range[lr + 1] - ls;
5038a25198SStefano Zampini     } else {
5138a25198SStefano Zampini       remote[i].rank  = lr;
5238a25198SStefano Zampini       remote[i].index = idx;
5338a25198SStefano Zampini     }
54b0c7db22SLisandro Dalcin   }
559566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER));
56b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
57b0c7db22SLisandro Dalcin }
58b0c7db22SLisandro Dalcin 
59eb58282aSVaclav Hapla /*@C
60cab54364SBarry Smith    PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe this star forest
61eb58282aSVaclav Hapla 
62eb58282aSVaclav Hapla    Collective
63eb58282aSVaclav Hapla 
64eb58282aSVaclav Hapla    Input Parameter:
65eb58282aSVaclav Hapla .  sf - star forest
66eb58282aSVaclav Hapla 
67eb58282aSVaclav Hapla    Output Parameters:
68cab54364SBarry Smith +  layout - `PetscLayout` defining the global space for roots
69eb58282aSVaclav Hapla .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
70eb58282aSVaclav Hapla .  ilocal - locations of leaves in leafdata buffers, or NULL for contiguous storage
71eb58282aSVaclav Hapla -  gremote - root vertices in global numbering corresponding to leaves in ilocal
72eb58282aSVaclav Hapla 
73eb58282aSVaclav Hapla    Level: intermediate
74eb58282aSVaclav Hapla 
75eb58282aSVaclav Hapla    Notes:
76eb58282aSVaclav Hapla    The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest.
77eb58282aSVaclav Hapla    The outputs layout and gremote are freshly created each time this function is called,
78eb58282aSVaclav Hapla    so they need to be freed by user and cannot be qualified as const.
79eb58282aSVaclav Hapla 
80cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
81eb58282aSVaclav Hapla @*/
82d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[])
83d71ae5a4SJacob Faibussowitsch {
84eb58282aSVaclav Hapla   PetscInt           nr, nl;
85eb58282aSVaclav Hapla   const PetscSFNode *ir;
86eb58282aSVaclav Hapla   PetscLayout        lt;
87eb58282aSVaclav Hapla 
88eb58282aSVaclav Hapla   PetscFunctionBegin;
89eb58282aSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir));
90eb58282aSVaclav Hapla   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, &lt));
91eb58282aSVaclav Hapla   if (gremote) {
92eb58282aSVaclav Hapla     PetscInt        i;
93eb58282aSVaclav Hapla     const PetscInt *range;
94eb58282aSVaclav Hapla     PetscInt       *gr;
95eb58282aSVaclav Hapla 
96eb58282aSVaclav Hapla     PetscCall(PetscLayoutGetRanges(lt, &range));
97eb58282aSVaclav Hapla     PetscCall(PetscMalloc1(nl, &gr));
98eb58282aSVaclav Hapla     for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index;
99eb58282aSVaclav Hapla     *gremote = gr;
100eb58282aSVaclav Hapla   }
101eb58282aSVaclav Hapla   if (nleaves) *nleaves = nl;
102eb58282aSVaclav Hapla   if (layout) *layout = lt;
103eb58282aSVaclav Hapla   else PetscCall(PetscLayoutDestroy(&lt));
104eb58282aSVaclav Hapla   PetscFunctionReturn(0);
105eb58282aSVaclav Hapla }
106eb58282aSVaclav Hapla 
107b0c7db22SLisandro Dalcin /*@
108cab54364SBarry Smith   PetscSFSetGraphSection - Sets the `PetscSF` graph encoding the parallel dof overlap based upon the `PetscSection` describing the data layout.
109b0c7db22SLisandro Dalcin 
110b0c7db22SLisandro Dalcin   Input Parameters:
111cab54364SBarry Smith + sf - The `PetscSF`
112cab54364SBarry Smith . localSection - `PetscSection` describing the local data layout
113cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout
114b0c7db22SLisandro Dalcin 
115b0c7db22SLisandro Dalcin   Level: developer
116b0c7db22SLisandro Dalcin 
117cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()`
118b0c7db22SLisandro Dalcin @*/
119d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
120d71ae5a4SJacob Faibussowitsch {
121b0c7db22SLisandro Dalcin   MPI_Comm        comm;
122b0c7db22SLisandro Dalcin   PetscLayout     layout;
123b0c7db22SLisandro Dalcin   const PetscInt *ranges;
124b0c7db22SLisandro Dalcin   PetscInt       *local;
125b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
126b0c7db22SLisandro Dalcin   PetscInt        pStart, pEnd, p, nroots, nleaves = 0, l;
127b0c7db22SLisandro Dalcin   PetscMPIInt     size, rank;
128b0c7db22SLisandro Dalcin 
129b0c7db22SLisandro Dalcin   PetscFunctionBegin;
130b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
131b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
132b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
133b0c7db22SLisandro Dalcin 
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
1359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1379566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
1389566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
1399566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &layout));
1409566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(layout, 1));
1419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
1429566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
1439566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRanges(layout, &ranges));
144b0c7db22SLisandro Dalcin   for (p = pStart; p < pEnd; ++p) {
145b0c7db22SLisandro Dalcin     PetscInt gdof, gcdof;
146b0c7db22SLisandro Dalcin 
1479566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
1489566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
149c9cc58a2SBarry 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));
150b0c7db22SLisandro Dalcin     nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
151b0c7db22SLisandro Dalcin   }
1529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &local));
1539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &remote));
154b0c7db22SLisandro Dalcin   for (p = pStart, l = 0; p < pEnd; ++p) {
155b0c7db22SLisandro Dalcin     const PetscInt *cind;
156b0c7db22SLisandro Dalcin     PetscInt        dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
157b0c7db22SLisandro Dalcin 
1589566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(localSection, p, &dof));
1599566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(localSection, p, &off));
1609566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
1619566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind));
1629566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
1639566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
1649566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
165b0c7db22SLisandro Dalcin     if (!gdof) continue; /* Censored point */
166b0c7db22SLisandro Dalcin     gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
167b0c7db22SLisandro Dalcin     if (gsize != dof - cdof) {
16808401ef6SPierre 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);
169b0c7db22SLisandro Dalcin       cdof = 0; /* Ignore constraints */
170b0c7db22SLisandro Dalcin     }
171b0c7db22SLisandro Dalcin     for (d = 0, c = 0; d < dof; ++d) {
1729371c9d4SSatish Balay       if ((c < cdof) && (cind[c] == d)) {
1739371c9d4SSatish Balay         ++c;
1749371c9d4SSatish Balay         continue;
1759371c9d4SSatish Balay       }
176b0c7db22SLisandro Dalcin       local[l + d - c] = off + d;
177b0c7db22SLisandro Dalcin     }
17808401ef6SPierre 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);
179b0c7db22SLisandro Dalcin     if (gdof < 0) {
180b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
181b0c7db22SLisandro Dalcin         PetscInt offset = -(goff + 1) + d, r;
182b0c7db22SLisandro Dalcin 
1839566063dSJacob Faibussowitsch         PetscCall(PetscFindInt(offset, size + 1, ranges, &r));
184b0c7db22SLisandro Dalcin         if (r < 0) r = -(r + 2);
18508401ef6SPierre 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);
186b0c7db22SLisandro Dalcin         remote[l].rank  = r;
187b0c7db22SLisandro Dalcin         remote[l].index = offset - ranges[r];
188b0c7db22SLisandro Dalcin       }
189b0c7db22SLisandro Dalcin     } else {
190b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
191b0c7db22SLisandro Dalcin         remote[l].rank  = rank;
192b0c7db22SLisandro Dalcin         remote[l].index = goff + d - ranges[rank];
193b0c7db22SLisandro Dalcin       }
194b0c7db22SLisandro Dalcin     }
195b0c7db22SLisandro Dalcin   }
19608401ef6SPierre Jolivet   PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
1979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
1989566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
199b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
200b0c7db22SLisandro Dalcin }
201b0c7db22SLisandro Dalcin 
202b0c7db22SLisandro Dalcin /*@C
203cab54364SBarry Smith   PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF`
204b0c7db22SLisandro Dalcin 
205b0c7db22SLisandro Dalcin   Collective on sf
206b0c7db22SLisandro Dalcin 
207b0c7db22SLisandro Dalcin   Input Parameters:
208cab54364SBarry Smith + sf - The `PetscSF`
209b0c7db22SLisandro Dalcin - rootSection - Section defined on root space
210b0c7db22SLisandro Dalcin 
211b0c7db22SLisandro Dalcin   Output Parameters:
212b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL
213b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space
214b0c7db22SLisandro Dalcin 
215b0c7db22SLisandro Dalcin   Level: advanced
216b0c7db22SLisandro Dalcin 
217cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
218b0c7db22SLisandro Dalcin @*/
219d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
220d71ae5a4SJacob Faibussowitsch {
221b0c7db22SLisandro Dalcin   PetscSF         embedSF;
222b0c7db22SLisandro Dalcin   const PetscInt *indices;
223b0c7db22SLisandro Dalcin   IS              selected;
224b0c7db22SLisandro Dalcin   PetscInt        numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c;
225b0c7db22SLisandro Dalcin   PetscBool      *sub, hasc;
226b0c7db22SLisandro Dalcin 
227b0c7db22SLisandro Dalcin   PetscFunctionBegin;
2289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0));
2299566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(rootSection, &numFields));
230029e8818Sksagiyam   if (numFields) {
231029e8818Sksagiyam     IS perm;
232029e8818Sksagiyam 
233029e8818Sksagiyam     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
234029e8818Sksagiyam        leafSection->perm. To keep this permutation set by the user, we grab
235029e8818Sksagiyam        the reference before calling PetscSectionSetNumFields() and set it
236029e8818Sksagiyam        back after. */
2379566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetPermutation(leafSection, &perm));
2389566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)perm));
2399566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetNumFields(leafSection, numFields));
2409566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetPermutation(leafSection, perm));
2419566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&perm));
242029e8818Sksagiyam   }
2439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numFields + 2, &sub));
244b0c7db22SLisandro Dalcin   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
245b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
2462ee82556SMatthew G. Knepley     PetscSectionSym sym, dsym = NULL;
247b0c7db22SLisandro Dalcin     const char     *name    = NULL;
248b0c7db22SLisandro Dalcin     PetscInt        numComp = 0;
249b0c7db22SLisandro Dalcin 
250b0c7db22SLisandro Dalcin     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
2519566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp));
2529566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(rootSection, f, &name));
2539566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym));
2549566063dSJacob Faibussowitsch     if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym));
2559566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp));
2569566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(leafSection, f, name));
2579566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym));
2589566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymDestroy(&dsym));
259b0c7db22SLisandro Dalcin     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
2609566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name));
2619566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(leafSection, f, c, name));
262b0c7db22SLisandro Dalcin     }
263b0c7db22SLisandro Dalcin   }
2649566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
2659566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
266b0c7db22SLisandro Dalcin   rpEnd = PetscMin(rpEnd, nroots);
267b0c7db22SLisandro Dalcin   rpEnd = PetscMax(rpStart, rpEnd);
268b0c7db22SLisandro Dalcin   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
269b0c7db22SLisandro Dalcin   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
2701c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf)));
271b0c7db22SLisandro Dalcin   if (sub[0]) {
2729566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
2739566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(selected, &indices));
2749566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
2759566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(selected, &indices));
2769566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&selected));
277b0c7db22SLisandro Dalcin   } else {
2789566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)sf));
279b0c7db22SLisandro Dalcin     embedSF = sf;
280b0c7db22SLisandro Dalcin   }
2819566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd));
282b0c7db22SLisandro Dalcin   lpEnd++;
283b0c7db22SLisandro Dalcin 
2849566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd));
285b0c7db22SLisandro Dalcin 
286b0c7db22SLisandro Dalcin   /* Constrained dof section */
287b0c7db22SLisandro Dalcin   hasc = sub[1];
288b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]);
289b0c7db22SLisandro Dalcin 
290b0c7db22SLisandro Dalcin   /* Could fuse these at the cost of copies and extra allocation */
2919566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE));
2929566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart], MPI_REPLACE));
293b0c7db22SLisandro Dalcin   if (sub[1]) {
2947c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(rootSection));
2957c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(leafSection));
2969566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
2979566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
298b0c7db22SLisandro Dalcin   }
299b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
3009566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE));
3019566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart], MPI_REPLACE));
302b0c7db22SLisandro Dalcin     if (sub[2 + f]) {
3037c0883d5SVaclav Hapla       PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f]));
3047c0883d5SVaclav Hapla       PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f]));
3059566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
3069566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
307b0c7db22SLisandro Dalcin     }
308b0c7db22SLisandro Dalcin   }
309b0c7db22SLisandro Dalcin   if (remoteOffsets) {
3109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets));
3119566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
3129566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
313b0c7db22SLisandro Dalcin   }
31469c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection));
3159566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSection));
316b0c7db22SLisandro Dalcin   if (hasc) { /* need to communicate bcIndices */
317b0c7db22SLisandro Dalcin     PetscSF   bcSF;
318b0c7db22SLisandro Dalcin     PetscInt *rOffBc;
319b0c7db22SLisandro Dalcin 
3209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc));
321b0c7db22SLisandro Dalcin     if (sub[1]) {
3229566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3239566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3249566063dSJacob Faibussowitsch       PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF));
3259566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
3269566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
3279566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&bcSF));
328b0c7db22SLisandro Dalcin     }
329b0c7db22SLisandro Dalcin     for (f = 0; f < numFields; ++f) {
330b0c7db22SLisandro Dalcin       if (sub[2 + f]) {
3319566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3329566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3339566063dSJacob Faibussowitsch         PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF));
3349566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
3359566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
3369566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&bcSF));
337b0c7db22SLisandro Dalcin       }
338b0c7db22SLisandro Dalcin     }
3399566063dSJacob Faibussowitsch     PetscCall(PetscFree(rOffBc));
340b0c7db22SLisandro Dalcin   }
3419566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&embedSF));
3429566063dSJacob Faibussowitsch   PetscCall(PetscFree(sub));
3439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0));
344b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
345b0c7db22SLisandro Dalcin }
346b0c7db22SLisandro Dalcin 
347b0c7db22SLisandro Dalcin /*@C
348b0c7db22SLisandro Dalcin   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
349b0c7db22SLisandro Dalcin 
350b0c7db22SLisandro Dalcin   Collective on sf
351b0c7db22SLisandro Dalcin 
352b0c7db22SLisandro Dalcin   Input Parameters:
353cab54364SBarry Smith + sf          - The `PetscSF`
354b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
355b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
356b0c7db22SLisandro Dalcin 
357b0c7db22SLisandro Dalcin   Output Parameter:
358b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
359b0c7db22SLisandro Dalcin 
360b0c7db22SLisandro Dalcin   Level: developer
361b0c7db22SLisandro Dalcin 
362cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
363b0c7db22SLisandro Dalcin @*/
364d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
365d71ae5a4SJacob Faibussowitsch {
366b0c7db22SLisandro Dalcin   PetscSF         embedSF;
367b0c7db22SLisandro Dalcin   const PetscInt *indices;
368b0c7db22SLisandro Dalcin   IS              selected;
369b0c7db22SLisandro Dalcin   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
370b0c7db22SLisandro Dalcin 
371b0c7db22SLisandro Dalcin   PetscFunctionBegin;
372b0c7db22SLisandro Dalcin   *remoteOffsets = NULL;
3739566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL));
374b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
3759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0));
3769566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
3779566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
3789566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
3799566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected, &indices));
3809566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
3819566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(selected, &indices));
3829566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selected));
3839566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets));
3849566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
3859566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart], MPI_REPLACE));
3869566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&embedSF));
3879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0));
388b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
389b0c7db22SLisandro Dalcin }
390b0c7db22SLisandro Dalcin 
391b0c7db22SLisandro Dalcin /*@C
392cab54364SBarry Smith   PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points
393b0c7db22SLisandro Dalcin 
394b0c7db22SLisandro Dalcin   Collective on sf
395b0c7db22SLisandro Dalcin 
396b0c7db22SLisandro Dalcin   Input Parameters:
397cab54364SBarry Smith + sf - The `PetscSF`
398b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
399b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
400b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is the distributed section)
401b0c7db22SLisandro Dalcin 
402b0c7db22SLisandro Dalcin   Output Parameters:
403cab54364SBarry Smith - sectionSF - The new `PetscSF`
404b0c7db22SLisandro Dalcin 
405b0c7db22SLisandro Dalcin   Level: advanced
406b0c7db22SLisandro Dalcin 
407cab54364SBarry Smith   Note:
408cab54364SBarry Smith   Either rootSection or remoteOffsets can be specified
409cab54364SBarry Smith 
410cab54364SBarry Smith .seealso:  `PetscSF`, `PetscSFCreate()`
411b0c7db22SLisandro Dalcin @*/
412d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
413d71ae5a4SJacob Faibussowitsch {
414b0c7db22SLisandro Dalcin   MPI_Comm           comm;
415b0c7db22SLisandro Dalcin   const PetscInt    *localPoints;
416b0c7db22SLisandro Dalcin   const PetscSFNode *remotePoints;
417b0c7db22SLisandro Dalcin   PetscInt           lpStart, lpEnd;
418b0c7db22SLisandro Dalcin   PetscInt           numRoots, numSectionRoots, numPoints, numIndices = 0;
419b0c7db22SLisandro Dalcin   PetscInt          *localIndices;
420b0c7db22SLisandro Dalcin   PetscSFNode       *remoteIndices;
421b0c7db22SLisandro Dalcin   PetscInt           i, ind;
422b0c7db22SLisandro Dalcin 
423b0c7db22SLisandro Dalcin   PetscFunctionBegin;
424b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
425b0c7db22SLisandro Dalcin   PetscValidPointer(rootSection, 2);
426b0c7db22SLisandro Dalcin   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
427b0c7db22SLisandro Dalcin   PetscValidPointer(leafSection, 4);
428b0c7db22SLisandro Dalcin   PetscValidPointer(sectionSF, 5);
4299566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
4309566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, sectionSF));
4319566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
4329566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots));
4339566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints));
434b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
4359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0));
436b0c7db22SLisandro Dalcin   for (i = 0; i < numPoints; ++i) {
437b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
438b0c7db22SLisandro Dalcin     PetscInt dof;
439b0c7db22SLisandro Dalcin 
440b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
4419566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
44222a18c9fSMatthew G. Knepley       numIndices += dof < 0 ? 0 : dof;
443b0c7db22SLisandro Dalcin     }
444b0c7db22SLisandro Dalcin   }
4459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numIndices, &localIndices));
4469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numIndices, &remoteIndices));
447b0c7db22SLisandro Dalcin   /* Create new index graph */
448b0c7db22SLisandro Dalcin   for (i = 0, ind = 0; i < numPoints; ++i) {
449b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
450b0c7db22SLisandro Dalcin     PetscInt rank       = remotePoints[i].rank;
451b0c7db22SLisandro Dalcin 
452b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
453b0c7db22SLisandro Dalcin       PetscInt remoteOffset = remoteOffsets[localPoint - lpStart];
454b0c7db22SLisandro Dalcin       PetscInt loff, dof, d;
455b0c7db22SLisandro Dalcin 
4569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff));
4579566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
458b0c7db22SLisandro Dalcin       for (d = 0; d < dof; ++d, ++ind) {
459b0c7db22SLisandro Dalcin         localIndices[ind]        = loff + d;
460b0c7db22SLisandro Dalcin         remoteIndices[ind].rank  = rank;
461b0c7db22SLisandro Dalcin         remoteIndices[ind].index = remoteOffset + d;
462b0c7db22SLisandro Dalcin       }
463b0c7db22SLisandro Dalcin     }
464b0c7db22SLisandro Dalcin   }
46508401ef6SPierre Jolivet   PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
4669566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER));
4679566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(*sectionSF));
4689566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0));
469b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
470b0c7db22SLisandro Dalcin }
4718fa9e22eSVaclav Hapla 
4728fa9e22eSVaclav Hapla /*@C
473cab54364SBarry Smith    PetscSFCreateFromLayouts - Creates a parallel star forest mapping two `PetscLayout` objects
4748fa9e22eSVaclav Hapla 
4758fa9e22eSVaclav Hapla    Collective
4768fa9e22eSVaclav Hapla 
4774165533cSJose E. Roman    Input Parameters:
478cab54364SBarry Smith +  rmap - `PetscLayout` defining the global root space
479cab54364SBarry Smith -  lmap - `PetscLayout` defining the global leaf space
4808fa9e22eSVaclav Hapla 
4814165533cSJose E. Roman    Output Parameter:
4828fa9e22eSVaclav Hapla .  sf - The parallel star forest
4838fa9e22eSVaclav Hapla 
4848fa9e22eSVaclav Hapla    Level: intermediate
4858fa9e22eSVaclav Hapla 
486cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()`
4878fa9e22eSVaclav Hapla @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf)
489d71ae5a4SJacob Faibussowitsch {
4908fa9e22eSVaclav Hapla   PetscInt     i, nroots, nleaves = 0;
4918fa9e22eSVaclav Hapla   PetscInt     rN, lst, len;
4928fa9e22eSVaclav Hapla   PetscMPIInt  owner = -1;
4938fa9e22eSVaclav Hapla   PetscSFNode *remote;
4948fa9e22eSVaclav Hapla   MPI_Comm     rcomm = rmap->comm;
4958fa9e22eSVaclav Hapla   MPI_Comm     lcomm = lmap->comm;
4968fa9e22eSVaclav Hapla   PetscMPIInt  flg;
4978fa9e22eSVaclav Hapla 
4988fa9e22eSVaclav Hapla   PetscFunctionBegin;
4998fa9e22eSVaclav Hapla   PetscValidPointer(sf, 3);
50028b400f6SJacob Faibussowitsch   PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup");
50128b400f6SJacob Faibussowitsch   PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup");
5029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg));
503c9cc58a2SBarry Smith   PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators");
5049566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(rcomm, sf));
5059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(rmap, &nroots));
5069566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(rmap, &rN));
5079566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(lmap, &lst, &len));
5089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len - lst, &remote));
5098fa9e22eSVaclav Hapla   for (i = lst; i < len && i < rN; i++) {
51048a46eb9SPierre Jolivet     if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner));
5118fa9e22eSVaclav Hapla     remote[nleaves].rank  = owner;
5128fa9e22eSVaclav Hapla     remote[nleaves].index = i - rmap->range[owner];
5138fa9e22eSVaclav Hapla     nleaves++;
5148fa9e22eSVaclav Hapla   }
5159566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFree(remote));
5178fa9e22eSVaclav Hapla   PetscFunctionReturn(0);
5188fa9e22eSVaclav Hapla }
5198fa9e22eSVaclav Hapla 
5208fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
521d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt **oidxs, PetscInt **ogidxs)
522d71ae5a4SJacob Faibussowitsch {
5238fa9e22eSVaclav Hapla   PetscInt    *owners = map->range;
5248fa9e22eSVaclav Hapla   PetscInt     n      = map->n;
5258fa9e22eSVaclav Hapla   PetscSF      sf;
5268fa9e22eSVaclav Hapla   PetscInt    *lidxs, *work = NULL;
5278fa9e22eSVaclav Hapla   PetscSFNode *ridxs;
5288fa9e22eSVaclav Hapla   PetscMPIInt  rank, p = 0;
5298fa9e22eSVaclav Hapla   PetscInt     r, len  = 0;
5308fa9e22eSVaclav Hapla 
5318fa9e22eSVaclav Hapla   PetscFunctionBegin;
5328fa9e22eSVaclav Hapla   if (on) *on = 0; /* squelch -Wmaybe-uninitialized */
5338fa9e22eSVaclav Hapla   /* Create SF where leaves are input idxs and roots are owned idxs */
5349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
5359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &lidxs));
5368fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r) lidxs[r] = -1;
5379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N, &ridxs));
5388fa9e22eSVaclav Hapla   for (r = 0; r < N; ++r) {
5398fa9e22eSVaclav Hapla     const PetscInt idx = idxs[r];
540c9cc58a2SBarry 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);
5418fa9e22eSVaclav Hapla     if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */
5429566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(map, idx, &p));
5438fa9e22eSVaclav Hapla     }
5448fa9e22eSVaclav Hapla     ridxs[r].rank  = p;
5458fa9e22eSVaclav Hapla     ridxs[r].index = idxs[r] - owners[p];
5468fa9e22eSVaclav Hapla   }
5479566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(map->comm, &sf));
5489566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER));
5499566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
5509566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
5518fa9e22eSVaclav Hapla   if (ogidxs) { /* communicate global idxs */
5528fa9e22eSVaclav Hapla     PetscInt cum = 0, start, *work2;
5538fa9e22eSVaclav Hapla 
5549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &work));
5559566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(N, &work2));
5569371c9d4SSatish Balay     for (r = 0; r < N; ++r)
5579371c9d4SSatish Balay       if (idxs[r] >= 0) cum++;
5589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm));
5598fa9e22eSVaclav Hapla     start -= cum;
5608fa9e22eSVaclav Hapla     cum = 0;
5619371c9d4SSatish Balay     for (r = 0; r < N; ++r)
5629371c9d4SSatish Balay       if (idxs[r] >= 0) work2[r] = start + cum++;
5639566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE));
5649566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE));
5659566063dSJacob Faibussowitsch     PetscCall(PetscFree(work2));
5668fa9e22eSVaclav Hapla   }
5679566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5688fa9e22eSVaclav Hapla   /* Compress and put in indices */
5698fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r)
5708fa9e22eSVaclav Hapla     if (lidxs[r] >= 0) {
5718fa9e22eSVaclav Hapla       if (work) work[len] = work[r];
5728fa9e22eSVaclav Hapla       lidxs[len++] = r;
5738fa9e22eSVaclav Hapla     }
5748fa9e22eSVaclav Hapla   if (on) *on = len;
5758fa9e22eSVaclav Hapla   if (oidxs) *oidxs = lidxs;
5768fa9e22eSVaclav Hapla   if (ogidxs) *ogidxs = work;
5778fa9e22eSVaclav Hapla   PetscFunctionReturn(0);
5788fa9e22eSVaclav Hapla }
579deffd5ebSksagiyam 
580deffd5ebSksagiyam /*@
581cab54364SBarry Smith   PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices
582deffd5ebSksagiyam 
583deffd5ebSksagiyam   Collective
584deffd5ebSksagiyam 
585deffd5ebSksagiyam   Input Parameters:
586cab54364SBarry Smith + layout - `PetscLayout` defining the global index space and the rank that brokers each index
587deffd5ebSksagiyam . numRootIndices - size of rootIndices
588cab54364SBarry Smith . rootIndices - `PetscInt` array of global indices of which this process requests ownership
589deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation)
590deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices
591deffd5ebSksagiyam . numLeafIndices - size of leafIndices
592cab54364SBarry Smith . leafIndices - `PetscInt` array of global indices with which this process requires data associated
593deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation)
594deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices
595deffd5ebSksagiyam 
596d8d19677SJose E. Roman   Output Parameters:
597deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed)
598deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space
599deffd5ebSksagiyam 
600cab54364SBarry Smith   Level: advanced
601cab54364SBarry Smith 
602deffd5ebSksagiyam   Example 1:
603cab54364SBarry Smith .vb
604cab54364SBarry Smith   rank             : 0            1            2
605cab54364SBarry Smith   rootIndices      : [1 0 2]      [3]          [3]
606cab54364SBarry Smith   rootLocalOffset  : 100          200          300
607cab54364SBarry Smith   layout           : [0 1]        [2]          [3]
608cab54364SBarry Smith   leafIndices      : [0]          [2]          [0 3]
609cab54364SBarry Smith   leafLocalOffset  : 400          500          600
610cab54364SBarry Smith 
611cab54364SBarry Smith would build the following PetscSF
612cab54364SBarry Smith 
613cab54364SBarry Smith   [0] 400 <- (0,101)
614cab54364SBarry Smith   [1] 500 <- (0,102)
615cab54364SBarry Smith   [2] 600 <- (0,101)
616cab54364SBarry Smith   [2] 601 <- (2,300)
617cab54364SBarry Smith .ve
618cab54364SBarry Smith 
619deffd5ebSksagiyam   Example 2:
620cab54364SBarry Smith .vb
621cab54364SBarry Smith   rank             : 0               1               2
622cab54364SBarry Smith   rootIndices      : [1 0 2]         [3]             [3]
623cab54364SBarry Smith   rootLocalOffset  : 100             200             300
624cab54364SBarry Smith   layout           : [0 1]           [2]             [3]
625cab54364SBarry Smith   leafIndices      : rootIndices     rootIndices     rootIndices
626cab54364SBarry Smith   leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset
627cab54364SBarry Smith 
628cab54364SBarry Smith would build the following PetscSF
629cab54364SBarry Smith 
630cab54364SBarry Smith   [1] 200 <- (2,300)
631cab54364SBarry Smith .ve
632cab54364SBarry Smith 
633deffd5ebSksagiyam   Example 3:
634cab54364SBarry Smith .vb
635cab54364SBarry Smith   No process requests ownership of global index 1, but no process needs it.
636cab54364SBarry Smith 
637cab54364SBarry Smith   rank             : 0            1            2
638cab54364SBarry Smith   numRootIndices   : 2            1            1
639cab54364SBarry Smith   rootIndices      : [0 2]        [3]          [3]
640cab54364SBarry Smith   rootLocalOffset  : 100          200          300
641cab54364SBarry Smith   layout           : [0 1]        [2]          [3]
642cab54364SBarry Smith   numLeafIndices   : 1            1            2
643cab54364SBarry Smith   leafIndices      : [0]          [2]          [0 3]
644cab54364SBarry Smith   leafLocalOffset  : 400          500          600
645cab54364SBarry Smith 
646cab54364SBarry Smith would build the following PetscSF
647cab54364SBarry Smith 
648cab54364SBarry Smith   [0] 400 <- (0,100)
649cab54364SBarry Smith   [1] 500 <- (0,101)
650cab54364SBarry Smith   [2] 600 <- (0,100)
651cab54364SBarry Smith   [2] 601 <- (2,300)
652cab54364SBarry Smith .ve
653deffd5ebSksagiyam 
654deffd5ebSksagiyam   Notes:
655deffd5ebSksagiyam   The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its
656cab54364SBarry Smith   local size can be set to `PETSC_DECIDE`.
657cab54364SBarry Smith 
658deffd5ebSksagiyam   If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests
659deffd5ebSksagiyam   ownership of x and sends its own rank and the local index of x to process i.
660deffd5ebSksagiyam   If multiple processes request ownership of x, the one with the highest rank is to own x.
661deffd5ebSksagiyam   Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the
662deffd5ebSksagiyam   ownership information of x.
663deffd5ebSksagiyam   The output sf is constructed by associating each leaf point to a root point in this way.
664deffd5ebSksagiyam 
665deffd5ebSksagiyam   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
666cab54364SBarry Smith   The optional output `PetscSF` object sfA can be used to push such data to leaf points.
667deffd5ebSksagiyam 
668deffd5ebSksagiyam   All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices
669deffd5ebSksagiyam   must cover that of leafIndices, but need not cover the entire layout.
670deffd5ebSksagiyam 
671deffd5ebSksagiyam   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
672deffd5ebSksagiyam   star forest is almost identity, so will only include non-trivial part of the map.
673deffd5ebSksagiyam 
674cab54364SBarry Smith   Developer Note:
675deffd5ebSksagiyam   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
676deffd5ebSksagiyam   hash(rank, root_local_index) as the bid for the ownership determination.
677deffd5ebSksagiyam 
678cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
679deffd5ebSksagiyam @*/
680d71ae5a4SJacob 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)
681d71ae5a4SJacob Faibussowitsch {
682deffd5ebSksagiyam   MPI_Comm     comm = layout->comm;
683deffd5ebSksagiyam   PetscMPIInt  size, rank;
684deffd5ebSksagiyam   PetscSF      sf1;
685deffd5ebSksagiyam   PetscSFNode *owners, *buffer, *iremote;
686deffd5ebSksagiyam   PetscInt    *ilocal, nleaves, N, n, i;
687deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
688deffd5ebSksagiyam   PetscInt N1;
689deffd5ebSksagiyam #endif
690deffd5ebSksagiyam   PetscBool flag;
691deffd5ebSksagiyam 
692deffd5ebSksagiyam   PetscFunctionBegin;
693deffd5ebSksagiyam   if (rootIndices) PetscValidIntPointer(rootIndices, 3);
694deffd5ebSksagiyam   if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices, 4);
695deffd5ebSksagiyam   if (leafIndices) PetscValidIntPointer(leafIndices, 7);
696deffd5ebSksagiyam   if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices, 8);
697deffd5ebSksagiyam   if (sfA) PetscValidPointer(sfA, 10);
698deffd5ebSksagiyam   PetscValidPointer(sf, 11);
69908401ef6SPierre Jolivet   PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
70008401ef6SPierre Jolivet   PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
7019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
7049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(layout, &N));
7059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &n));
706deffd5ebSksagiyam   flag = (PetscBool)(leafIndices == rootIndices);
7071c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm));
708c9cc58a2SBarry Smith   PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
709deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
710deffd5ebSksagiyam   N1 = PETSC_MIN_INT;
7119371c9d4SSatish Balay   for (i = 0; i < numRootIndices; i++)
7129371c9d4SSatish Balay     if (rootIndices[i] > N1) N1 = rootIndices[i];
7139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
71408401ef6SPierre 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);
715deffd5ebSksagiyam   if (!flag) {
716deffd5ebSksagiyam     N1 = PETSC_MIN_INT;
7179371c9d4SSatish Balay     for (i = 0; i < numLeafIndices; i++)
7189371c9d4SSatish Balay       if (leafIndices[i] > N1) N1 = leafIndices[i];
7199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
72008401ef6SPierre 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);
721deffd5ebSksagiyam   }
722deffd5ebSksagiyam #endif
723deffd5ebSksagiyam   /* Reduce: owners -> buffer */
7249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &buffer));
7259566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, &sf1));
7269566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf1));
7279566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices));
7289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootIndices, &owners));
729deffd5ebSksagiyam   for (i = 0; i < numRootIndices; ++i) {
730deffd5ebSksagiyam     owners[i].rank  = rank;
731deffd5ebSksagiyam     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
732deffd5ebSksagiyam   }
733deffd5ebSksagiyam   for (i = 0; i < n; ++i) {
734deffd5ebSksagiyam     buffer[i].index = -1;
735deffd5ebSksagiyam     buffer[i].rank  = -1;
736deffd5ebSksagiyam   }
7379566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
7389566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
739deffd5ebSksagiyam   /* Bcast: buffer -> owners */
740deffd5ebSksagiyam   if (!flag) {
741deffd5ebSksagiyam     /* leafIndices is different from rootIndices */
7429566063dSJacob Faibussowitsch     PetscCall(PetscFree(owners));
7439566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices));
7449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeafIndices, &owners));
745deffd5ebSksagiyam   }
7469566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
7479566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
74808401ef6SPierre 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]);
7499566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
7509371c9d4SSatish Balay   if (sfA) {
7519371c9d4SSatish Balay     *sfA = sf1;
7529371c9d4SSatish Balay   } else PetscCall(PetscSFDestroy(&sf1));
753deffd5ebSksagiyam   /* Create sf */
754deffd5ebSksagiyam   if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
755deffd5ebSksagiyam     /* leaf space == root space */
7569371c9d4SSatish Balay     for (i = 0, nleaves = 0; i < numLeafIndices; ++i)
7579371c9d4SSatish Balay       if (owners[i].rank != rank) ++nleaves;
7589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &ilocal));
7599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &iremote));
760deffd5ebSksagiyam     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
761deffd5ebSksagiyam       if (owners[i].rank != rank) {
762deffd5ebSksagiyam         ilocal[nleaves]        = leafLocalOffset + i;
763deffd5ebSksagiyam         iremote[nleaves].rank  = owners[i].rank;
764deffd5ebSksagiyam         iremote[nleaves].index = owners[i].index;
765deffd5ebSksagiyam         ++nleaves;
766deffd5ebSksagiyam       }
767deffd5ebSksagiyam     }
7689566063dSJacob Faibussowitsch     PetscCall(PetscFree(owners));
769deffd5ebSksagiyam   } else {
770deffd5ebSksagiyam     nleaves = numLeafIndices;
7719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &ilocal));
772ad540459SPierre Jolivet     for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);
773deffd5ebSksagiyam     iremote = owners;
774deffd5ebSksagiyam   }
7759566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, sf));
7769566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sf));
7779566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
778deffd5ebSksagiyam   PetscFunctionReturn(0);
779deffd5ebSksagiyam }
780