Actual source code: sfutils.c

  1: #include <petsc/private/sfimpl.h>
  2: #include <petsc/private/sectionimpl.h>

  4: /*@
  5:   PetscSFSetGraphLayout - Set a `PetscSF` communication pattern using global indices and a `PetscLayout`

  7:   Collective

  9:   Input Parameters:
 10: + sf        - star forest
 11: . layout    - `PetscLayout` defining the global space for roots, i.e. which roots are owned by each MPI process
 12: . nleaves   - number of leaf vertices on the current process, each of these references a root on any MPI process
 13: . ilocal    - locations of leaves in leafdata buffers, pass `NULL` for contiguous storage, that is the locations are in [0,`nleaves`)
 14: . localmode - copy mode for `ilocal`
 15: - gremote   - root vertices in global numbering corresponding to the leaves

 17:   Level: intermediate

 19:   Note:
 20:   Global indices must lie in [0, N) where N is the global size of `layout`.
 21:   Leaf indices in `ilocal` get sorted; this means the user-provided array gets sorted if localmode is `PETSC_OWN_POINTER`.

 23:   Developer Notes:
 24:   Local indices which are the identity permutation in the range [0,`nleaves`) are discarded as they
 25:   encode contiguous storage. In such case, if localmode is `PETSC_OWN_POINTER`, the memory is deallocated as it is not
 26:   needed

 28: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFGetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
 29: @*/
 30: PetscErrorCode PetscSFSetGraphLayout(PetscSF sf, PetscLayout layout, PetscInt nleaves, PetscInt ilocal[], PetscCopyMode localmode, const PetscInt gremote[])
 31: {
 32:   const PetscInt *range;
 33:   PetscInt        i, nroots, ls = -1, ln = -1;
 34:   PetscMPIInt     lr = -1;
 35:   PetscSFNode    *remote;

 37:   PetscFunctionBegin;
 39:   PetscAssertPointer(layout, 2);
 40:   if (nleaves > 0 && ilocal) PetscAssertPointer(ilocal, 4);
 41:   if (nleaves > 0) PetscAssertPointer(gremote, 6);
 42:   PetscCall(PetscLayoutSetUp(layout));
 43:   PetscCall(PetscLayoutGetLocalSize(layout, &nroots));
 44:   PetscCall(PetscLayoutGetRanges(layout, &range));
 45:   PetscCall(PetscMalloc1(nleaves, &remote));
 46:   if (nleaves) ls = gremote[0] + 1;
 47:   for (i = 0; i < nleaves; i++) {
 48:     const PetscInt idx = gremote[i] - ls;
 49:     if (idx < 0 || idx >= ln) { /* short-circuit the search */
 50:       PetscCall(PetscLayoutFindOwnerIndex(layout, gremote[i], &lr, &remote[i].index));
 51:       remote[i].rank = lr;
 52:       ls             = range[lr];
 53:       ln             = range[lr + 1] - ls;
 54:     } else {
 55:       remote[i].rank  = lr;
 56:       remote[i].index = idx;
 57:     }
 58:   }
 59:   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, localmode, remote, PETSC_OWN_POINTER));
 60:   PetscFunctionReturn(PETSC_SUCCESS);
 61: }

 63: /*@C
 64:   PetscSFGetGraphLayout - Get the global indices and `PetscLayout` that describe a `PetscSF`

 66:   Collective

 68:   Input Parameter:
 69: . sf - star forest

 71:   Output Parameters:
 72: + layout  - `PetscLayout` defining the global space for roots
 73: . nleaves - number of leaf vertices on the current process, each of these references a root on any process
 74: . ilocal  - locations of leaves in leafdata buffers, or `NULL` for contiguous storage
 75: - gremote - root vertices in global numbering corresponding to the leaves

 77:   Level: intermediate

 79:   Notes:
 80:   The outputs are such that passing them as inputs to `PetscSFSetGraphLayout()` would lead to the same star forest.
 81:   The outputs `layout` and `gremote` are freshly created each time this function is called,
 82:   so they need to be freed (with `PetscLayoutDestroy()` and `PetscFree()`) by the user.

 84: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFSetGraphLayout()`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
 85: @*/
 86: PetscErrorCode PetscSFGetGraphLayout(PetscSF sf, PetscLayout *layout, PetscInt *nleaves, const PetscInt *ilocal[], PetscInt *gremote[])
 87: {
 88:   PetscInt           nr, nl;
 89:   const PetscSFNode *ir;
 90:   PetscLayout        lt;

 92:   PetscFunctionBegin;
 94:   if (layout) PetscAssertPointer(layout, 2);
 95:   if (nleaves) PetscAssertPointer(nleaves, 3);
 96:   if (ilocal) PetscAssertPointer(ilocal, 4);
 97:   if (gremote) PetscAssertPointer(gremote, 5);
 98:   PetscCall(PetscSFGetGraph(sf, &nr, &nl, ilocal, &ir));
 99:   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sf), nr, PETSC_DECIDE, 1, &lt));
100:   if (gremote) {
101:     const PetscInt *range;
102:     PetscInt       *gr;

104:     PetscCall(PetscLayoutGetRanges(lt, &range));
105:     PetscCall(PetscMalloc1(nl, &gr));
106:     for (PetscInt i = 0; i < nl; i++) gr[i] = range[ir[i].rank] + ir[i].index;
107:     *gremote = gr;
108:   }
109:   if (nleaves) *nleaves = nl;
110:   if (layout) *layout = lt;
111:   else PetscCall(PetscLayoutDestroy(&lt));
112:   PetscFunctionReturn(PETSC_SUCCESS);
113: }

