xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision 2f04c52277ae86d0bd99bd90d9d5574dfa2d51e6)
1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h> /*I  "petscsf.h"   I*/
2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h>
3b0c7db22SLisandro Dalcin 
4eb58282aSVaclav Hapla /*@
5*2f04c522SBarry Smith   PetscSFSetGraphLayout - Set a `PetscSF` communication pattern using global indices and a `PetscLayout`
6b0c7db22SLisandro Dalcin 
7b0c7db22SLisandro Dalcin   Collective
8b0c7db22SLisandro Dalcin 
94165533cSJose E. Roman   Input Parameters:
10b0c7db22SLisandro Dalcin + sf        - star forest
11*2f04c522SBarry Smith . layout    - `PetscLayout` defining the global space for roots, i.e. which roots are owned by each MPI process
12*2f04c522SBarry Smith . nleaves   - number of leaf vertices on the current process, each of these references a root on any MPI process
13*2f04c522SBarry Smith . ilocal    - locations of leaves in leafdata buffers, pass `NULL` for contiguous storage, that is the locations are in [0,`nleaves`)
14*2f04c522SBarry Smith . localmode - copy mode for `ilocal`
15*2f04c522SBarry Smith - gremote   - root vertices in global numbering corresponding to the leaves
16b0c7db22SLisandro Dalcin 
17b0c7db22SLisandro Dalcin   Level: intermediate
18b0c7db22SLisandro Dalcin 
19cab54364SBarry Smith   Note:
20*2f04c522SBarry Smith   Global indices must lie in [0, N) where N is the global size of `layout`.
21*2f04c522SBarry Smith   Leaf indices in `ilocal` get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`.
2286c56f52SVaclav Hapla 
2338b5cf2dSJacob Faibussowitsch   Developer Notes:
24*2f04c522SBarry Smith   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 
28*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
29b0c7db22SLisandro Dalcin @*/
30cf84bf9eSBarry Smith 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;
38da0802e2SStefano Zampini   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
39b114984aSVaclav Hapla   PetscCall(PetscLayoutSetUp(layout));
409566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &nroots));
419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRanges(layout, &range));
429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &remote));
43eb58282aSVaclav Hapla   if (nleaves) ls = gremote[0] + 1;
44b0c7db22SLisandro Dalcin   for (i = 0; i < nleaves; i++) {
45eb58282aSVaclav Hapla     const PetscInt idx = gremote[i] - ls;
4638a25198SStefano Zampini     if (idx < 0 || idx >= ln) { /* short-circuit the search */
47eb58282aSVaclav Hapla       PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index));
4838a25198SStefano Zampini       remote[i].rank = lr;
4938a25198SStefano Zampini       ls             = range[lr];
5038a25198SStefano Zampini       ln             = range[lr + 1] - ls;
5138a25198SStefano Zampini     } else {
5238a25198SStefano Zampini       remote[i].rank  = lr;
5338a25198SStefano Zampini       remote[i].index = idx;
5438a25198SStefano Zampini     }
55b0c7db22SLisandro Dalcin   }
569566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58b0c7db22SLisandro Dalcin }
59b0c7db22SLisandro Dalcin 
60eb58282aSVaclav Hapla /*@C
61*2f04c522SBarry Smith   PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe a `PetscSF`
62eb58282aSVaclav Hapla 
63eb58282aSVaclav Hapla   Collective
64eb58282aSVaclav Hapla 
65eb58282aSVaclav Hapla   Input Parameter:
66eb58282aSVaclav Hapla . sf - star forest
67eb58282aSVaclav Hapla 
68eb58282aSVaclav Hapla   Output Parameters:
69cab54364SBarry Smith + layout  - `PetscLayout` defining the global space for roots
70eb58282aSVaclav Hapla . nleaves - number of leaf vertices on the current process, each of these references a root on any process
71f13dfd9eSBarry Smith . ilocal  - locations of leaves in leafdata buffers, or `NULL` for contiguous storage
72*2f04c522SBarry Smith - gremote - root vertices in global numbering corresponding to the leaves
73eb58282aSVaclav Hapla 
74eb58282aSVaclav Hapla   Level: intermediate
75eb58282aSVaclav Hapla 
76eb58282aSVaclav Hapla   Notes:
77eb58282aSVaclav Hapla   The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest.
78f13dfd9eSBarry Smith   The outputs `layout` and `gremote` are freshly created each time this function is called,
79f13dfd9eSBarry Smith   so they need to be freed (with `PetscLayoutDestroy()` and `PetscFree()`) by the user.
80eb58282aSVaclav Hapla 
81*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
82eb58282aSVaclav Hapla @*/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[])
84d71ae5a4SJacob Faibussowitsch {
85eb58282aSVaclav Hapla   PetscInt           nr, nl;
86eb58282aSVaclav Hapla   const PetscSFNode *ir;
87eb58282aSVaclav Hapla   PetscLayout        lt;
88eb58282aSVaclav Hapla 
89eb58282aSVaclav Hapla   PetscFunctionBegin;
90eb58282aSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir));
91eb58282aSVaclav Hapla   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, &lt));
92eb58282aSVaclav Hapla   if (gremote) {
93eb58282aSVaclav Hapla     PetscInt        i;
94eb58282aSVaclav Hapla     const PetscInt *range;
95eb58282aSVaclav Hapla     PetscInt       *gr;
96eb58282aSVaclav Hapla 
97eb58282aSVaclav Hapla     PetscCall(PetscLayoutGetRanges(lt, &range));
98eb58282aSVaclav Hapla     PetscCall(PetscMalloc1(nl, &gr));
99eb58282aSVaclav Hapla     for (i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index;
100eb58282aSVaclav Hapla     *gremote = gr;
101eb58282aSVaclav Hapla   }
102eb58282aSVaclav Hapla   if (nleaves) *nleaves = nl;
103eb58282aSVaclav Hapla   if (layout) *layout = lt;
104eb58282aSVaclav Hapla   else PetscCall(PetscLayoutDestroy(&lt));
1053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106eb58282aSVaclav Hapla }
107eb58282aSVaclav Hapla 
108b0c7db22SLisandro Dalcin /*@
109*2f04c522SBarry Smith   PetscSFSetGraphSection - Sets the `PetscSF` graph (communication pattern) encoding the parallel dof overlap based upon the `PetscSection` describing the data layout.
110b0c7db22SLisandro Dalcin 
111b0c7db22SLisandro Dalcin   Input Parameters:
112cab54364SBarry Smith + sf            - The `PetscSF`
113cab54364SBarry Smith . localSection  - `PetscSection` describing the local data layout
114cab54364SBarry Smith - globalSection - `PetscSection` describing the global data layout
115b0c7db22SLisandro Dalcin 
116b0c7db22SLisandro Dalcin   Level: developer
117b0c7db22SLisandro Dalcin 
118*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()`
119b0c7db22SLisandro Dalcin @*/
120d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
121d71ae5a4SJacob Faibussowitsch {
122b0c7db22SLisandro Dalcin   MPI_Comm        comm;
123b0c7db22SLisandro Dalcin   PetscLayout     layout;
124b0c7db22SLisandro Dalcin   const PetscInt *ranges;
125b0c7db22SLisandro Dalcin   PetscInt       *local;
126b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
127b0c7db22SLisandro Dalcin   PetscInt        pStart, pEnd, p, nroots, nleaves = 0, l;
128b0c7db22SLisandro Dalcin   PetscMPIInt     size, rank;
129b0c7db22SLisandro Dalcin 
130b0c7db22SLisandro Dalcin   PetscFunctionBegin;
131b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
132b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
133b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
134b0c7db22SLisandro Dalcin 
1359566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
1369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1389566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
1399566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
1409566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &layout));
1419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(layout, 1));
1429566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
1439566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
1449566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRanges(layout, &ranges));
145b0c7db22SLisandro Dalcin   for (p = pStart; p < pEnd; ++p) {
146b0c7db22SLisandro Dalcin     PetscInt gdof, gcdof;
147b0c7db22SLisandro Dalcin 
1489566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
1499566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
15057508eceSPierre Jolivet     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);
151b0c7db22SLisandro Dalcin     nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
152b0c7db22SLisandro Dalcin   }
1539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &local));
1549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &remote));
155b0c7db22SLisandro Dalcin   for (p = pStart, l = 0; p < pEnd; ++p) {
156b0c7db22SLisandro Dalcin     const PetscInt *cind;
157b0c7db22SLisandro Dalcin     PetscInt        dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
158b0c7db22SLisandro Dalcin 
1599566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(localSection, p, &dof));
1609566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(localSection, p, &off));
1619566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
1629566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind));
1639566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
1649566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
1659566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
166b0c7db22SLisandro Dalcin     if (!gdof) continue; /* Censored point */
167b0c7db22SLisandro Dalcin     gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
168b0c7db22SLisandro Dalcin     if (gsize != dof - cdof) {
16908401ef6SPierre 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);
170b0c7db22SLisandro Dalcin       cdof = 0; /* Ignore constraints */
171b0c7db22SLisandro Dalcin     }
172b0c7db22SLisandro Dalcin     for (d = 0, c = 0; d < dof; ++d) {
1739371c9d4SSatish Balay       if ((c < cdof) && (cind[c] == d)) {
1749371c9d4SSatish Balay         ++c;
1759371c9d4SSatish Balay         continue;
1769371c9d4SSatish Balay       }
177b0c7db22SLisandro Dalcin       local[l + d - c] = off + d;
178b0c7db22SLisandro Dalcin     }
17908401ef6SPierre 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);
180b0c7db22SLisandro Dalcin     if (gdof < 0) {
181b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
1826497c311SBarry Smith         PetscInt    offset = -(goff + 1) + d, ir;
1836497c311SBarry Smith         PetscMPIInt r;
184b0c7db22SLisandro Dalcin 
1856497c311SBarry Smith         PetscCall(PetscFindInt(offset, size + 1, ranges, &ir));
1866497c311SBarry Smith         PetscCall(PetscMPIIntCast(ir, &r));
187b0c7db22SLisandro Dalcin         if (r < 0) r = -(r + 2);
1886497c311SBarry Smith         PetscCheck(!(r < 0) && !(r >= size), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %d (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff);
189b0c7db22SLisandro Dalcin         remote[l].rank  = r;
190b0c7db22SLisandro Dalcin         remote[l].index = offset - ranges[r];
191b0c7db22SLisandro Dalcin       }
192b0c7db22SLisandro Dalcin     } else {
193b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
194b0c7db22SLisandro Dalcin         remote[l].rank  = rank;
195b0c7db22SLisandro Dalcin         remote[l].index = goff + d - ranges[rank];
196b0c7db22SLisandro Dalcin       }
197b0c7db22SLisandro Dalcin     }
198b0c7db22SLisandro Dalcin   }
19908401ef6SPierre Jolivet   PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
2009566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
2019566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
203b0c7db22SLisandro Dalcin }
204b0c7db22SLisandro Dalcin 
205b0c7db22SLisandro Dalcin /*@C
206cab54364SBarry Smith   PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF`
207b0c7db22SLisandro Dalcin 
208c3339decSBarry Smith   Collective
209b0c7db22SLisandro Dalcin 
210b0c7db22SLisandro Dalcin   Input Parameters:
211cab54364SBarry Smith + sf          - The `PetscSF`
212b0c7db22SLisandro Dalcin - rootSection - Section defined on root space
213b0c7db22SLisandro Dalcin 
214b0c7db22SLisandro Dalcin   Output Parameters:
215ce78bad3SBarry Smith + remoteOffsets - root offsets in leaf storage, or `NULL`, its length will be the size of the chart of `leafSection`
216b0c7db22SLisandro Dalcin - leafSection   - Section defined on the leaf space
217b0c7db22SLisandro Dalcin 
218b0c7db22SLisandro Dalcin   Level: advanced
219b0c7db22SLisandro Dalcin 
220ce78bad3SBarry Smith   Note:
221ce78bad3SBarry Smith   Caller must `PetscFree()` `remoteOffsets` if it was requested
222ce78bad3SBarry Smith 
223ce78bad3SBarry Smith   Fortran Note:
224ce78bad3SBarry Smith   Use `PetscSFDestroyRemoteOffsets()` when `remoteOffsets` is no longer needed.
22523e9620eSJunchao Zhang 
226*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`
227b0c7db22SLisandro Dalcin @*/
228ce78bad3SBarry Smith PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt *remoteOffsets[], PetscSection leafSection)
229d71ae5a4SJacob Faibussowitsch {
230b0c7db22SLisandro Dalcin   PetscSF         embedSF;
231b0c7db22SLisandro Dalcin   const PetscInt *indices;
232b0c7db22SLisandro Dalcin   IS              selected;
2331690c2aeSBarry Smith   PetscInt        numFields, nroots, rpStart, rpEnd, lpStart = PETSC_INT_MAX, lpEnd = -1, f, c;
234b0c7db22SLisandro Dalcin   PetscBool      *sub, hasc;
235b0c7db22SLisandro Dalcin 
236b0c7db22SLisandro Dalcin   PetscFunctionBegin;
2379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0));
2389566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(rootSection, &numFields));
239029e8818Sksagiyam   if (numFields) {
240029e8818Sksagiyam     IS perm;
241029e8818Sksagiyam 
242029e8818Sksagiyam     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
243029e8818Sksagiyam        leafSection->perm. To keep this permutation set by the user, we grab
244029e8818Sksagiyam        the reference before calling PetscSectionSetNumFields() and set it
245029e8818Sksagiyam        back after. */
2469566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetPermutation(leafSection, &perm));
2479566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)perm));
2489566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetNumFields(leafSection, numFields));
2499566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetPermutation(leafSection, perm));
2509566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&perm));
251029e8818Sksagiyam   }
2529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numFields + 2, &sub));
253b0c7db22SLisandro Dalcin   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
254b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
2552ee82556SMatthew G. Knepley     PetscSectionSym sym, dsym = NULL;
256b0c7db22SLisandro Dalcin     const char     *name    = NULL;
257b0c7db22SLisandro Dalcin     PetscInt        numComp = 0;
258b0c7db22SLisandro Dalcin 
259b0c7db22SLisandro Dalcin     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
2609566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp));
2619566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldName(rootSection, f, &name));
2629566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym));
2639566063dSJacob Faibussowitsch     if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym));
2649566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp));
2659566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldName(leafSection, f, name));
2669566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym));
2679566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymDestroy(&dsym));
268b0c7db22SLisandro Dalcin     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
2699566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name));
2709566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetComponentName(leafSection, f, c, name));
271b0c7db22SLisandro Dalcin     }
272b0c7db22SLisandro Dalcin   }
2739566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
2749566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
275b0c7db22SLisandro Dalcin   rpEnd = PetscMin(rpEnd, nroots);
276b0c7db22SLisandro Dalcin   rpEnd = PetscMax(rpStart, rpEnd);
277b0c7db22SLisandro Dalcin   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
278b0c7db22SLisandro Dalcin   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
2795440e5dcSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf)));
280b0c7db22SLisandro Dalcin   if (sub[0]) {
2819566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
2829566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(selected, &indices));
2839566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
2849566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(selected, &indices));
2859566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&selected));
286b0c7db22SLisandro Dalcin   } else {
2879566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)sf));
288b0c7db22SLisandro Dalcin     embedSF = sf;
289b0c7db22SLisandro Dalcin   }
2909566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd));
291b0c7db22SLisandro Dalcin   lpEnd++;
292b0c7db22SLisandro Dalcin 
2939566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd));
294b0c7db22SLisandro Dalcin 
295b0c7db22SLisandro Dalcin   /* Constrained dof section */
296b0c7db22SLisandro Dalcin   hasc = sub[1];
297b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]);
298b0c7db22SLisandro Dalcin 
299b0c7db22SLisandro Dalcin   /* Could fuse these at the cost of copies and extra allocation */
30016cd844bSPierre Jolivet   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE));
30116cd844bSPierre Jolivet   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE));
302b0c7db22SLisandro Dalcin   if (sub[1]) {
3037c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(rootSection));
3047c0883d5SVaclav Hapla     PetscCall(PetscSectionCheckConstraints_Private(leafSection));
3059566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
3069566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
307b0c7db22SLisandro Dalcin   }
308b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
30916cd844bSPierre Jolivet     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE));
31016cd844bSPierre Jolivet     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE));
311b0c7db22SLisandro Dalcin     if (sub[2 + f]) {
3127c0883d5SVaclav Hapla       PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f]));
3137c0883d5SVaclav Hapla       PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f]));
3149566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
3159566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
316b0c7db22SLisandro Dalcin     }
317b0c7db22SLisandro Dalcin   }
318b0c7db22SLisandro Dalcin   if (remoteOffsets) {
3199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets));
32016cd844bSPierre Jolivet     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
32116cd844bSPierre Jolivet     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
322b0c7db22SLisandro Dalcin   }
32369c11d05SVaclav Hapla   PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection));
3249566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSection));
325b0c7db22SLisandro Dalcin   if (hasc) { /* need to communicate bcIndices */
326b0c7db22SLisandro Dalcin     PetscSF   bcSF;
327b0c7db22SLisandro Dalcin     PetscInt *rOffBc;
328b0c7db22SLisandro Dalcin 
3299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc));
330b0c7db22SLisandro Dalcin     if (sub[1]) {
3319566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3329566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3339566063dSJacob Faibussowitsch       PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF));
3349566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
3359566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
3369566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&bcSF));
337b0c7db22SLisandro Dalcin     }
338b0c7db22SLisandro Dalcin     for (f = 0; f < numFields; ++f) {
339b0c7db22SLisandro Dalcin       if (sub[2 + f]) {
3409566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3419566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
3429566063dSJacob Faibussowitsch         PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF));
3439566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
3449566063dSJacob Faibussowitsch         PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
3459566063dSJacob Faibussowitsch         PetscCall(PetscSFDestroy(&bcSF));
346b0c7db22SLisandro Dalcin       }
347b0c7db22SLisandro Dalcin     }
3489566063dSJacob Faibussowitsch     PetscCall(PetscFree(rOffBc));
349b0c7db22SLisandro Dalcin   }
3509566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&embedSF));
3519566063dSJacob Faibussowitsch   PetscCall(PetscFree(sub));
3529566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354b0c7db22SLisandro Dalcin }
355b0c7db22SLisandro Dalcin 
356b0c7db22SLisandro Dalcin /*@C
357b0c7db22SLisandro Dalcin   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
358b0c7db22SLisandro Dalcin 
359c3339decSBarry Smith   Collective
360b0c7db22SLisandro Dalcin 
361b0c7db22SLisandro Dalcin   Input Parameters:
362cab54364SBarry Smith + sf          - The `PetscSF`
363ce78bad3SBarry Smith . rootSection - Data layout of remote points for outgoing data (this is layout for roots)
364ce78bad3SBarry Smith - leafSection - Data layout of local points for incoming data  (this is layout for leaves)
365b0c7db22SLisandro Dalcin 
366b0c7db22SLisandro Dalcin   Output Parameter:
367ce78bad3SBarry Smith . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or `NULL`
368b0c7db22SLisandro Dalcin 
369b0c7db22SLisandro Dalcin   Level: developer
370b0c7db22SLisandro Dalcin 
371ce78bad3SBarry Smith   Note:
372ce78bad3SBarry Smith   Caller must `PetscFree()` `remoteOffsets` if it was requested
373ce78bad3SBarry Smith 
374ce78bad3SBarry Smith   Fortran Note:
375ce78bad3SBarry Smith   Use `PetscSFDestroyRemoteOffsets()` when `remoteOffsets` is no longer needed.
37623e9620eSJunchao Zhang 
377*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`
378b0c7db22SLisandro Dalcin @*/
379ce78bad3SBarry Smith PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt *remoteOffsets[])
380d71ae5a4SJacob Faibussowitsch {
381b0c7db22SLisandro Dalcin   PetscSF         embedSF;
382b0c7db22SLisandro Dalcin   const PetscInt *indices;
383b0c7db22SLisandro Dalcin   IS              selected;
384b0c7db22SLisandro Dalcin   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
385b0c7db22SLisandro Dalcin 
386b0c7db22SLisandro Dalcin   PetscFunctionBegin;
387b0c7db22SLisandro Dalcin   *remoteOffsets = NULL;
3889566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL));
3893ba16761SJacob Faibussowitsch   if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS);
3909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0));
3919566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
3929566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
3939566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
3949566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(selected, &indices));
3959566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
3969566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(selected, &indices));
3979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&selected));
3989566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets));
3998e3a54c0SPierre Jolivet   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
4008e3a54c0SPierre Jolivet   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
4019566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&embedSF));
4029566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0));
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
404b0c7db22SLisandro Dalcin }
405b0c7db22SLisandro Dalcin 
406ce78bad3SBarry Smith /*@
407cab54364SBarry Smith   PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points
408b0c7db22SLisandro Dalcin 
409c3339decSBarry Smith   Collective
410b0c7db22SLisandro Dalcin 
411b0c7db22SLisandro Dalcin   Input Parameters:
412cab54364SBarry Smith + sf            - The `PetscSF`
413b0c7db22SLisandro Dalcin . rootSection   - Data layout of remote points for outgoing data (this is usually the serial section)
414*2f04c522SBarry Smith . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or `NULL`
415b0c7db22SLisandro Dalcin - leafSection   - Data layout of local points for incoming data  (this is the distributed section)
416b0c7db22SLisandro Dalcin 
4172fe279fdSBarry Smith   Output Parameter:
4182fe279fdSBarry Smith . sectionSF - The new `PetscSF`
419b0c7db22SLisandro Dalcin 
420b0c7db22SLisandro Dalcin   Level: advanced
421b0c7db22SLisandro Dalcin 
42223e9620eSJunchao Zhang   Notes:
423ce78bad3SBarry Smith   `remoteOffsets` can be `NULL` if `sf` does not reference any points in leafSection
42423e9620eSJunchao Zhang 
425*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`
426b0c7db22SLisandro Dalcin @*/
427d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
428d71ae5a4SJacob Faibussowitsch {
429b0c7db22SLisandro Dalcin   MPI_Comm           comm;
430b0c7db22SLisandro Dalcin   const PetscInt    *localPoints;
431b0c7db22SLisandro Dalcin   const PetscSFNode *remotePoints;
432b0c7db22SLisandro Dalcin   PetscInt           lpStart, lpEnd;
433b0c7db22SLisandro Dalcin   PetscInt           numRoots, numSectionRoots, numPoints, numIndices = 0;
434b0c7db22SLisandro Dalcin   PetscInt          *localIndices;
435b0c7db22SLisandro Dalcin   PetscSFNode       *remoteIndices;
436b0c7db22SLisandro Dalcin   PetscInt           i, ind;
437b0c7db22SLisandro Dalcin 
438b0c7db22SLisandro Dalcin   PetscFunctionBegin;
439b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
4404f572ea9SToby Isaac   PetscAssertPointer(rootSection, 2);
4414f572ea9SToby Isaac   /* Cannot check PetscAssertPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
4424f572ea9SToby Isaac   PetscAssertPointer(leafSection, 4);
4434f572ea9SToby Isaac   PetscAssertPointer(sectionSF, 5);
4449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
4459566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, sectionSF));
4469566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
4479566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots));
4489566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints));
4493ba16761SJacob Faibussowitsch   if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS);
4509566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0));
451b0c7db22SLisandro Dalcin   for (i = 0; i < numPoints; ++i) {
452b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
453b0c7db22SLisandro Dalcin     PetscInt dof;
454b0c7db22SLisandro Dalcin 
455b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
4569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
45722a18c9fSMatthew G. Knepley       numIndices += dof < 0 ? 0 : dof;
458b0c7db22SLisandro Dalcin     }
459b0c7db22SLisandro Dalcin   }
4609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numIndices, &localIndices));
4619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numIndices, &remoteIndices));
462b0c7db22SLisandro Dalcin   /* Create new index graph */
463b0c7db22SLisandro Dalcin   for (i = 0, ind = 0; i < numPoints; ++i) {
464b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
465835f2295SStefano Zampini     PetscInt rank       = remotePoints[i].rank;
466b0c7db22SLisandro Dalcin 
467b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
468b0c7db22SLisandro Dalcin       PetscInt remoteOffset = remoteOffsets[localPoint - lpStart];
469b0c7db22SLisandro Dalcin       PetscInt loff, dof, d;
470b0c7db22SLisandro Dalcin 
4719566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff));
4729566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
473b0c7db22SLisandro Dalcin       for (d = 0; d < dof; ++d, ++ind) {
474b0c7db22SLisandro Dalcin         localIndices[ind]        = loff + d;
475b0c7db22SLisandro Dalcin         remoteIndices[ind].rank  = rank;
476b0c7db22SLisandro Dalcin         remoteIndices[ind].index = remoteOffset + d;
477b0c7db22SLisandro Dalcin       }
478b0c7db22SLisandro Dalcin     }
479b0c7db22SLisandro Dalcin   }
48008401ef6SPierre Jolivet   PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
4819566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER));
4829566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(*sectionSF));
4839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0));
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
485b0c7db22SLisandro Dalcin }
4868fa9e22eSVaclav Hapla 
4878fa9e22eSVaclav Hapla /*@C
488*2f04c522SBarry Smith   PetscSFCreateFromLayouts - Creates a parallel star forest mapping between two `PetscLayout` objects
4898fa9e22eSVaclav Hapla 
4908fa9e22eSVaclav Hapla   Collective
4918fa9e22eSVaclav Hapla 
4924165533cSJose E. Roman   Input Parameters:
493cab54364SBarry Smith + rmap - `PetscLayout` defining the global root space
494cab54364SBarry Smith - lmap - `PetscLayout` defining the global leaf space
4958fa9e22eSVaclav Hapla 
4964165533cSJose E. Roman   Output Parameter:
4978fa9e22eSVaclav Hapla . sf - The parallel star forest
4988fa9e22eSVaclav Hapla 
4998fa9e22eSVaclav Hapla   Level: intermediate
5008fa9e22eSVaclav Hapla 
501*2f04c522SBarry Smith   Notes:
502*2f04c522SBarry Smith   If the global length of `lmap` differs from the global length of `rmap` then the excess entries are ignored.
503*2f04c522SBarry Smith 
504*2f04c522SBarry Smith   The resulting `sf` used with `PetscSFBcastBegin()` and `PetscSFBcastEnd()` merely copies the array entries of `rootdata` to
505*2f04c522SBarry Smith   `leafdata`; moving them between MPI processes if needed. For example,
506*2f04c522SBarry Smith   if rmap is [0, 3, 5) and lmap is [0, 2, 6) and `rootdata` is (1, 2, 3) on MPI rank 0 and (4, 5) on MPI rank 1 then the
507*2f04c522SBarry Smith   `leafdata` would become (1, 2) on MPI rank 0 and (3, 4, 5, x) on MPI rank 1.
508*2f04c522SBarry Smith 
509*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscLayout`, `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()`
5108fa9e22eSVaclav Hapla @*/
511d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf)
512d71ae5a4SJacob Faibussowitsch {
5138fa9e22eSVaclav Hapla   PetscInt     i, nroots, nleaves = 0;
5148fa9e22eSVaclav Hapla   PetscInt     rN, lst, len;
5158fa9e22eSVaclav Hapla   PetscMPIInt  owner = -1;
5168fa9e22eSVaclav Hapla   PetscSFNode *remote;
5178fa9e22eSVaclav Hapla   MPI_Comm     rcomm = rmap->comm;
5188fa9e22eSVaclav Hapla   MPI_Comm     lcomm = lmap->comm;
5198fa9e22eSVaclav Hapla   PetscMPIInt  flg;
5208fa9e22eSVaclav Hapla 
5218fa9e22eSVaclav Hapla   PetscFunctionBegin;
5224f572ea9SToby Isaac   PetscAssertPointer(sf, 3);
52328b400f6SJacob Faibussowitsch   PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup");
52428b400f6SJacob Faibussowitsch   PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup");
5259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg));
526c9cc58a2SBarry Smith   PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators");
5279566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(rcomm, sf));
5289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(rmap, &nroots));
5299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(rmap, &rN));
5309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(lmap, &lst, &len));
5319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len - lst, &remote));
5328fa9e22eSVaclav Hapla   for (i = lst; i < len && i < rN; i++) {
53348a46eb9SPierre Jolivet     if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner));
5348fa9e22eSVaclav Hapla     remote[nleaves].rank  = owner;
5358fa9e22eSVaclav Hapla     remote[nleaves].index = i - rmap->range[owner];
5368fa9e22eSVaclav Hapla     nleaves++;
5378fa9e22eSVaclav Hapla   }
5389566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES));
5399566063dSJacob Faibussowitsch   PetscCall(PetscFree(remote));
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5418fa9e22eSVaclav Hapla }
5428fa9e22eSVaclav Hapla 
5438fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
544ce78bad3SBarry Smith PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt *oidxs[], PetscInt *ogidxs[])
545d71ae5a4SJacob Faibussowitsch {
5468fa9e22eSVaclav Hapla   PetscInt    *owners = map->range;
5478fa9e22eSVaclav Hapla   PetscInt     n      = map->n;
5488fa9e22eSVaclav Hapla   PetscSF      sf;
549d61846d3SStefano Zampini   PetscInt    *lidxs, *work = NULL, *ilocal;
5508fa9e22eSVaclav Hapla   PetscSFNode *ridxs;
5518fa9e22eSVaclav Hapla   PetscMPIInt  rank, p = 0;
552d61846d3SStefano Zampini   PetscInt     r, len = 0, nleaves = 0;
5538fa9e22eSVaclav Hapla 
5548fa9e22eSVaclav Hapla   PetscFunctionBegin;
5558fa9e22eSVaclav Hapla   if (on) *on = 0; /* squelch -Wmaybe-uninitialized */
5568fa9e22eSVaclav Hapla   /* Create SF where leaves are input idxs and roots are owned idxs */
5579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
5589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &lidxs));
5598fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r) lidxs[r] = -1;
5609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N, &ridxs));
561d61846d3SStefano Zampini   PetscCall(PetscMalloc1(N, &ilocal));
5628fa9e22eSVaclav Hapla   for (r = 0; r < N; ++r) {
5638fa9e22eSVaclav Hapla     const PetscInt idx = idxs[r];
564d61846d3SStefano Zampini 
565d61846d3SStefano Zampini     if (idx < 0) continue;
566d61846d3SStefano Zampini     PetscCheck(idx < map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")", idx, map->N);
5678fa9e22eSVaclav Hapla     if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */
5689566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(map, idx, &p));
5698fa9e22eSVaclav Hapla     }
570d61846d3SStefano Zampini     ridxs[nleaves].rank  = p;
571d61846d3SStefano Zampini     ridxs[nleaves].index = idxs[r] - owners[p];
572d61846d3SStefano Zampini     ilocal[nleaves]      = r;
573d61846d3SStefano Zampini     nleaves++;
5748fa9e22eSVaclav Hapla   }
5759566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(map->comm, &sf));
576d61846d3SStefano Zampini   PetscCall(PetscSFSetGraph(sf, n, nleaves, ilocal, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER));
5779566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
5789566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
5798fa9e22eSVaclav Hapla   if (ogidxs) { /* communicate global idxs */
5808fa9e22eSVaclav Hapla     PetscInt cum = 0, start, *work2;
5818fa9e22eSVaclav Hapla 
5829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &work));
5839566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(N, &work2));
5849371c9d4SSatish Balay     for (r = 0; r < N; ++r)
5859371c9d4SSatish Balay       if (idxs[r] >= 0) cum++;
5869566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm));
5878fa9e22eSVaclav Hapla     start -= cum;
5888fa9e22eSVaclav Hapla     cum = 0;
5899371c9d4SSatish Balay     for (r = 0; r < N; ++r)
5909371c9d4SSatish Balay       if (idxs[r] >= 0) work2[r] = start + cum++;
5919566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE));
5929566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE));
5939566063dSJacob Faibussowitsch     PetscCall(PetscFree(work2));
5948fa9e22eSVaclav Hapla   }
5959566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5968fa9e22eSVaclav Hapla   /* Compress and put in indices */
5978fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r)
5988fa9e22eSVaclav Hapla     if (lidxs[r] >= 0) {
5998fa9e22eSVaclav Hapla       if (work) work[len] = work[r];
6008fa9e22eSVaclav Hapla       lidxs[len++] = r;
6018fa9e22eSVaclav Hapla     }
6028fa9e22eSVaclav Hapla   if (on) *on = len;
6038fa9e22eSVaclav Hapla   if (oidxs) *oidxs = lidxs;
6048fa9e22eSVaclav Hapla   if (ogidxs) *ogidxs = work;
6053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6068fa9e22eSVaclav Hapla }
607deffd5ebSksagiyam 
608deffd5ebSksagiyam /*@
609cab54364SBarry Smith   PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices
610deffd5ebSksagiyam 
611deffd5ebSksagiyam   Collective
612deffd5ebSksagiyam 
613deffd5ebSksagiyam   Input Parameters:
614cf84bf9eSBarry Smith + layout           - `PetscLayout` defining the global index space and the MPI rank that brokers each index
615cf84bf9eSBarry Smith . numRootIndices   - size of `rootIndices`
616cf84bf9eSBarry Smith . rootIndices      - array of global indices of which this process requests ownership
617cf84bf9eSBarry Smith . rootLocalIndices - root local index permutation (`NULL` if no permutation)
618cf84bf9eSBarry Smith . rootLocalOffset  - offset to be added to `rootLocalIndices`
619cf84bf9eSBarry Smith . numLeafIndices   - size of `leafIndices`
620cf84bf9eSBarry Smith . leafIndices      - array of global indices with which this process requires data associated
621cf84bf9eSBarry Smith . leafLocalIndices - leaf local index permutation (`NULL` if no permutation)
622cf84bf9eSBarry Smith - leafLocalOffset  - offset to be added to `leafLocalIndices`
623deffd5ebSksagiyam 
624d8d19677SJose E. Roman   Output Parameters:
625cf84bf9eSBarry Smith + sfA - star forest representing the communication pattern from the layout space to the leaf space (`NULL` if not needed)
626deffd5ebSksagiyam - sf  - star forest representing the communication pattern from the root space to the leaf space
627deffd5ebSksagiyam 
628cab54364SBarry Smith   Level: advanced
629cab54364SBarry Smith 
630deffd5ebSksagiyam   Example 1:
631cab54364SBarry Smith .vb
632cab54364SBarry Smith   rank             : 0            1            2
633cab54364SBarry Smith   rootIndices      : [1 0 2]      [3]          [3]
634cab54364SBarry Smith   rootLocalOffset  : 100          200          300
635cab54364SBarry Smith   layout           : [0 1]        [2]          [3]
636cab54364SBarry Smith   leafIndices      : [0]          [2]          [0 3]
637cab54364SBarry Smith   leafLocalOffset  : 400          500          600
638cab54364SBarry Smith 
639cab54364SBarry Smith would build the following PetscSF
640cab54364SBarry Smith 
641cab54364SBarry Smith   [0] 400 <- (0,101)
642cab54364SBarry Smith   [1] 500 <- (0,102)
643cab54364SBarry Smith   [2] 600 <- (0,101)
644cab54364SBarry Smith   [2] 601 <- (2,300)
645cab54364SBarry Smith .ve
646cab54364SBarry Smith 
647deffd5ebSksagiyam   Example 2:
648cab54364SBarry Smith .vb
649cab54364SBarry Smith   rank             : 0               1               2
650cab54364SBarry Smith   rootIndices      : [1 0 2]         [3]             [3]
651cab54364SBarry Smith   rootLocalOffset  : 100             200             300
652cab54364SBarry Smith   layout           : [0 1]           [2]             [3]
653cab54364SBarry Smith   leafIndices      : rootIndices     rootIndices     rootIndices
654cab54364SBarry Smith   leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset
655cab54364SBarry Smith 
656cab54364SBarry Smith would build the following PetscSF
657cab54364SBarry Smith 
658cab54364SBarry Smith   [1] 200 <- (2,300)
659cab54364SBarry Smith .ve
660cab54364SBarry Smith 
661deffd5ebSksagiyam   Example 3:
662cab54364SBarry Smith .vb
663cab54364SBarry Smith   No process requests ownership of global index 1, but no process needs it.
664cab54364SBarry Smith 
665cab54364SBarry Smith   rank             : 0            1            2
666cab54364SBarry Smith   numRootIndices   : 2            1            1
667cab54364SBarry Smith   rootIndices      : [0 2]        [3]          [3]
668cab54364SBarry Smith   rootLocalOffset  : 100          200          300
669cab54364SBarry Smith   layout           : [0 1]        [2]          [3]
670cab54364SBarry Smith   numLeafIndices   : 1            1            2
671cab54364SBarry Smith   leafIndices      : [0]          [2]          [0 3]
672cab54364SBarry Smith   leafLocalOffset  : 400          500          600
673cab54364SBarry Smith 
674cab54364SBarry Smith would build the following PetscSF
675cab54364SBarry Smith 
676cab54364SBarry Smith   [0] 400 <- (0,100)
677cab54364SBarry Smith   [1] 500 <- (0,101)
678cab54364SBarry Smith   [2] 600 <- (0,100)
679cab54364SBarry Smith   [2] 601 <- (2,300)
680cab54364SBarry Smith .ve
681deffd5ebSksagiyam 
682deffd5ebSksagiyam   Notes:
683*2f04c522SBarry Smith   `layout` represents any partitioning of [0, N), where N is the total number of global indices, and its
684cab54364SBarry Smith   local size can be set to `PETSC_DECIDE`.
685cab54364SBarry Smith 
686*2f04c522SBarry Smith   If a global index x lies in the partition owned by process i, each process whose `rootIndices` contains x requests
687deffd5ebSksagiyam   ownership of x and sends its own rank and the local index of x to process i.
688deffd5ebSksagiyam   If multiple processes request ownership of x, the one with the highest rank is to own x.
689*2f04c522SBarry Smith   Process i then broadcasts the ownership information, so that each process whose `leafIndices` contains x knows the
690deffd5ebSksagiyam   ownership information of x.
691*2f04c522SBarry Smith   The output `sf` is constructed by associating each leaf point to a root point in this way.
692deffd5ebSksagiyam 
693deffd5ebSksagiyam   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
694*2f04c522SBarry Smith   The optional output `sfA` can be used to push such data to leaf points.
695deffd5ebSksagiyam 
696*2f04c522SBarry Smith   All indices in `rootIndices` and `leafIndices` must lie in the layout range. The union (over all processes) of `rootIndices`
697*2f04c522SBarry Smith   must cover that of `leafIndices`, but need not cover the entire layout.
698deffd5ebSksagiyam 
699deffd5ebSksagiyam   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
700deffd5ebSksagiyam   star forest is almost identity, so will only include non-trivial part of the map.
701deffd5ebSksagiyam 
70238b5cf2dSJacob Faibussowitsch   Developer Notes:
703deffd5ebSksagiyam   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
704deffd5ebSksagiyam   hash(rank, root_local_index) as the bid for the ownership determination.
705deffd5ebSksagiyam 
706*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`
707deffd5ebSksagiyam @*/
708cf84bf9eSBarry Smith 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)
709d71ae5a4SJacob Faibussowitsch {
710deffd5ebSksagiyam   MPI_Comm     comm = layout->comm;
711deffd5ebSksagiyam   PetscMPIInt  size, rank;
712deffd5ebSksagiyam   PetscSF      sf1;
713deffd5ebSksagiyam   PetscSFNode *owners, *buffer, *iremote;
714deffd5ebSksagiyam   PetscInt    *ilocal, nleaves, N, n, i;
715deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
716deffd5ebSksagiyam   PetscInt N1;
717deffd5ebSksagiyam #endif
718deffd5ebSksagiyam   PetscBool flag;
719deffd5ebSksagiyam 
720deffd5ebSksagiyam   PetscFunctionBegin;
7214f572ea9SToby Isaac   if (rootIndices) PetscAssertPointer(rootIndices, 3);
7224f572ea9SToby Isaac   if (rootLocalIndices) PetscAssertPointer(rootLocalIndices, 4);
7234f572ea9SToby Isaac   if (leafIndices) PetscAssertPointer(leafIndices, 7);
7244f572ea9SToby Isaac   if (leafLocalIndices) PetscAssertPointer(leafLocalIndices, 8);
7254f572ea9SToby Isaac   if (sfA) PetscAssertPointer(sfA, 10);
7264f572ea9SToby Isaac   PetscAssertPointer(sf, 11);
72708401ef6SPierre Jolivet   PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
72808401ef6SPierre Jolivet   PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
7299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7319566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(layout));
7329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(layout, &N));
7339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(layout, &n));
734deffd5ebSksagiyam   flag = (PetscBool)(leafIndices == rootIndices);
7355440e5dcSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPI_C_BOOL, MPI_LAND, comm));
736c9cc58a2SBarry Smith   PetscCheck(!flag || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
737deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
7381690c2aeSBarry Smith   N1 = PETSC_INT_MIN;
7399371c9d4SSatish Balay   for (i = 0; i < numRootIndices; i++)
7409371c9d4SSatish Balay     if (rootIndices[i] > N1) N1 = rootIndices[i];
741462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
74208401ef6SPierre 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);
743deffd5ebSksagiyam   if (!flag) {
7441690c2aeSBarry Smith     N1 = PETSC_INT_MIN;
7459371c9d4SSatish Balay     for (i = 0; i < numLeafIndices; i++)
7469371c9d4SSatish Balay       if (leafIndices[i] > N1) N1 = leafIndices[i];
747462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
74808401ef6SPierre 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);
749deffd5ebSksagiyam   }
750deffd5ebSksagiyam #endif
751deffd5ebSksagiyam   /* Reduce: owners -> buffer */
7529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &buffer));
7539566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, &sf1));
7549566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf1));
7559566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices));
7569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootIndices, &owners));
757deffd5ebSksagiyam   for (i = 0; i < numRootIndices; ++i) {
758deffd5ebSksagiyam     owners[i].rank  = rank;
759deffd5ebSksagiyam     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
760deffd5ebSksagiyam   }
761deffd5ebSksagiyam   for (i = 0; i < n; ++i) {
762deffd5ebSksagiyam     buffer[i].index = -1;
763deffd5ebSksagiyam     buffer[i].rank  = -1;
764deffd5ebSksagiyam   }
7656497c311SBarry Smith   PetscCall(PetscSFReduceBegin(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC));
7666497c311SBarry Smith   PetscCall(PetscSFReduceEnd(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC));
767deffd5ebSksagiyam   /* Bcast: buffer -> owners */
768deffd5ebSksagiyam   if (!flag) {
769deffd5ebSksagiyam     /* leafIndices is different from rootIndices */
7709566063dSJacob Faibussowitsch     PetscCall(PetscFree(owners));
7719566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices));
7729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeafIndices, &owners));
773deffd5ebSksagiyam   }
7746497c311SBarry Smith   PetscCall(PetscSFBcastBegin(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE));
7756497c311SBarry Smith   PetscCall(PetscSFBcastEnd(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE));
77608401ef6SPierre 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]);
7779566063dSJacob Faibussowitsch   PetscCall(PetscFree(buffer));
7789371c9d4SSatish Balay   if (sfA) {
7799371c9d4SSatish Balay     *sfA = sf1;
7809371c9d4SSatish Balay   } else PetscCall(PetscSFDestroy(&sf1));
781deffd5ebSksagiyam   /* Create sf */
782deffd5ebSksagiyam   if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
783deffd5ebSksagiyam     /* leaf space == root space */
7849371c9d4SSatish Balay     for (i = 0, nleaves = 0; i < numLeafIndices; ++i)
7859371c9d4SSatish Balay       if (owners[i].rank != rank) ++nleaves;
7869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &ilocal));
7879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &iremote));
788deffd5ebSksagiyam     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
789deffd5ebSksagiyam       if (owners[i].rank != rank) {
790deffd5ebSksagiyam         ilocal[nleaves]        = leafLocalOffset + i;
791deffd5ebSksagiyam         iremote[nleaves].rank  = owners[i].rank;
792deffd5ebSksagiyam         iremote[nleaves].index = owners[i].index;
793deffd5ebSksagiyam         ++nleaves;
794deffd5ebSksagiyam       }
795deffd5ebSksagiyam     }
7969566063dSJacob Faibussowitsch     PetscCall(PetscFree(owners));
797deffd5ebSksagiyam   } else {
798deffd5ebSksagiyam     nleaves = numLeafIndices;
7999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &ilocal));
800ad540459SPierre Jolivet     for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);
801deffd5ebSksagiyam     iremote = owners;
802deffd5ebSksagiyam   }
8039566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, sf));
8049566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sf));
8059566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
8063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
807deffd5ebSksagiyam }
808fbc7033fSJed Brown 
809fbc7033fSJed Brown /*@
81053c0d4aeSBarry Smith   PetscSFMerge - append/merge indices of `sfb` into `sfa`, with preference for `sfb`
811fbc7033fSJed Brown 
812fbc7033fSJed Brown   Collective
813fbc7033fSJed Brown 
81438b5cf2dSJacob Faibussowitsch   Input Parameters:
815fbc7033fSJed Brown + sfa - default `PetscSF`
816*2f04c522SBarry Smith - sfb - additional edges to add/replace edges in `sfa`
817fbc7033fSJed Brown 
81838b5cf2dSJacob Faibussowitsch   Output Parameter:
819fbc7033fSJed Brown . merged - new `PetscSF` with combined edges
820fbc7033fSJed Brown 
82153c0d4aeSBarry Smith   Level: intermediate
82253c0d4aeSBarry Smith 
823*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCompose()`
824fbc7033fSJed Brown @*/
825fbc7033fSJed Brown PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged)
826fbc7033fSJed Brown {
827fbc7033fSJed Brown   PetscInt maxleaf;
828fbc7033fSJed Brown 
829fbc7033fSJed Brown   PetscFunctionBegin;
830fbc7033fSJed Brown   PetscValidHeaderSpecific(sfa, PETSCSF_CLASSID, 1);
831fbc7033fSJed Brown   PetscValidHeaderSpecific(sfb, PETSCSF_CLASSID, 2);
832fbc7033fSJed Brown   PetscCheckSameComm(sfa, 1, sfb, 2);
8334f572ea9SToby Isaac   PetscAssertPointer(merged, 3);
834fbc7033fSJed Brown   {
835fbc7033fSJed Brown     PetscInt aleaf, bleaf;
836fbc7033fSJed Brown     PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf));
837fbc7033fSJed Brown     PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf));
838fbc7033fSJed Brown     maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index
839fbc7033fSJed Brown   }
840fbc7033fSJed Brown   PetscInt          *clocal, aroots, aleaves, broots, bleaves;
841fbc7033fSJed Brown   PetscSFNode       *cremote;
842fbc7033fSJed Brown   const PetscInt    *alocal, *blocal;
843fbc7033fSJed Brown   const PetscSFNode *aremote, *bremote;
844fbc7033fSJed Brown   PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote));
845fbc7033fSJed Brown   for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1;
846fbc7033fSJed Brown   PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote));
847fbc7033fSJed Brown   PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote));
848fbc7033fSJed Brown   PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space");
849fbc7033fSJed Brown   for (PetscInt i = 0; i < aleaves; i++) {
850fbc7033fSJed Brown     PetscInt a = alocal ? alocal[i] : i;
851fbc7033fSJed Brown     clocal[a]  = a;
852fbc7033fSJed Brown     cremote[a] = aremote[i];
853fbc7033fSJed Brown   }
854fbc7033fSJed Brown   for (PetscInt i = 0; i < bleaves; i++) {
855fbc7033fSJed Brown     PetscInt b = blocal ? blocal[i] : i;
856fbc7033fSJed Brown     clocal[b]  = b;
857fbc7033fSJed Brown     cremote[b] = bremote[i];
858fbc7033fSJed Brown   }
859fbc7033fSJed Brown   PetscInt nleaves = 0;
860fbc7033fSJed Brown   for (PetscInt i = 0; i < maxleaf; i++) {
861fbc7033fSJed Brown     if (clocal[i] < 0) continue;
862fbc7033fSJed Brown     clocal[nleaves]  = clocal[i];
863fbc7033fSJed Brown     cremote[nleaves] = cremote[i];
864fbc7033fSJed Brown     nleaves++;
865fbc7033fSJed Brown   }
866fbc7033fSJed Brown   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged));
867fbc7033fSJed Brown   PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES));
868fbc7033fSJed Brown   PetscCall(PetscFree2(clocal, cremote));
8693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
870fbc7033fSJed Brown }
8710dd791a8SStefano Zampini 
8720dd791a8SStefano Zampini /*@
8730dd791a8SStefano Zampini   PetscSFCreateStridedSF - Create an `PetscSF` to communicate interleaved blocks of data
8740dd791a8SStefano Zampini 
8750dd791a8SStefano Zampini   Collective
8760dd791a8SStefano Zampini 
8770dd791a8SStefano Zampini   Input Parameters:
8780dd791a8SStefano Zampini + sf  - star forest
8790dd791a8SStefano Zampini . bs  - stride
8800dd791a8SStefano Zampini . ldr - leading dimension of root space
8810dd791a8SStefano Zampini - ldl - leading dimension of leaf space
8820dd791a8SStefano Zampini 
8830dd791a8SStefano Zampini   Output Parameter:
8840dd791a8SStefano Zampini . vsf - the new `PetscSF`
8850dd791a8SStefano Zampini 
8860dd791a8SStefano Zampini   Level: intermediate
8870dd791a8SStefano Zampini 
8880dd791a8SStefano Zampini   Notes:
889*2f04c522SBarry Smith   This can be useful to perform communications on multiple right-hand sides stored in a Fortran-style two dimensional array.
890*2f04c522SBarry Smith   For example, the calling sequence
8910dd791a8SStefano Zampini .vb
8920dd791a8SStefano Zampini   c_datatype *roots, *leaves;
8930dd791a8SStefano Zampini   for i in [0,bs) do
8940dd791a8SStefano Zampini     PetscSFBcastBegin(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op)
8950dd791a8SStefano Zampini     PetscSFBcastEnd(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op)
8960dd791a8SStefano Zampini .ve
8970dd791a8SStefano Zampini   is equivalent to
8980dd791a8SStefano Zampini .vb
8990dd791a8SStefano Zampini   c_datatype *roots, *leaves;
9000dd791a8SStefano Zampini   PetscSFCreateStridedSF(sf, bs, ldr, ldl, &vsf)
9010dd791a8SStefano Zampini   PetscSFBcastBegin(vsf, mpi_datatype, roots, leaves, op)
9020dd791a8SStefano Zampini   PetscSFBcastEnd(vsf, mpi_datatype, roots, leaves, op)
9030dd791a8SStefano Zampini .ve
9040dd791a8SStefano Zampini 
9050dd791a8SStefano Zampini   Developer Notes:
9060dd791a8SStefano Zampini   Should this functionality be handled with a new API instead of creating a new object?
9070dd791a8SStefano Zampini 
908*2f04c522SBarry Smith .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()`
9090dd791a8SStefano Zampini @*/
9100dd791a8SStefano Zampini PetscErrorCode PetscSFCreateStridedSF(PetscSF sf, PetscInt bs, PetscInt ldr, PetscInt ldl, PetscSF *vsf)
9110dd791a8SStefano Zampini {
9120dd791a8SStefano Zampini   PetscSF            rankssf;
9130dd791a8SStefano Zampini   const PetscSFNode *iremote, *sfrremote;
9140dd791a8SStefano Zampini   PetscSFNode       *viremote;
9150dd791a8SStefano Zampini   const PetscInt    *ilocal;
9160dd791a8SStefano Zampini   PetscInt          *vilocal = NULL, *ldrs;
9170dd791a8SStefano Zampini   PetscInt           nranks, nr, nl, vnr, vnl, maxl;
9180dd791a8SStefano Zampini   PetscMPIInt        rank;
9190dd791a8SStefano Zampini   MPI_Comm           comm;
9200dd791a8SStefano Zampini   PetscSFType        sftype;
9210dd791a8SStefano Zampini 
9220dd791a8SStefano Zampini   PetscFunctionBegin;
9230dd791a8SStefano Zampini   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
9240dd791a8SStefano Zampini   PetscValidLogicalCollectiveInt(sf, bs, 2);
9250dd791a8SStefano Zampini   PetscAssertPointer(vsf, 5);
9260dd791a8SStefano Zampini   if (bs == 1) {
9270dd791a8SStefano Zampini     PetscCall(PetscObjectReference((PetscObject)sf));
9280dd791a8SStefano Zampini     *vsf = sf;
9290dd791a8SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
9300dd791a8SStefano Zampini   }
9310dd791a8SStefano Zampini   PetscCall(PetscSFSetUp(sf));
9320dd791a8SStefano Zampini   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
9330dd791a8SStefano Zampini   PetscCallMPI(MPI_Comm_rank(comm, &rank));
9340dd791a8SStefano Zampini   PetscCall(PetscSFGetGraph(sf, &nr, &nl, &ilocal, &iremote));
9350dd791a8SStefano Zampini   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxl));
9360dd791a8SStefano Zampini   maxl += 1;
9370dd791a8SStefano Zampini   if (ldl == PETSC_DECIDE) ldl = maxl;
9380dd791a8SStefano Zampini   if (ldr == PETSC_DECIDE) ldr = nr;
9390dd791a8SStefano Zampini   PetscCheck(ldr >= nr, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid leading dimension %" PetscInt_FMT " must be smaller than number of roots %" PetscInt_FMT, ldr, nr);
9400dd791a8SStefano Zampini   PetscCheck(ldl >= maxl, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid leading dimension %" PetscInt_FMT " must be larger than leaf range %" PetscInt_FMT, ldl, maxl - 1);
9410dd791a8SStefano Zampini   vnr = nr * bs;
9420dd791a8SStefano Zampini   vnl = nl * bs;
9430dd791a8SStefano Zampini   PetscCall(PetscMalloc1(vnl, &viremote));
9440dd791a8SStefano Zampini   PetscCall(PetscMalloc1(vnl, &vilocal));
9450dd791a8SStefano Zampini 
9460dd791a8SStefano Zampini   /* Communicate root leading dimensions to leaf ranks */
9470dd791a8SStefano Zampini   PetscCall(PetscSFGetRanksSF(sf, &rankssf));
9480dd791a8SStefano Zampini   PetscCall(PetscSFGetGraph(rankssf, NULL, &nranks, NULL, &sfrremote));
9490dd791a8SStefano Zampini   PetscCall(PetscMalloc1(nranks, &ldrs));
9500dd791a8SStefano Zampini   PetscCall(PetscSFBcastBegin(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE));
9510dd791a8SStefano Zampini   PetscCall(PetscSFBcastEnd(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE));
9520dd791a8SStefano Zampini 
9530dd791a8SStefano Zampini   for (PetscInt i = 0, rold = -1, lda = -1; i < nl; i++) {
954835f2295SStefano Zampini     const PetscInt r  = iremote[i].rank;
9550dd791a8SStefano Zampini     const PetscInt ii = iremote[i].index;
9560dd791a8SStefano Zampini 
9570dd791a8SStefano Zampini     if (r == rank) lda = ldr;
9580dd791a8SStefano Zampini     else if (rold != r) {
9590dd791a8SStefano Zampini       PetscInt j;
9600dd791a8SStefano Zampini 
9610dd791a8SStefano Zampini       for (j = 0; j < nranks; j++)
9620dd791a8SStefano Zampini         if (sfrremote[j].rank == r) break;
963835f2295SStefano Zampini       PetscCheck(j < nranks, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to locate neighbor rank %" PetscInt_FMT, r);
9640dd791a8SStefano Zampini       lda = ldrs[j];
9650dd791a8SStefano Zampini     }
9660dd791a8SStefano Zampini     rold = r;
9670dd791a8SStefano Zampini     for (PetscInt v = 0; v < bs; v++) {
9680dd791a8SStefano Zampini       viremote[v * nl + i].rank  = r;
9690dd791a8SStefano Zampini       viremote[v * nl + i].index = v * lda + ii;
9700dd791a8SStefano Zampini       vilocal[v * nl + i]        = v * ldl + (ilocal ? ilocal[i] : i);
9710dd791a8SStefano Zampini     }
9720dd791a8SStefano Zampini   }
9730dd791a8SStefano Zampini   PetscCall(PetscFree(ldrs));
9740dd791a8SStefano Zampini   PetscCall(PetscSFCreate(comm, vsf));
9750dd791a8SStefano Zampini   PetscCall(PetscSFGetType(sf, &sftype));
9760dd791a8SStefano Zampini   PetscCall(PetscSFSetType(*vsf, sftype));
9770dd791a8SStefano Zampini   PetscCall(PetscSFSetGraph(*vsf, vnr, vnl, vilocal, PETSC_OWN_POINTER, viremote, PETSC_OWN_POINTER));
9780dd791a8SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
9790dd791a8SStefano Zampini }
980