115: /*@
116:   PetscSFSetGraphSection - Sets the `PetscSF` graph (communication pattern) encoding the parallel dof overlap based upon the `PetscSection` describing the data layout.

118:   Input Parameters:
119: + sf            - The `PetscSF`
120: . localSection  - `PetscSection` describing the local data layout
121: - globalSection - `PetscSection` describing the global data layout

123:   Level: developer

125: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFSetGraph()`, `PetscSFSetGraphLayout()`
126: @*/
127: PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
128: {
129:   MPI_Comm        comm;
130:   PetscLayout     layout;
131:   const PetscInt *ranges;
132:   PetscInt       *local;
133:   PetscSFNode    *remote;
134:   PetscInt        pStart, pEnd, p, nroots, nleaves = 0, l;
135:   PetscMPIInt     size, rank;

137:   PetscFunctionBegin;

142:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
143:   PetscCallMPI(MPI_Comm_size(comm, &size));
144:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
145:   PetscCall(PetscSectionGetChart(globalSection, &pStart, &pEnd));
146:   PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
147:   PetscCall(PetscLayoutCreate(comm, &layout));
148:   PetscCall(PetscLayoutSetBlockSize(layout, 1));
149:   PetscCall(PetscLayoutSetLocalSize(layout, nroots));
150:   PetscCall(PetscLayoutSetUp(layout));
151:   PetscCall(PetscLayoutGetRanges(layout, &ranges));
152:   for (p = pStart; p < pEnd; ++p) {
153:     PetscInt gdof, gcdof;

155:     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
156:     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
157:     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);
158:     nleaves += gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
159:   }
160:   PetscCall(PetscMalloc1(nleaves, &local));
161:   PetscCall(PetscMalloc1(nleaves, &remote));
162:   for (p = pStart, l = 0; p < pEnd; ++p) {
163:     const PetscInt *cind;
164:     PetscInt        dof, cdof, off, gdof, gcdof, goff, gsize, d, c;

166:     PetscCall(PetscSectionGetDof(localSection, p, &dof));
167:     PetscCall(PetscSectionGetOffset(localSection, p, &off));
168:     PetscCall(PetscSectionGetConstraintDof(localSection, p, &cdof));
169:     PetscCall(PetscSectionGetConstraintIndices(localSection, p, &cind));
170:     PetscCall(PetscSectionGetDof(globalSection, p, &gdof));
171:     PetscCall(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
172:     PetscCall(PetscSectionGetOffset(globalSection, p, &goff));
173:     if (!gdof) continue; /* Censored point */
174:     gsize = gdof < 0 ? -(gdof + 1) - gcdof : gdof - gcdof;
175:     if (gsize != dof - cdof) {
176:       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);
177:       cdof = 0; /* Ignore constraints */
178:     }
179:     for (d = 0, c = 0; d < dof; ++d) {
180:       if ((c < cdof) && (cind[c] == d)) {
181:         ++c;
182:         continue;
183:       }
184:       local[l + d - c] = off + d;
185:     }
186:     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);
187:     if (gdof < 0) {
188:       for (d = 0; d < gsize; ++d, ++l) {
189:         PetscInt    offset = -(goff + 1) + d, ir;
190:         PetscMPIInt r;

192:         PetscCall(PetscFindInt(offset, size + 1, ranges, &ir));
193:         PetscCall(PetscMPIIntCast(ir, &r));
194:         if (r < 0) r = -(r + 2);
195:         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);
196:         remote[l].rank  = r;
197:         remote[l].index = offset - ranges[r];
198:       }
199:     } else {
200:       for (d = 0; d < gsize; ++d, ++l) {
201:         remote[l].rank  = rank;
202:         remote[l].index = goff + d - ranges[rank];
203:       }
204:     }
205:   }
206:   PetscCheck(l == nleaves, comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
207:   PetscCall(PetscLayoutDestroy(&layout));
208:   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: /*@C
213:   PetscSFDistributeSection - Create a new `PetscSection` reorganized, moving from the root to the leaves of the `PetscSF`

215:   Collective

217:   Input Parameters:
218: + sf          - The `PetscSF`
219: - rootSection - Section defined on root space

221:   Output Parameters:
222: + remoteOffsets - root offsets in leaf storage, or `NULL`, its length will be the size of the chart of `leafSection`
223: - leafSection   - Section defined on the leaf space

225:   Level: advanced

227:   Note:
228:   Caller must `PetscFree()` `remoteOffsets` if it was requested

230:   To distribute data from the `rootSection` to the `leafSection`, see  `PetscSFCreateSectionSF()` or `PetscSectionMigrateData()`.

232:   Fortran Note:
233:   Use `PetscSFDestroyRemoteOffsets()` when `remoteOffsets` is no longer needed.

235: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`, `PetscSFCreateSectionSF()`
236: @*/
237: PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt *remoteOffsets[], PetscSection leafSection)
238: {
239:   PetscSF         embedSF;
240:   const PetscInt *indices;
241:   IS              selected;
242:   PetscInt        numFields, nroots, rpStart, rpEnd, lpStart = PETSC_INT_MAX, lpEnd = -1, f, c;
243:   PetscBool      *sub, hasc;

245:   PetscFunctionBegin;
248:   if (remoteOffsets) PetscAssertPointer(remoteOffsets, 3);
250:   PetscCall(PetscLogEventBegin(PETSCSF_DistSect, sf, 0, 0, 0));
251:   PetscCall(PetscSectionGetNumFields(rootSection, &numFields));
252:   if (numFields) {
253:     IS perm;

255:     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
256:        leafSection->perm. To keep this permutation set by the user, we grab
257:        the reference before calling PetscSectionSetNumFields() and set it
258:        back after. */
259:     PetscCall(PetscSectionGetPermutation(leafSection, &perm));
260:     PetscCall(PetscObjectReference((PetscObject)perm));
261:     PetscCall(PetscSectionSetNumFields(leafSection, numFields));
262:     PetscCall(PetscSectionSetPermutation(leafSection, perm));
263:     PetscCall(ISDestroy(&perm));
264:   }
265:   PetscCall(PetscMalloc1(numFields + 2, &sub));
266:   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
267:   for (f = 0; f < numFields; ++f) {
268:     PetscSectionSym sym, dsym = NULL;
269:     const char     *name    = NULL;
270:     PetscInt        numComp = 0;

272:     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
273:     PetscCall(PetscSectionGetFieldComponents(rootSection, f, &numComp));
274:     PetscCall(PetscSectionGetFieldName(rootSection, f, &name));
275:     PetscCall(PetscSectionGetFieldSym(rootSection, f, &sym));
276:     if (sym) PetscCall(PetscSectionSymDistribute(sym, sf, &dsym));
277:     PetscCall(PetscSectionSetFieldComponents(leafSection, f, numComp));
278:     PetscCall(PetscSectionSetFieldName(leafSection, f, name));
279:     PetscCall(PetscSectionSetFieldSym(leafSection, f, dsym));
280:     PetscCall(PetscSectionSymDestroy(&dsym));
281:     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
282:       PetscCall(PetscSectionGetComponentName(rootSection, f, c, &name));
283:       PetscCall(PetscSectionSetComponentName(leafSection, f, c, name));
284:     }
285:   }
286:   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
287:   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
288:   rpEnd = PetscMin(rpEnd, nroots);
289:   rpEnd = PetscMax(rpStart, rpEnd);
290:   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
291:   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
292:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, sub, 2 + numFields, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf)));
293:   if (sub[0]) {
294:     PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
295:     PetscCall(ISGetIndices(selected, &indices));
296:     PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
297:     PetscCall(ISRestoreIndices(selected, &indices));
298:     PetscCall(ISDestroy(&selected));
299:   } else {
300:     PetscCall(PetscObjectReference((PetscObject)sf));
301:     embedSF = sf;
302:   }
303:   PetscCall(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd));
304:   lpEnd++;

306:   PetscCall(PetscSectionSetChart(leafSection, lpStart, lpEnd));

308:   /* Constrained dof section */
309:   hasc = sub[1];
310:   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2 + f]);

312:   /* Could fuse these at the cost of copies and extra allocation */
313:   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE));
314:   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->atlasDof, -lpStart), MPI_REPLACE));
315:   if (sub[1]) {
316:     PetscCall(PetscSectionCheckConstraints_Private(rootSection));
317:     PetscCall(PetscSectionCheckConstraints_Private(leafSection));
318:     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
319:     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart], MPI_REPLACE));
320:   }
321:   for (f = 0; f < numFields; ++f) {
322:     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE));
323:     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->field[f]->atlasDof, -rpStart), PetscSafePointerPlusOffset(leafSection->field[f]->atlasDof, -lpStart), MPI_REPLACE));
324:     if (sub[2 + f]) {
325:       PetscCall(PetscSectionCheckConstraints_Private(rootSection->field[f]));
326:       PetscCall(PetscSectionCheckConstraints_Private(leafSection->field[f]));
327:       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
328:       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart], MPI_REPLACE));
329:     }
330:   }
331:   if (remoteOffsets) {
332:     PetscCall(PetscMalloc1(lpEnd - lpStart, remoteOffsets));
333:     PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
334:     PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
335:   }
336:   PetscCall(PetscSectionInvalidateMaxDof_Internal(leafSection));
337:   PetscCall(PetscSectionSetUp(leafSection));
338:   if (hasc) { /* need to communicate bcIndices */
339:     PetscSF   bcSF;
340:     PetscInt *rOffBc;

342:     PetscCall(PetscMalloc1(lpEnd - lpStart, &rOffBc));
343:     if (sub[1]) {
344:       PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
345:       PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
346:       PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF));
347:       PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
348:       PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices, MPI_REPLACE));
349:       PetscCall(PetscSFDestroy(&bcSF));
350:     }
351:     for (f = 0; f < numFields; ++f) {
352:       if (sub[2 + f]) {
353:         PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
354:         PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart], MPI_REPLACE));
355:         PetscCall(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF));
356:         PetscCall(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
357:         PetscCall(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices, MPI_REPLACE));
358:         PetscCall(PetscSFDestroy(&bcSF));
359:       }
360:     }
361:     PetscCall(PetscFree(rOffBc));
362:   }
363:   PetscCall(PetscSFDestroy(&embedSF));
364:   PetscCall(PetscFree(sub));
365:   PetscCall(PetscLogEventEnd(PETSCSF_DistSect, sf, 0, 0, 0));
366:   PetscFunctionReturn(PETSC_SUCCESS);
367: }

369: /*@C
370:   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes

372:   Collective

374:   Input Parameters:
375: + sf          - The `PetscSF`
376: . rootSection - Data layout of remote points for outgoing data (this is layout for roots)
377: - leafSection - Data layout of local points for incoming data  (this is layout for leaves)

379:   Output Parameter:
380: . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or `NULL`

382:   Level: developer

384:   Note:
385:   Caller must `PetscFree()` `remoteOffsets` if it was requested

387:   Fortran Note:
388:   Use `PetscSFDestroyRemoteOffsets()` when `remoteOffsets` is no longer needed.

390: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`
391: @*/
392: PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt *remoteOffsets[])
393: {
394:   PetscSF         embedSF;
395:   const PetscInt *indices;
396:   IS              selected;
397:   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;

399:   PetscFunctionBegin;
403:   PetscAssertPointer(remoteOffsets, 4);
404:   *remoteOffsets = NULL;
405:   PetscCall(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL));
406:   if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS);
407:   PetscCall(PetscLogEventBegin(PETSCSF_RemoteOff, sf, 0, 0, 0));
408:   PetscCall(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
409:   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
410:   PetscCall(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
411:   PetscCall(ISGetIndices(selected, &indices));
412:   PetscCall(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
413:   PetscCall(ISRestoreIndices(selected, &indices));
414:   PetscCall(ISDestroy(&selected));
415:   PetscCall(PetscCalloc1(lpEnd - lpStart, remoteOffsets));
416:   PetscCall(PetscSFBcastBegin(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
417:   PetscCall(PetscSFBcastEnd(embedSF, MPIU_INT, PetscSafePointerPlusOffset(rootSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(*remoteOffsets, -lpStart), MPI_REPLACE));
418:   PetscCall(PetscSFDestroy(&embedSF));
419:   PetscCall(PetscLogEventEnd(PETSCSF_RemoteOff, sf, 0, 0, 0));
420:   PetscFunctionReturn(PETSC_SUCCESS);
421: }

423: /*@
424:   PetscSFCreateSectionSF - Create an expanded `PetscSF` of dofs, assuming the input `PetscSF` relates points

426:   Collective

428:   Input Parameters:
429: + sf            - The `PetscSF`
430: . rootSection   - Data layout of remote points for outgoing data (this is usually the serial section)
431: . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or `NULL`
432: - leafSection   - Data layout of local points for incoming data  (this is the distributed section)

434:   Output Parameter:
435: . sectionSF - The new `PetscSF`

437:   Level: advanced

439:   Notes:
440:   `remoteOffsets` can be `NULL` if `sf` does not reference any points in `leafSection`

442: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`, `PetscSFDistributeSection()`
443: @*/
444: PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
445: {
446:   MPI_Comm           comm;
447:   const PetscInt    *localPoints;
448:   const PetscSFNode *remotePoints;
449:   PetscInt           lpStart, lpEnd;
450:   PetscInt           numRoots, numSectionRoots, numPoints, numIndices = 0;
451:   PetscInt          *localIndices;
452:   PetscSFNode       *remoteIndices;
453:   PetscInt           i, ind;

455:   PetscFunctionBegin;
457:   PetscAssertPointer(rootSection, 2);
458:   /* Cannot check PetscAssertPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
459:   PetscAssertPointer(leafSection, 4);
460:   PetscAssertPointer(sectionSF, 5);
461:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
462:   PetscCall(PetscSFCreate(comm, sectionSF));
463:   PetscCall(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
464:   PetscCall(PetscSectionGetStorageSize(rootSection, &numSectionRoots));
465:   PetscCall(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints));
466:   if (numRoots < 0) PetscFunctionReturn(PETSC_SUCCESS);
467:   PetscCall(PetscLogEventBegin(PETSCSF_SectSF, sf, 0, 0, 0));
468:   for (i = 0; i < numPoints; ++i) {
469:     PetscInt localPoint = localPoints ? localPoints[i] : i;
470:     PetscInt dof;

472:     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
473:       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
474:       numIndices += dof < 0 ? 0 : dof;
475:     }
476:   }
477:   PetscCall(PetscMalloc1(numIndices, &localIndices));
478:   PetscCall(PetscMalloc1(numIndices, &remoteIndices));
479:   /* Create new index graph */
480:   for (i = 0, ind = 0; i < numPoints; ++i) {
481:     PetscInt localPoint = localPoints ? localPoints[i] : i;
482:     PetscInt rank       = remotePoints[i].rank;

484:     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
485:       PetscInt remoteOffset = remoteOffsets[localPoint - lpStart];
486:       PetscInt loff, dof, d;

488:       PetscCall(PetscSectionGetOffset(leafSection, localPoint, &loff));
489:       PetscCall(PetscSectionGetDof(leafSection, localPoint, &dof));
490:       for (d = 0; d < dof; ++d, ++ind) {
491:         localIndices[ind]        = loff + d;
492:         remoteIndices[ind].rank  = rank;
493:         remoteIndices[ind].index = remoteOffset + d;
494:       }
495:     }
496:   }
497:   PetscCheck(numIndices == ind, comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
498:   PetscCall(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER));
499:   PetscCall(PetscSFSetUp(*sectionSF));
500:   PetscCall(PetscLogEventEnd(PETSCSF_SectSF, sf, 0, 0, 0));
501:   PetscFunctionReturn(PETSC_SUCCESS);
502: }

504: /*@C
505:   PetscSFCreateFromLayouts - Creates a parallel star forest mapping between two `PetscLayout` objects

507:   Collective

509:   Input Parameters:
510: + rmap - `PetscLayout` defining the global root space
511: - lmap - `PetscLayout` defining the global leaf space

513:   Output Parameter:
514: . sf - The parallel star forest

516:   Level: intermediate

518:   Notes:
519:   If the global length of `lmap` differs from the global length of `rmap` then the excess entries are ignored.

521:   The resulting `sf` used with `PetscSFBcastBegin()` and `PetscSFBcastEnd()` merely copies the array entries of `rootdata` to
522:   `leafdata`; moving them between MPI processes if needed. For example,
523:   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
524:   `leafdata` would become (1, 2) on MPI rank 0 and (3, 4, 5, x) on MPI rank 1.

526: .seealso: [](sec_petscsf), `PetscSF`, `PetscLayout`, `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscLayoutCreate()`, `PetscSFSetGraphLayout()`
527: @*/
528: PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF *sf)
529: {
530:   PetscInt     i, nroots, nleaves = 0;
531:   PetscInt     rN, lst, len;
532:   PetscMPIInt  owner = -1;
533:   PetscSFNode *remote;
534:   MPI_Comm     rcomm = rmap->comm;
535:   MPI_Comm     lcomm = lmap->comm;
536:   PetscMPIInt  flg;

538:   PetscFunctionBegin;
539:   PetscAssertPointer(rmap, 1);
540:   PetscAssertPointer(lmap, 2);
541:   PetscAssertPointer(sf, 3);
542:   PetscCheck(rmap->setupcalled, rcomm, PETSC_ERR_ARG_WRONGSTATE, "Root layout not setup");
543:   PetscCheck(lmap->setupcalled, lcomm, PETSC_ERR_ARG_WRONGSTATE, "Leaf layout not setup");
544:   PetscCallMPI(MPI_Comm_compare(rcomm, lcomm, &flg));
545:   PetscCheck(flg == MPI_CONGRUENT || flg == MPI_IDENT, rcomm, PETSC_ERR_SUP, "cannot map two layouts with non-matching communicators");
546:   PetscCall(PetscSFCreate(rcomm, sf));
547:   PetscCall(PetscLayoutGetLocalSize(rmap, &nroots));
548:   PetscCall(PetscLayoutGetSize(rmap, &rN));
549:   PetscCall(PetscLayoutGetRange(lmap, &lst, &len));
550:   PetscCall(PetscMalloc1(len - lst, &remote));
551:   for (i = lst; i < len && i < rN; i++) {
552:     if (owner < -1 || i >= rmap->range[owner + 1]) PetscCall(PetscLayoutFindOwner(rmap, i, &owner));
553:     remote[nleaves].rank  = owner;
554:     remote[nleaves].index = i - rmap->range[owner];
555:     nleaves++;
556:   }
557:   PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, remote, PETSC_COPY_VALUES));
558:   PetscCall(PetscFree(remote));
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

562: /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
563: PetscErrorCode PetscLayoutMapLocal(PetscLayout map, PetscInt N, const PetscInt idxs[], PetscInt *on, PetscInt *oidxs[], PetscInt *ogidxs[])
564: {
565:   PetscInt    *owners = map->range;
566:   PetscInt     n      = map->n;
567:   PetscSF      sf;
568:   PetscInt    *lidxs, *work = NULL, *ilocal;
569:   PetscSFNode *ridxs;
570:   PetscMPIInt  rank, p = 0;
571:   PetscInt     r, len = 0, nleaves = 0;

573:   PetscFunctionBegin;
574:   if (on) *on = 0; /* squelch -Wmaybe-uninitialized */
575:   /* Create SF where leaves are input idxs and roots are owned idxs */
576:   PetscCallMPI(MPI_Comm_rank(map->comm, &rank));
577:   PetscCall(PetscMalloc1(n, &lidxs));
578:   for (r = 0; r < n; ++r) lidxs[r] = -1;
579:   PetscCall(PetscMalloc1(N, &ridxs));
580:   PetscCall(PetscMalloc1(N, &ilocal));
581:   for (r = 0; r < N; ++r) {
582:     const PetscInt idx = idxs[r];

584:     if (idx < 0) continue;
585:     PetscCheck(idx < map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")", idx, map->N);
586:     if (idx < owners[p] || owners[p + 1] <= idx) { /* short-circuit the search if the last p owns this idx too */
587:       PetscCall(PetscLayoutFindOwner(map, idx, &p));
588:     }
589:     ridxs[nleaves].rank  = p;
590:     ridxs[nleaves].index = idxs[r] - owners[p];
591:     ilocal[nleaves]      = r;
592:     nleaves++;
593:   }
594:   PetscCall(PetscSFCreate(map->comm, &sf));
595:   PetscCall(PetscSFSetGraph(sf, n, nleaves, ilocal, PETSC_OWN_POINTER, ridxs, PETSC_OWN_POINTER));
596:   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
597:   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *)idxs, lidxs, MPI_LOR));
598:   if (ogidxs) { /* communicate global idxs */
599:     PetscInt cum = 0, start, *work2;

601:     PetscCall(PetscMalloc1(n, &work));
602:     PetscCall(PetscCalloc1(N, &work2));
603:     for (r = 0; r < N; ++r)
604:       if (idxs[r] >= 0) cum++;
605:     PetscCallMPI(MPI_Scan(&cum, &start, 1, MPIU_INT, MPI_SUM, map->comm));
606:     start -= cum;
607:     cum = 0;
608:     for (r = 0; r < N; ++r)
609:       if (idxs[r] >= 0) work2[r] = start + cum++;
610:     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, work2, work, MPI_REPLACE));
611:     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, work2, work, MPI_REPLACE));
612:     PetscCall(PetscFree(work2));
613:   }
614:   PetscCall(PetscSFDestroy(&sf));
615:   /* Compress and put in indices */
616:   for (r = 0; r < n; ++r)
617:     if (lidxs[r] >= 0) {
618:       if (work) work[len] = work[r];
619:       lidxs[len++] = r;
620:     }
621:   if (on) *on = len;
622:   if (oidxs) *oidxs = lidxs;
623:   if (ogidxs) *ogidxs = work;
624:   PetscFunctionReturn(PETSC_SUCCESS);
625: }

627: /*@
628:   PetscSFCreateByMatchingIndices - Create `PetscSF` by matching root and leaf indices

630:   Collective

632:   Input Parameters:
633: + layout           - `PetscLayout` defining the global index space and the MPI rank that brokers each index
634: . numRootIndices   - size of `rootIndices`
635: . rootIndices      - array of global indices of which this process requests ownership
636: . rootLocalIndices - root local index permutation (`NULL` if no permutation)
637: . rootLocalOffset  - offset to be added to `rootLocalIndices`
638: . numLeafIndices   - size of `leafIndices`
639: . leafIndices      - array of global indices with which this process requires data associated
640: . leafLocalIndices - leaf local index permutation (`NULL` if no permutation)
641: - leafLocalOffset  - offset to be added to `leafLocalIndices`

643:   Output Parameters:
644: + sfA - star forest representing the communication pattern from the layout space to the leaf space (`NULL` if not needed)
645: - sf  - star forest representing the communication pattern from the root space to the leaf space

647:   Level: advanced

649:   Example 1:
650: .vb
651:   rank             : 0            1            2
652:   rootIndices      : [1 0 2]      [3]          [3]
653:   rootLocalOffset  : 100          200          300
654:   layout           : [0 1]        [2]          [3]
655:   leafIndices      : [0]          [2]          [0 3]
656:   leafLocalOffset  : 400          500          600

658: would build the following PetscSF

660:   [0] 400 <- (0,101)
661:   [1] 500 <- (0,102)
662:   [2] 600 <- (0,101)
663:   [2] 601 <- (2,300)
664: .ve

666:   Example 2:
667: .vb
668:   rank             : 0               1               2
669:   rootIndices      : [1 0 2]         [3]             [3]
670:   rootLocalOffset  : 100             200             300
671:   layout           : [0 1]           [2]             [3]
672:   leafIndices      : rootIndices     rootIndices     rootIndices
673:   leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset

675: would build the following PetscSF

677:   [1] 200 <- (2,300)
678: .ve

680:   Example 3:
681: .vb
682:   No process requests ownership of global index 1, but no process needs it.

684:   rank             : 0            1            2
685:   numRootIndices   : 2            1            1
686:   rootIndices      : [0 2]        [3]          [3]
687:   rootLocalOffset  : 100          200          300
688:   layout           : [0 1]        [2]          [3]
689:   numLeafIndices   : 1            1            2
690:   leafIndices      : [0]          [2]          [0 3]
691:   leafLocalOffset  : 400          500          600

693: would build the following PetscSF

695:   [0] 400 <- (0,100)
696:   [1] 500 <- (0,101)
697:   [2] 600 <- (0,100)
698:   [2] 601 <- (2,300)
699: .ve

701:   Notes:
702:   `layout` represents any partitioning of [0, N), where N is the total number of global indices, and its
703:   local size can be set to `PETSC_DECIDE`.

705:   If a global index x lies in the partition owned by process i, each process whose `rootIndices` contains x requests
706:   ownership of x and sends its own rank and the local index of x to process i.
707:   If multiple processes request ownership of x, the one with the highest rank is to own x.
708:   Process i then broadcasts the ownership information, so that each process whose `leafIndices` contains x knows the
709:   ownership information of x.
710:   The output `sf` is constructed by associating each leaf point to a root point in this way.

712:   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
713:   The optional output `sfA` can be used to push such data to leaf points.

715:   All indices in `rootIndices` and `leafIndices` must lie in the layout range. The union (over all processes) of `rootIndices`
716:   must cover that of `leafIndices`, but need not cover the entire layout.

718:   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
719:   star forest is almost identity, so will only include non-trivial part of the map.

721:   Developer Notes:
722:   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
723:   hash(rank, root_local_index) as the bid for the ownership determination.

725: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`
726: @*/
727: 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)
728: {
729:   MPI_Comm     comm = layout->comm;
730:   PetscMPIInt  rank;
731:   PetscSF      sf1;
732:   PetscSFNode *owners, *buffer, *iremote;
733:   PetscInt    *ilocal, nleaves, N, n, i;
734:   PetscBool    areIndicesSame;

736:   PetscFunctionBegin;
737:   PetscAssertPointer(layout, 1);
738:   if (rootIndices) PetscAssertPointer(rootIndices, 3);
739:   if (rootLocalIndices) PetscAssertPointer(rootLocalIndices, 4);
740:   if (leafIndices) PetscAssertPointer(leafIndices, 7);
741:   if (leafLocalIndices) PetscAssertPointer(leafLocalIndices, 8);
742:   if (sfA) PetscAssertPointer(sfA, 10);
743:   PetscAssertPointer(sf, 11);
744:   PetscCheck(numRootIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
745:   PetscCheck(numLeafIndices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
746:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
747:   PetscCall(PetscLayoutSetUp(layout));
748:   PetscCall(PetscLayoutGetSize(layout, &N));
749:   PetscCall(PetscLayoutGetLocalSize(layout, &n));
750:   areIndicesSame = (PetscBool)(leafIndices == rootIndices);
751:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &areIndicesSame, 1, MPI_C_BOOL, MPI_LAND, comm));
752:   PetscCheck(!areIndicesSame || numLeafIndices == numRootIndices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
753:   if (PetscDefined(USE_DEBUG)) {
754:     PetscInt N1 = PETSC_INT_MIN;
755:     for (i = 0; i < numRootIndices; i++)
756:       if (rootIndices[i] > N1) N1 = rootIndices[i];
757:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
758:     PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
759:     if (!areIndicesSame) {
760:       N1 = PETSC_INT_MIN;
761:       for (i = 0; i < numLeafIndices; i++)
762:         if (leafIndices[i] > N1) N1 = leafIndices[i];
763:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
764:       PetscCheck(N1 < N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
765:     }
766:   }

768:   /* Reduce: owners -> buffer */
769:   PetscCall(PetscMalloc1(n, &buffer));
770:   PetscCall(PetscSFCreate(comm, &sf1));
771:   PetscCall(PetscSFSetFromOptions(sf1));
772:   PetscCall(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices));
773:   PetscCall(PetscMalloc1(numRootIndices, &owners));
774:   for (i = 0; i < numRootIndices; ++i) {
775:     owners[i].rank  = rank;
776:     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
777:   }
778:   for (i = 0; i < n; ++i) {
779:     buffer[i].index = -1;
780:     buffer[i].rank  = -1;
781:   }
782:   PetscCall(PetscSFReduceBegin(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC));
783:   PetscCall(PetscSFReduceEnd(sf1, MPIU_SF_NODE, owners, buffer, MPI_MAXLOC));
784:   /* Bcast: buffer -> owners */
785:   if (!areIndicesSame) {
786:     PetscCall(PetscFree(owners));
787:     PetscCall(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices));
788:     PetscCall(PetscMalloc1(numLeafIndices, &owners));
789:   }
790:   PetscCall(PetscSFBcastBegin(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE));
791:   PetscCall(PetscSFBcastEnd(sf1, MPIU_SF_NODE, buffer, owners, MPI_REPLACE));
792:   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]);
793:   PetscCall(PetscFree(buffer));
794:   if (sfA) *sfA = sf1;
795:   else PetscCall(PetscSFDestroy(&sf1));
796:   /* Create sf */
797:   if (areIndicesSame && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
798:     /* leaf space == root space */
799:     for (i = 0, nleaves = 0; i < numLeafIndices; ++i)
800:       if (owners[i].rank != rank) ++nleaves;
801:     PetscCall(PetscMalloc1(nleaves, &ilocal));
802:     PetscCall(PetscMalloc1(nleaves, &iremote));
803:     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
804:       if (owners[i].rank != rank) {
805:         ilocal[nleaves]        = leafLocalOffset + i;
806:         iremote[nleaves].rank  = owners[i].rank;
807:         iremote[nleaves].index = owners[i].index;
808:         ++nleaves;
809:       }
810:     }
811:     PetscCall(PetscFree(owners));
812:   } else {
813:     nleaves = numLeafIndices;
814:     PetscCall(PetscMalloc1(nleaves, &ilocal));
815:     for (i = 0; i < nleaves; ++i) ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);
816:     iremote = owners;
817:   }
818:   PetscCall(PetscSFCreate(comm, sf));
819:   PetscCall(PetscSFSetFromOptions(*sf));
820:   PetscCall(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
821:   PetscFunctionReturn(PETSC_SUCCESS);
822: }

824: /*@
825:   PetscSFMerge - append/merge indices of `sfb` into `sfa`, with preference for `sfb`

827:   Collective

829:   Input Parameters:
830: + sfa - default `PetscSF`
831: - sfb - additional edges to add/replace edges in `sfa`

833:   Output Parameter:
834: . merged - new `PetscSF` with combined edges

836:   Level: intermediate

838: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCompose()`
839: @*/
840: PetscErrorCode PetscSFMerge(PetscSF sfa, PetscSF sfb, PetscSF *merged)
841: {
842:   PetscInt maxleaf;

844:   PetscFunctionBegin;
847:   PetscCheckSameComm(sfa, 1, sfb, 2);
848:   PetscAssertPointer(merged, 3);
849:   {
850:     PetscInt aleaf, bleaf;
851:     PetscCall(PetscSFGetLeafRange(sfa, NULL, &aleaf));
852:     PetscCall(PetscSFGetLeafRange(sfb, NULL, &bleaf));
853:     maxleaf = PetscMax(aleaf, bleaf) + 1; // One more than the last index
854:   }
855:   PetscInt          *clocal, aroots, aleaves, broots, bleaves;
856:   PetscSFNode       *cremote;
857:   const PetscInt    *alocal, *blocal;
858:   const PetscSFNode *aremote, *bremote;
859:   PetscCall(PetscMalloc2(maxleaf, &clocal, maxleaf, &cremote));
860:   for (PetscInt i = 0; i < maxleaf; i++) clocal[i] = -1;
861:   PetscCall(PetscSFGetGraph(sfa, &aroots, &aleaves, &alocal, &aremote));
862:   PetscCall(PetscSFGetGraph(sfb, &broots, &bleaves, &blocal, &bremote));
863:   PetscCheck(aroots == broots, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Both sfa and sfb must have the same root space");
864:   for (PetscInt i = 0; i < aleaves; i++) {
865:     PetscInt a = alocal ? alocal[i] : i;
866:     clocal[a]  = a;
867:     cremote[a] = aremote[i];
868:   }
869:   for (PetscInt i = 0; i < bleaves; i++) {
870:     PetscInt b = blocal ? blocal[i] : i;
871:     clocal[b]  = b;
872:     cremote[b] = bremote[i];
873:   }
874:   PetscInt nleaves = 0;
875:   for (PetscInt i = 0; i < maxleaf; i++) {
876:     if (clocal[i] < 0) continue;
877:     clocal[nleaves]  = clocal[i];
878:     cremote[nleaves] = cremote[i];
879:     nleaves++;
880:   }
881:   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfa), merged));
882:   PetscCall(PetscSFSetGraph(*merged, aroots, nleaves, clocal, PETSC_COPY_VALUES, cremote, PETSC_COPY_VALUES));
883:   PetscCall(PetscFree2(clocal, cremote));
884:   PetscFunctionReturn(PETSC_SUCCESS);
885: }

887: /*@
888:   PetscSFCreateStridedSF - Create an `PetscSF` to communicate interleaved blocks of data

890:   Collective

892:   Input Parameters:
893: + sf  - star forest
894: . bs  - stride
895: . ldr - leading dimension of root space
896: - ldl - leading dimension of leaf space

898:   Output Parameter:
899: . vsf - the new `PetscSF`

901:   Level: intermediate

903:   Notes:
904:   This can be useful to perform communications on multiple right-hand sides stored in a Fortran-style two dimensional array.
905:   For example, the calling sequence
906: .vb
907:   c_datatype *roots, *leaves;
908:   for i in [0,bs) do
909:     PetscSFBcastBegin(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op)
910:     PetscSFBcastEnd(sf, mpi_datatype, roots + i*ldr, leaves + i*ldl, op)
911: .ve
912:   is equivalent to
913: .vb
914:   c_datatype *roots, *leaves;
915:   PetscSFCreateStridedSF(sf, bs, ldr, ldl, &vsf)
916:   PetscSFBcastBegin(vsf, mpi_datatype, roots, leaves, op)
917:   PetscSFBcastEnd(vsf, mpi_datatype, roots, leaves, op)
918: .ve

920:   Developer Notes:
921:   Should this functionality be handled with a new API instead of creating a new object?

923: .seealso: [](sec_petscsf), `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()`
924: @*/
925: PetscErrorCode PetscSFCreateStridedSF(PetscSF sf, PetscInt bs, PetscInt ldr, PetscInt ldl, PetscSF *vsf)
926: {
927:   PetscSF            rankssf;
928:   const PetscSFNode *iremote, *sfrremote;
929:   PetscSFNode       *viremote;
930:   const PetscInt    *ilocal;
931:   PetscInt          *vilocal = NULL, *ldrs;
932:   PetscInt           nranks, nr, nl, vnr, vnl, maxl;
933:   PetscMPIInt        rank;
934:   MPI_Comm           comm;
935:   PetscSFType        sftype;

937:   PetscFunctionBegin;
940:   PetscAssertPointer(vsf, 5);
941:   if (bs == 1) {
942:     PetscCall(PetscObjectReference((PetscObject)sf));
943:     *vsf = sf;
944:     PetscFunctionReturn(PETSC_SUCCESS);
945:   }
946:   PetscCall(PetscSFSetUp(sf));
947:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
948:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
949:   PetscCall(PetscSFGetGraph(sf, &nr, &nl, &ilocal, &iremote));
950:   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxl));
951:   maxl += 1;
952:   if (ldl == PETSC_DECIDE) ldl = maxl;
953:   if (ldr == PETSC_DECIDE) ldr = nr;
954:   ldl /= PetscMax(1, sf->vscat.bs); // SFs created from VecScatterCreate() may have a nonzero block size. If not 0, we need to scale ldl and ldr
955:   ldr /= PetscMax(1, sf->vscat.bs);
956:   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);
957:   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);
958:   vnr = nr * bs;
959:   vnl = nl * bs;
960:   PetscCall(PetscMalloc1(vnl, &viremote));
961:   PetscCall(PetscMalloc1(vnl, &vilocal));

963:   /* Communicate root leading dimensions to leaf ranks */
964:   PetscCall(PetscSFGetRanksSF(sf, &rankssf));
965:   PetscCall(PetscSFGetGraph(rankssf, NULL, &nranks, NULL, &sfrremote));
966:   PetscCall(PetscMalloc1(nranks, &ldrs));
967:   PetscCall(PetscSFBcastBegin(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE));
968:   PetscCall(PetscSFBcastEnd(rankssf, MPIU_INT, &ldr, ldrs, MPI_REPLACE));

970:   for (PetscInt i = 0, rold = -1, lda = -1; i < nl; i++) {
971:     const PetscInt r  = iremote[i].rank;
972:     const PetscInt ii = iremote[i].index;

974:     if (r == rank) lda = ldr;
975:     else if (rold != r) {
976:       PetscInt j;

978:       for (j = 0; j < nranks; j++)
979:         if (sfrremote[j].rank == r) break;
980:       PetscCheck(j < nranks, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unable to locate neighbor rank %" PetscInt_FMT, r);
981:       lda = ldrs[j];
982:     }
983:     rold = r;
984:     for (PetscInt v = 0; v < bs; v++) {
985:       viremote[v * nl + i].rank  = r;
986:       viremote[v * nl + i].index = v * lda + ii;
987:       vilocal[v * nl + i]        = v * ldl + (ilocal ? ilocal[i] : i);
988:     }
989:   }
990:   PetscCall(PetscFree(ldrs));
991:   PetscCall(PetscSFCreate(comm, vsf));
992:   if (sf->vscat.bs > 1) {
993:     (*vsf)->vscat.bs = sf->vscat.bs;
994:     PetscCallMPI(MPI_Type_dup(sf->vscat.unit, &(*vsf)->vscat.unit));
995:     (*vsf)->vscat.to_n   = bs * sf->vscat.to_n;
996:     (*vsf)->vscat.from_n = bs * sf->vscat.from_n;
997:   }
998:   PetscCall(PetscSFGetType(sf, &sftype));
999:   PetscCall(PetscSFSetType(*vsf, sftype));
1000:   PetscCall(PetscSFSetGraph(*vsf, vnr, vnl, vilocal, PETSC_OWN_POINTER, viremote, PETSC_OWN_POINTER));
1001:   PetscFunctionReturn(PETSC_SUCCESS);
1002: }