Actual source code: plex.c

  1: #include <petsc/private/dmpleximpl.h>
  2: #include <petsc/private/dmlabelimpl.h>
  3: #include <petsc/private/isimpl.h>
  4: #include <petsc/private/vecimpl.h>
  5: #include <petsc/private/glvisvecimpl.h>
  6: #include <petscsf.h>
  7: #include <petscds.h>
  8: #include <petscdraw.h>
  9: #include <petscdmfield.h>
 10: #include <petscdmplextransform.h>
 11: #include <petscblaslapack.h>

 13: /* Logging support */
 14: PetscLogEvent DMPLEX_Interpolate, DMPLEX_Partition, DMPLEX_Distribute, DMPLEX_DistributeCones, DMPLEX_DistributeLabels, DMPLEX_DistributeSF, DMPLEX_DistributeOverlap, DMPLEX_DistributeField, DMPLEX_DistributeData, DMPLEX_Migrate, DMPLEX_InterpolateSF, DMPLEX_GlobalToNaturalBegin, DMPLEX_GlobalToNaturalEnd, DMPLEX_NaturalToGlobalBegin, DMPLEX_NaturalToGlobalEnd, DMPLEX_Stratify, DMPLEX_Symmetrize, DMPLEX_Preallocate, DMPLEX_ResidualFEM, DMPLEX_JacobianFEM, DMPLEX_InterpolatorFEM, DMPLEX_InjectorFEM, DMPLEX_IntegralFEM, DMPLEX_CreateGmsh, DMPLEX_CreateBoxSFC, DMPLEX_RebalanceSharedPoints, DMPLEX_PartSelf, DMPLEX_PartLabelInvert, DMPLEX_PartLabelCreateSF, DMPLEX_PartStratSF, DMPLEX_CreatePointSF, DMPLEX_LocatePoints, DMPLEX_TopologyView, DMPLEX_LabelsView, DMPLEX_CoordinatesView, DMPLEX_SectionView, DMPLEX_GlobalVectorView, DMPLEX_LocalVectorView, DMPLEX_TopologyLoad, DMPLEX_LabelsLoad, DMPLEX_CoordinatesLoad, DMPLEX_SectionLoad, DMPLEX_GlobalVectorLoad, DMPLEX_LocalVectorLoad;
 15: PetscLogEvent DMPLEX_RebalBuildGraph, DMPLEX_RebalRewriteSF, DMPLEX_RebalGatherGraph, DMPLEX_RebalPartition, DMPLEX_RebalScatterPart, DMPLEX_Generate, DMPLEX_Transform, DMPLEX_GetLocalOffsets, DMPLEX_Uninterpolate;

 17: PetscBool  Plexcite       = PETSC_FALSE;
 18: const char PlexCitation[] = "@article{LangeMitchellKnepleyGorman2015,\n"
 19:                             "title     = {Efficient mesh management in {Firedrake} using {PETSc-DMPlex}},\n"
 20:                             "author    = {Michael Lange and Lawrence Mitchell and Matthew G. Knepley and Gerard J. Gorman},\n"
 21:                             "journal   = {SIAM Journal on Scientific Computing},\n"
 22:                             "volume    = {38},\n"
 23:                             "number    = {5},\n"
 24:                             "pages     = {S143--S155},\n"
 25:                             "eprint    = {http://arxiv.org/abs/1506.07749},\n"
 26:                             "doi       = {10.1137/15M1026092},\n"
 27:                             "year      = {2016},\n"
 28:                             "petsc_uses={DMPlex},\n}\n";

 30: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);

 32: /*@
 33:   DMPlexIsSimplex - Is the first cell in this mesh a simplex?

 35:   Input Parameter:
 36: . dm - The `DMPLEX` object

 38:   Output Parameter:
 39: . simplex - Flag checking for a simplex

 41:   Level: intermediate

 43:   Note:
 44:   This just gives the first range of cells found. If the mesh has several cell types, it will only give the first.
 45:   If the mesh has no cells, this returns `PETSC_FALSE`.

 47: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetSimplexOrBoxCells()`, `DMPlexGetCellType()`, `DMPlexGetHeightStratum()`, `DMPolytopeTypeGetNumVertices()`
 48: @*/
 49: PetscErrorCode DMPlexIsSimplex(DM dm, PetscBool *simplex)
 50: {
 51:   DMPolytopeType ct;
 52:   PetscInt       cStart, cEnd;

 54:   PetscFunctionBegin;
 55:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
 56:   if (cEnd <= cStart) {
 57:     *simplex = PETSC_FALSE;
 58:     PetscFunctionReturn(PETSC_SUCCESS);
 59:   }
 60:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
 61:   *simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
 62:   PetscFunctionReturn(PETSC_SUCCESS);
 63: }

 65: /*@
 66:   DMPlexGetSimplexOrBoxCells - Get the range of cells which are neither prisms nor ghost FV cells

 68:   Input Parameters:
 69: + dm     - The `DMPLEX` object
 70: - height - The cell height in the Plex, 0 is the default

 72:   Output Parameters:
 73: + cStart - The first "normal" cell
 74: - cEnd   - The upper bound on "normal" cells

 76:   Level: developer

 78:   Note:
 79:   This function requires that tensor cells are ordered last.

 81: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexConstructGhostCells()`, `DMPlexGetCellTypeStratum()`
 82: @*/
 83: PetscErrorCode DMPlexGetSimplexOrBoxCells(DM dm, PetscInt height, PetscInt *cStart, PetscInt *cEnd)
 84: {
 85:   DMLabel         ctLabel;
 86:   IS              valueIS;
 87:   const PetscInt *ctypes;
 88:   PetscBool       found = PETSC_FALSE;
 89:   PetscInt        Nct, cS = PETSC_INT_MAX, cE = 0;

 91:   PetscFunctionBegin;
 92:   PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel));
 93:   PetscCall(DMLabelGetValueIS(ctLabel, &valueIS));
 94:   PetscCall(ISGetLocalSize(valueIS, &Nct));
 95:   PetscCall(ISGetIndices(valueIS, &ctypes));
 96:   for (PetscInt t = 0; t < Nct; ++t) {
 97:     const DMPolytopeType ct = (DMPolytopeType)ctypes[t];
 98:     PetscInt             ctS, ctE, ht;

100:     if (ct == DM_POLYTOPE_UNKNOWN) {
101:       // If any cells are not typed, just use all cells
102:       PetscCall(DMPlexGetHeightStratum(dm, PetscMax(height, 0), cStart, cEnd));
103:       break;
104:     }
105:     if (DMPolytopeTypeIsHybrid(ct) || ct == DM_POLYTOPE_FV_GHOST) continue;
106:     PetscCall(DMLabelGetStratumBounds(ctLabel, ct, &ctS, &ctE));
107:     if (ctS >= ctE) continue;
108:     // Check that a point has the right height
109:     PetscCall(DMPlexGetPointHeight(dm, ctS, &ht));
110:     if (ht != height) continue;
111:     cS    = PetscMin(cS, ctS);
112:     cE    = PetscMax(cE, ctE);
113:     found = PETSC_TRUE;
114:   }
115:   if (!Nct || !found) cS = cE = 0;
116:   PetscCall(ISDestroy(&valueIS));
117:   // Reset label for fast lookup
118:   PetscCall(DMLabelMakeAllInvalid_Internal(ctLabel));
119:   if (cStart) *cStart = cS;
120:   if (cEnd) *cEnd = cE;
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }

124: PetscErrorCode DMPlexGetFieldTypes_Internal(DM dm, PetscSection section, PetscInt field, PetscInt *types, PetscInt **ssStart, PetscInt **ssEnd, PetscViewerVTKFieldType **sft)
125: {
126:   PetscInt                 cdim, pStart, pEnd, vStart, vEnd, cStart, cEnd, c, depth, cellHeight, t;
127:   PetscInt                *sStart, *sEnd;
128:   PetscViewerVTKFieldType *ft;
129:   PetscInt                 vcdof[DM_NUM_POLYTOPES + 1], globalvcdof[DM_NUM_POLYTOPES + 1];
130:   DMLabel                  depthLabel, ctLabel;

132:   PetscFunctionBegin;
133:   /* the vcdof and globalvcdof are sized to allow every polytope type and simple vertex at DM_NUM_POLYTOPES */
134:   PetscCall(PetscArrayzero(vcdof, DM_NUM_POLYTOPES + 1));
135:   PetscCall(DMGetCoordinateDim(dm, &cdim));
136:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
137:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
138:   if (field >= 0) {
139:     if ((vStart >= pStart) && (vStart < pEnd)) PetscCall(PetscSectionGetFieldDof(section, vStart, field, &vcdof[DM_NUM_POLYTOPES]));
140:   } else {
141:     if ((vStart >= pStart) && (vStart < pEnd)) PetscCall(PetscSectionGetDof(section, vStart, &vcdof[DM_NUM_POLYTOPES]));
142:   }

144:   PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
145:   PetscCall(DMPlexGetDepth(dm, &depth));
146:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
147:   PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel));
148:   for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
149:     const DMPolytopeType ict = (DMPolytopeType)c;
150:     PetscInt             dep;

152:     if (ict == DM_POLYTOPE_FV_GHOST) continue;
153:     PetscCall(DMLabelGetStratumBounds(ctLabel, ict, &cStart, &cEnd));
154:     if (pStart >= 0) {
155:       PetscCall(DMLabelGetValue(depthLabel, cStart, &dep));
156:       if (dep != depth - cellHeight) continue;
157:     }
158:     if (field >= 0) {
159:       if ((cStart >= pStart) && (cStart < pEnd)) PetscCall(PetscSectionGetFieldDof(section, cStart, field, &vcdof[c]));
160:     } else {
161:       if ((cStart >= pStart) && (cStart < pEnd)) PetscCall(PetscSectionGetDof(section, cStart, &vcdof[c]));
162:     }
163:   }

165:   PetscCallMPI(MPIU_Allreduce(vcdof, globalvcdof, DM_NUM_POLYTOPES + 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
166:   *types = 0;

168:   for (c = 0; c < DM_NUM_POLYTOPES + 1; ++c) {
169:     if (globalvcdof[c]) ++(*types);
170:   }

172:   PetscCall(PetscMalloc3(*types, &sStart, *types, &sEnd, *types, &ft));
173:   t = 0;
174:   if (globalvcdof[DM_NUM_POLYTOPES]) {
175:     sStart[t] = vStart;
176:     sEnd[t]   = vEnd;
177:     ft[t]     = (globalvcdof[t] == cdim) ? PETSC_VTK_POINT_VECTOR_FIELD : PETSC_VTK_POINT_FIELD;
178:     ++t;
179:   }

181:   for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
182:     if (globalvcdof[c]) {
183:       const DMPolytopeType ict = (DMPolytopeType)c;

185:       PetscCall(DMLabelGetStratumBounds(ctLabel, ict, &cStart, &cEnd));
186:       sStart[t] = cStart;
187:       sEnd[t]   = cEnd;
188:       ft[t]     = (globalvcdof[c] == cdim) ? PETSC_VTK_CELL_VECTOR_FIELD : PETSC_VTK_CELL_FIELD;
189:       ++t;
190:     }
191:   }

193:   if (!*types) {
194:     if (field >= 0) {
195:       const char *fieldname;

197:       PetscCall(PetscSectionGetFieldName(section, field, &fieldname));
198:       PetscCall(PetscInfo((PetscObject)dm, "Could not classify VTK output type of section field %" PetscInt_FMT " \"%s\"\n", field, fieldname));
199:     } else {
200:       PetscCall(PetscInfo((PetscObject)dm, "Could not classify VTK output type of section\n"));
201:     }
202:   }

204:   *ssStart = sStart;
205:   *ssEnd   = sEnd;
206:   *sft     = ft;
207:   PetscFunctionReturn(PETSC_SUCCESS);
208: }

210: PetscErrorCode DMPlexRestoreFieldTypes_Internal(DM dm, PetscSection section, PetscInt field, PetscInt *types, PetscInt **sStart, PetscInt **sEnd, PetscViewerVTKFieldType **ft)
211: {
212:   PetscFunctionBegin;
213:   PetscCall(PetscFree3(*sStart, *sEnd, *ft));
214:   PetscFunctionReturn(PETSC_SUCCESS);
215: }

217: PetscErrorCode DMPlexGetFieldType_Internal(DM dm, PetscSection section, PetscInt field, PetscInt *sStart, PetscInt *sEnd, PetscViewerVTKFieldType *ft)
218: {
219:   PetscInt cdim, pStart, pEnd, vStart, vEnd, cStart, cEnd;
220:   PetscInt vcdof[2] = {0, 0}, globalvcdof[2];

222:   PetscFunctionBegin;
223:   *ft = PETSC_VTK_INVALID;
224:   PetscCall(DMGetCoordinateDim(dm, &cdim));
225:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
226:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
227:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
228:   if (field >= 0) {
229:     if ((vStart >= pStart) && (vStart < pEnd)) PetscCall(PetscSectionGetFieldDof(section, vStart, field, &vcdof[0]));
230:     if ((cStart >= pStart) && (cStart < pEnd)) PetscCall(PetscSectionGetFieldDof(section, cStart, field, &vcdof[1]));
231:   } else {
232:     if ((vStart >= pStart) && (vStart < pEnd)) PetscCall(PetscSectionGetDof(section, vStart, &vcdof[0]));
233:     if ((cStart >= pStart) && (cStart < pEnd)) PetscCall(PetscSectionGetDof(section, cStart, &vcdof[1]));
234:   }
235:   PetscCallMPI(MPIU_Allreduce(vcdof, globalvcdof, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
236:   if (globalvcdof[0]) {
237:     *sStart = vStart;
238:     *sEnd   = vEnd;
239:     if (globalvcdof[0] == cdim) *ft = PETSC_VTK_POINT_VECTOR_FIELD;
240:     else *ft = PETSC_VTK_POINT_FIELD;
241:   } else if (globalvcdof[1]) {
242:     *sStart = cStart;
243:     *sEnd   = cEnd;
244:     if (globalvcdof[1] == cdim) *ft = PETSC_VTK_CELL_VECTOR_FIELD;
245:     else *ft = PETSC_VTK_CELL_FIELD;
246:   } else {
247:     if (field >= 0) {
248:       const char *fieldname;

250:       PetscCall(PetscSectionGetFieldName(section, field, &fieldname));
251:       PetscCall(PetscInfo(dm, "Could not classify VTK output type of section field %" PetscInt_FMT " \"%s\"\n", field, fieldname));
252:     } else {
253:       PetscCall(PetscInfo(dm, "Could not classify VTK output type of section\n"));
254:     }
255:   }
256:   PetscFunctionReturn(PETSC_SUCCESS);
257: }

259: /*@
260:   DMPlexVecView1D - Plot many 1D solutions on the same line graph

262:   Collective

264:   Input Parameters:
265: + dm     - The `DMPLEX` object
266: . n      - The number of vectors
267: . u      - The array of local vectors
268: - viewer - The `PetscViewer`

270:   Level: advanced

272: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `VecViewFromOptions()`, `VecView()`
273: @*/
274: PetscErrorCode DMPlexVecView1D(DM dm, PetscInt n, Vec u[], PetscViewer viewer)
275: {
276:   DM                 cdm;
277:   PetscDS            ds;
278:   PetscDraw          draw = NULL;
279:   PetscDrawLG        lg;
280:   Vec                coordinates;
281:   const PetscScalar *coords, **sol;
282:   PetscReal         *vals;
283:   PetscInt          *Nc;
284:   PetscInt           Nf, Nl, vStart, vEnd, eStart, eEnd;
285:   char             **names;

287:   PetscFunctionBegin;
288:   PetscCall(DMGetCoordinateDM(dm, &cdm));
289:   PetscCall(DMGetDS(dm, &ds));
290:   PetscCall(PetscDSGetNumFields(ds, &Nf));
291:   PetscCall(PetscDSGetTotalComponents(ds, &Nl));
292:   PetscCall(PetscDSGetComponents(ds, &Nc));

294:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
295:   if (!draw) PetscFunctionReturn(PETSC_SUCCESS);
296:   PetscCall(PetscDrawLGCreate(draw, n * Nl, &lg));

298:   PetscCall(PetscMalloc3(n, &sol, n * Nl, &names, n * Nl, &vals));
299:   for (PetscInt i = 0, l = 0; i < n; ++i) {
300:     const char *vname;

302:     PetscCall(PetscObjectGetName((PetscObject)u[i], &vname));
303:     for (PetscInt f = 0; f < Nf; ++f) {
304:       PetscObject disc;
305:       const char *fname;
306:       char        tmpname[PETSC_MAX_PATH_LEN];

308:       PetscCall(PetscDSGetDiscretization(ds, f, &disc));
309:       /* TODO Create names for components */
310:       for (PetscInt c = 0; c < Nc[f]; ++c, ++l) {
311:         PetscCall(PetscObjectGetName(disc, &fname));
312:         PetscCall(PetscStrncpy(tmpname, vname, sizeof(tmpname)));
313:         PetscCall(PetscStrlcat(tmpname, ":", sizeof(tmpname)));
314:         PetscCall(PetscStrlcat(tmpname, fname, sizeof(tmpname)));
315:         PetscCall(PetscStrallocpy(tmpname, &names[l]));
316:       }
317:     }
318:   }
319:   PetscCall(PetscDrawLGSetLegend(lg, (const char *const *)names));
320:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
321:   PetscCall(VecGetArrayRead(coordinates, &coords));
322:   for (PetscInt i = 0; i < n; ++i) PetscCall(VecGetArrayRead(u[i], &sol[i]));
323:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
324:   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
325:   PetscSection s;
326:   PetscInt     cdof, vdof;

328:   PetscCall(DMGetLocalSection(dm, &s));
329:   PetscCall(PetscSectionGetDof(s, eStart, &cdof));
330:   PetscCall(PetscSectionGetDof(s, vStart, &vdof));
331:   if (cdof) {
332:     if (vdof) {
333:       // P_2
334:       PetscInt vFirst = -1;

336:       for (PetscInt e = eStart; e < eEnd; ++e) {
337:         PetscScalar    *xa, *xb, *svals;
338:         const PetscInt *cone;

340:         PetscCall(DMPlexGetCone(dm, e, &cone));
341:         PetscCall(DMPlexPointLocalRead(cdm, cone[0], coords, &xa));
342:         PetscCall(DMPlexPointLocalRead(cdm, cone[1], coords, &xb));
343:         if (e == eStart) vFirst = cone[0];
344:         for (PetscInt i = 0; i < n; ++i) {
345:           PetscCall(DMPlexPointLocalRead(dm, cone[0], sol[i], &svals));
346:           for (PetscInt l = 0; l < Nl; ++l) vals[i * Nl + l] = PetscRealPart(svals[l]);
347:         }
348:         PetscCall(PetscDrawLGAddCommonPoint(lg, PetscRealPart(xa[0]), vals));
349:         if (e == eEnd - 1 && cone[1] != vFirst) {
350:           for (PetscInt i = 0; i < n; ++i) {
351:             PetscCall(DMPlexPointLocalRead(dm, e, sol[i], &svals));
352:             for (PetscInt l = 0; l < Nl; ++l) vals[i * Nl + l] = PetscRealPart(svals[l]);
353:           }
354:           PetscCall(PetscDrawLGAddCommonPoint(lg, 0.5 * (PetscRealPart(xa[0]) + PetscRealPart(xb[0])), vals));
355:           for (PetscInt i = 0; i < n; ++i) {
356:             PetscCall(DMPlexPointLocalRead(dm, cone[1], sol[i], &svals));
357:             for (PetscInt l = 0; l < Nl; ++l) vals[i * Nl + l] = PetscRealPart(svals[l]);
358:           }
359:           PetscCall(PetscDrawLGAddCommonPoint(lg, PetscRealPart(xb[0]), vals));
360:         }
361:       }
362:     } else {
363:       // P_0
364:       for (PetscInt e = eStart; e < eEnd; ++e) {
365:         PetscScalar    *xa, *xb, *svals;
366:         const PetscInt *cone;

368:         PetscCall(DMPlexGetCone(dm, e, &cone));
369:         PetscCall(DMPlexPointLocalRead(cdm, cone[0], coords, &xa));
370:         PetscCall(DMPlexPointLocalRead(cdm, cone[1], coords, &xb));
371:         for (PetscInt i = 0; i < n; ++i) {
372:           PetscCall(DMPlexPointLocalRead(dm, e, sol[i], &svals));
373:           for (PetscInt l = 0; l < Nl; ++l) vals[i * Nl + l] = PetscRealPart(svals[l]);
374:         }
375:         PetscCall(PetscDrawLGAddCommonPoint(lg, 0.5 * (PetscRealPart(xa[0]) + PetscRealPart(xb[0])), vals));
376:       }
377:     }
378:   } else if (vdof) {
379:     // P_1
380:     for (PetscInt v = vStart; v < vEnd; ++v) {
381:       PetscScalar *x, *svals;

383:       PetscCall(DMPlexPointLocalRead(cdm, v, coords, &x));
384:       for (PetscInt i = 0; i < n; ++i) {
385:         PetscCall(DMPlexPointLocalRead(dm, v, sol[i], &svals));
386:         for (PetscInt l = 0; l < Nl; ++l) vals[i * Nl + l] = PetscRealPart(svals[l]);
387:       }
388:       PetscCall(PetscDrawLGAddCommonPoint(lg, PetscRealPart(x[0]), vals));
389:     }
390:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Discretization not supported");
391:   PetscCall(VecRestoreArrayRead(coordinates, &coords));
392:   for (PetscInt i = 0; i < n; ++i) PetscCall(VecRestoreArrayRead(u[i], &sol[i]));
393:   for (PetscInt l = 0; l < n * Nl; ++l) PetscCall(PetscFree(names[l]));
394:   PetscCall(PetscFree3(sol, names, vals));

396:   PetscCall(PetscDrawLGDraw(lg));
397:   PetscCall(PetscDrawLGDestroy(&lg));
398:   PetscFunctionReturn(PETSC_SUCCESS);
399: }

401: static PetscErrorCode VecView_Plex_Local_Draw_1D(Vec u, PetscViewer viewer)
402: {
403:   DM dm;

405:   PetscFunctionBegin;
406:   PetscCall(VecGetDM(u, &dm));
407:   PetscCall(DMPlexVecView1D(dm, 1, &u, viewer));
408:   PetscFunctionReturn(PETSC_SUCCESS);
409: }

411: static PetscErrorCode VecView_Plex_Local_Draw_2D(Vec v, PetscViewer viewer)
412: {
413:   DM                 dm;
414:   PetscSection       s;
415:   PetscDraw          draw, popup;
416:   DM                 cdm;
417:   PetscSection       coordSection;
418:   Vec                coordinates;
419:   const PetscScalar *array;
420:   PetscReal          lbound[3], ubound[3];
421:   PetscReal          vbound[2], time;
422:   PetscBool          flg;
423:   PetscInt           dim, Nf, f, Nc, comp, vStart, vEnd, cStart, cEnd, c, N, level, step, w = 0;
424:   const char        *name;
425:   char               title[PETSC_MAX_PATH_LEN];

427:   PetscFunctionBegin;
428:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
429:   PetscCall(VecGetDM(v, &dm));
430:   PetscCall(DMGetCoordinateDim(dm, &dim));
431:   PetscCall(DMGetLocalSection(dm, &s));
432:   PetscCall(PetscSectionGetNumFields(s, &Nf));
433:   PetscCall(DMGetCoarsenLevel(dm, &level));
434:   PetscCall(DMGetCoordinateDM(dm, &cdm));
435:   PetscCall(DMGetLocalSection(cdm, &coordSection));
436:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
437:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
438:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));

440:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
441:   PetscCall(DMGetOutputSequenceNumber(dm, &step, &time));

443:   PetscCall(VecGetLocalSize(coordinates, &N));
444:   PetscCall(DMGetBoundingBox(dm, lbound, ubound));
445:   PetscCall(PetscDrawClear(draw));

447:   /* Could implement something like DMDASelectFields() */
448:   for (f = 0; f < Nf; ++f) {
449:     DM          fdm = dm;
450:     Vec         fv  = v;
451:     IS          fis;
452:     char        prefix[PETSC_MAX_PATH_LEN];
453:     const char *fname;

455:     PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
456:     PetscCall(PetscSectionGetFieldName(s, f, &fname));

458:     if (v->hdr.prefix) PetscCall(PetscStrncpy(prefix, v->hdr.prefix, sizeof(prefix)));
459:     else prefix[0] = '\0';
460:     if (Nf > 1) {
461:       PetscCall(DMCreateSubDM(dm, 1, &f, &fis, &fdm));
462:       PetscCall(VecGetSubVector(v, fis, &fv));
463:       PetscCall(PetscStrlcat(prefix, fname, sizeof(prefix)));
464:       PetscCall(PetscStrlcat(prefix, "_", sizeof(prefix)));
465:     }
466:     for (comp = 0; comp < Nc; ++comp, ++w) {
467:       PetscInt nmax = 2;

469:       PetscCall(PetscViewerDrawGetDraw(viewer, w, &draw));
470:       if (Nc > 1) PetscCall(PetscSNPrintf(title, sizeof(title), "%s:%s_%" PetscInt_FMT " Step: %" PetscInt_FMT " Time: %.4g", name, fname, comp, step, (double)time));
471:       else PetscCall(PetscSNPrintf(title, sizeof(title), "%s:%s Step: %" PetscInt_FMT " Time: %.4g", name, fname, step, (double)time));
472:       PetscCall(PetscDrawSetTitle(draw, title));

474:       /* TODO Get max and min only for this component */
475:       PetscCall(PetscOptionsGetRealArray(NULL, prefix, "-vec_view_bounds", vbound, &nmax, &flg));
476:       if (!flg) {
477:         PetscCall(VecMin(fv, NULL, &vbound[0]));
478:         PetscCall(VecMax(fv, NULL, &vbound[1]));
479:         if (vbound[1] <= vbound[0]) vbound[1] = vbound[0] + 1.0;
480:       }

482:       PetscCall(PetscDrawGetPopup(draw, &popup));
483:       PetscCall(PetscDrawScalePopup(popup, vbound[0], vbound[1]));
484:       PetscCall(PetscDrawSetCoordinates(draw, lbound[0], lbound[1], ubound[0], ubound[1]));
485:       PetscCall(VecGetArrayRead(fv, &array));
486:       for (c = cStart; c < cEnd; ++c) {
487:         DMPolytopeType     ct;
488:         PetscScalar       *coords = NULL, *a = NULL;
489:         const PetscScalar *coords_arr;
490:         PetscBool          isDG;
491:         PetscInt           numCoords;
492:         int                color[4] = {-1, -1, -1, -1};

494:         PetscCall(DMPlexGetCellType(dm, c, &ct));
495:         PetscCall(DMPlexPointLocalRead(fdm, c, array, &a));
496:         if (a) {
497:           color[0] = PetscDrawRealToColor(PetscRealPart(a[comp]), vbound[0], vbound[1]);
498:           color[1] = color[2] = color[3] = color[0];
499:         } else {
500:           PetscScalar *vals = NULL;
501:           PetscInt     numVals, va;

503:           PetscCall(DMPlexVecGetClosure(fdm, NULL, fv, c, &numVals, &vals));
504:           if (!numVals) {
505:             PetscCall(DMPlexVecRestoreClosure(fdm, NULL, fv, c, &numVals, &vals));
506:             continue;
507:           }
508:           PetscCheck(numVals % Nc == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of components %" PetscInt_FMT " does not divide the number of values in the closure %" PetscInt_FMT, Nc, numVals);
509:           switch (numVals / Nc) {
510:           case 1: /* P1 Clamped Segment Prism */
511:           case 2: /* P1 Segment Prism, P2 Clamped Segment Prism */
512:             PetscCheck(ct == DM_POLYTOPE_SEG_PRISM_TENSOR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell should be a tensor segment, but it is a %s", DMPolytopeTypes[ct]);
513:             for (va = 0; va < numVals / Nc; ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va * Nc + comp]), vbound[0], vbound[1]);
514:             break;
515:           case 3: /* P1 Triangle */
516:           case 4: /* P1 Quadrangle */
517:             PetscCheck(ct == DM_POLYTOPE_TRIANGLE || ct == DM_POLYTOPE_QUADRILATERAL || ct == DM_POLYTOPE_SEG_PRISM_TENSOR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell should be a triangle or quad, but it is a %s", DMPolytopeTypes[ct]);
518:             for (va = 0; va < numVals / Nc; ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va * Nc + comp]), vbound[0], vbound[1]);
519:             break;
520:           case 6: /* P2 Triangle */
521:           case 8: /* P2 Quadrangle */
522:             PetscCheck(ct == DM_POLYTOPE_TRIANGLE || ct == DM_POLYTOPE_QUADRILATERAL || ct == DM_POLYTOPE_SEG_PRISM_TENSOR, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell should be a triangle or quad, but it is a %s", DMPolytopeTypes[ct]);
523:             for (va = 0; va < numVals / (Nc * 2); ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va * Nc + comp + numVals / (Nc * 2)]), vbound[0], vbound[1]);
524:             break;
525:           default:
526:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of values for cell closure %" PetscInt_FMT " cannot be handled", numVals / Nc);
527:           }
528:           PetscCall(DMPlexVecRestoreClosure(fdm, NULL, fv, c, &numVals, &vals));
529:         }
530:         PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
531:         switch (numCoords) {
532:         case 6:
533:         case 12: /* Localized triangle */
534:           PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), color[0], color[1], color[2]));
535:           break;
536:         case 8:
537:         case 16: /* Localized quadrilateral */
538:           if (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) {
539:             PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscMax(color[0], color[1])));
540:           } else {
541:             PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), color[0], color[1], color[2]));
542:             PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), color[2], color[3], color[0]));
543:           }
544:           break;
545:         default:
546:           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells with %" PetscInt_FMT " coordinates", numCoords);
547:         }
548:         PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
549:       }
550:       PetscCall(VecRestoreArrayRead(fv, &array));
551:       PetscCall(PetscDrawFlush(draw));
552:       PetscCall(PetscDrawPause(draw));
553:       PetscCall(PetscDrawSave(draw));
554:     }
555:     if (Nf > 1) {
556:       PetscCall(VecRestoreSubVector(v, fis, &fv));
557:       PetscCall(ISDestroy(&fis));
558:       PetscCall(DMDestroy(&fdm));
559:     }
560:   }
561:   PetscFunctionReturn(PETSC_SUCCESS);
562: }

564: static PetscErrorCode VecView_Plex_Local_Draw(Vec v, PetscViewer viewer)
565: {
566:   DM        dm;
567:   PetscDraw draw;
568:   PetscInt  dim;
569:   PetscBool isnull;

571:   PetscFunctionBegin;
572:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
573:   PetscCall(PetscDrawIsNull(draw, &isnull));
574:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

576:   PetscCall(VecGetDM(v, &dm));
577:   PetscCall(DMGetCoordinateDim(dm, &dim));
578:   switch (dim) {
579:   case 1:
580:     PetscCall(VecView_Plex_Local_Draw_1D(v, viewer));
581:     break;
582:   case 2:
583:     PetscCall(VecView_Plex_Local_Draw_2D(v, viewer));
584:     break;
585:   default:
586:     SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot draw meshes of dimension %" PetscInt_FMT ". Try PETSCVIEWERGLVIS", dim);
587:   }
588:   PetscFunctionReturn(PETSC_SUCCESS);
589: }

591: static PetscErrorCode VecView_Plex_Local_VTK(Vec v, PetscViewer viewer)
592: {
593:   DM                      dm;
594:   Vec                     locv;
595:   const char             *name;
596:   PetscSection            section;
597:   PetscInt                pStart, pEnd;
598:   PetscInt                numFields;
599:   PetscViewerVTKFieldType ft;

601:   PetscFunctionBegin;
602:   PetscCall(VecGetDM(v, &dm));
603:   PetscCall(DMCreateLocalVector(dm, &locv)); /* VTK viewer requires exclusive ownership of the vector */
604:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
605:   PetscCall(PetscObjectSetName((PetscObject)locv, name));
606:   PetscCall(VecCopy(v, locv));
607:   PetscCall(DMGetLocalSection(dm, &section));
608:   PetscCall(PetscSectionGetNumFields(section, &numFields));
609:   if (!numFields) {
610:     PetscCall(DMPlexGetFieldType_Internal(dm, section, PETSC_DETERMINE, &pStart, &pEnd, &ft));
611:     PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)dm, DMPlexVTKWriteAll, PETSC_DEFAULT, ft, PETSC_TRUE, (PetscObject)locv));
612:   } else {
613:     PetscInt f;

615:     for (f = 0; f < numFields; f++) {
616:       PetscCall(DMPlexGetFieldType_Internal(dm, section, f, &pStart, &pEnd, &ft));
617:       if (ft == PETSC_VTK_INVALID) continue;
618:       PetscCall(PetscObjectReference((PetscObject)locv));
619:       PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)dm, DMPlexVTKWriteAll, f, ft, PETSC_TRUE, (PetscObject)locv));
620:     }
621:     PetscCall(VecDestroy(&locv));
622:   }
623:   PetscFunctionReturn(PETSC_SUCCESS);
624: }

626: PetscErrorCode VecView_Plex_Local(Vec v, PetscViewer viewer)
627: {
628:   DM        dm;
629:   PetscBool isvtk, ishdf5, isdraw, isglvis, iscgns, ispython;

631:   PetscFunctionBegin;
632:   PetscCall(VecGetDM(v, &dm));
633:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
634:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
635:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
636:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
637:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
638:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
639:   PetscCall(PetscObjectHasFunction((PetscObject)viewer, "PetscViewerPythonViewObject_C", &ispython));
640:   if (isvtk || ishdf5 || isdraw || isglvis || iscgns || ispython) {
641:     PetscInt    i, numFields;
642:     PetscObject fe;
643:     PetscBool   fem  = PETSC_FALSE;
644:     Vec         locv = v;
645:     const char *name;
646:     PetscInt    step;
647:     PetscReal   time;

649:     PetscCall(DMGetNumFields(dm, &numFields));
650:     for (i = 0; i < numFields; i++) {
651:       PetscCall(DMGetField(dm, i, NULL, &fe));
652:       if (fe->classid == PETSCFE_CLASSID) {
653:         fem = PETSC_TRUE;
654:         break;
655:       }
656:     }
657:     if (fem) {
658:       PetscObject isZero;

660:       PetscCall(DMGetLocalVector(dm, &locv));
661:       PetscCall(PetscObjectGetName((PetscObject)v, &name));
662:       PetscCall(PetscObjectSetName((PetscObject)locv, name));
663:       PetscCall(PetscObjectQuery((PetscObject)v, "__Vec_bc_zero__", &isZero));
664:       PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", isZero));
665:       PetscCall(VecCopy(v, locv));
666:       PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
667:       PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL));
668:     }
669:     if (isvtk) {
670:       PetscCall(VecView_Plex_Local_VTK(locv, viewer));
671:     } else if (ishdf5) {
672: #if defined(PETSC_HAVE_HDF5)
673:       PetscCall(VecView_Plex_Local_HDF5_Internal(locv, viewer));
674: #else
675:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
676: #endif
677:     } else if (isdraw) {
678:       PetscCall(VecView_Plex_Local_Draw(locv, viewer));
679:     } else if (ispython) {
680:       PetscCall(PetscViewerPythonViewObject(viewer, (PetscObject)locv));
681:     } else if (isglvis) {
682:       PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
683:       PetscCall(PetscViewerGLVisSetSnapId(viewer, step));
684:       PetscCall(VecView_GLVis(locv, viewer));
685:     } else if (iscgns) {
686: #if defined(PETSC_HAVE_CGNS)
687:       PetscCall(VecView_Plex_Local_CGNS(locv, viewer));
688: #else
689:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CGNS not supported in this build.\nPlease reconfigure using --download-cgns");
690: #endif
691:     }
692:     if (fem) {
693:       PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", NULL));
694:       PetscCall(DMRestoreLocalVector(dm, &locv));
695:     }
696:   } else {
697:     PetscBool isseq;

699:     PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
700:     if (isseq) PetscCall(VecView_Seq(v, viewer));
701:     else PetscCall(VecView_MPI(v, viewer));
702:   }
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: PetscErrorCode VecView_Plex(Vec v, PetscViewer viewer)
707: {
708:   DM        dm;
709:   PetscBool isvtk, ishdf5, isdraw, isglvis, isexodusii, iscgns, ispython;

711:   PetscFunctionBegin;
712:   PetscCall(VecGetDM(v, &dm));
713:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
714:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
715:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
716:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
717:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
718:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
719:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWEREXODUSII, &isexodusii));
720:   PetscCall(PetscObjectHasFunction((PetscObject)viewer, "PetscViewerPythonViewObject_C", &ispython));
721:   if (isvtk || isdraw || isglvis || iscgns || ispython) {
722:     Vec         locv;
723:     PetscObject isZero;
724:     const char *name;

726:     PetscCall(DMGetLocalVector(dm, &locv));
727:     PetscCall(PetscObjectGetName((PetscObject)v, &name));
728:     PetscCall(PetscObjectSetName((PetscObject)locv, name));
729:     PetscCall(DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv));
730:     PetscCall(DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv));
731:     PetscCall(PetscObjectQuery((PetscObject)v, "__Vec_bc_zero__", &isZero));
732:     PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", isZero));
733:     PetscCall(VecView_Plex_Local(locv, viewer));
734:     PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", NULL));
735:     PetscCall(DMRestoreLocalVector(dm, &locv));
736:   } else if (ishdf5) {
737: #if defined(PETSC_HAVE_HDF5)
738:     PetscCall(VecView_Plex_HDF5_Internal(v, viewer));
739: #else
740:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
741: #endif
742:   } else if (isexodusii) {
743: #if defined(PETSC_HAVE_EXODUSII)
744:     PetscCall(VecView_PlexExodusII_Internal(v, viewer));
745: #else
746:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "ExodusII not supported in this build.\nPlease reconfigure using --download-exodusii");
747: #endif
748:   } else {
749:     PetscBool isseq;

751:     PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
752:     if (isseq) PetscCall(VecView_Seq(v, viewer));
753:     else PetscCall(VecView_MPI(v, viewer));
754:   }
755:   PetscFunctionReturn(PETSC_SUCCESS);
756: }

758: PetscErrorCode VecView_Plex_Native(Vec originalv, PetscViewer viewer)
759: {
760:   DM                dm;
761:   MPI_Comm          comm;
762:   PetscViewerFormat format;
763:   Vec               v;
764:   PetscBool         isvtk, ishdf5;

766:   PetscFunctionBegin;
767:   PetscCall(VecGetDM(originalv, &dm));
768:   PetscCall(PetscObjectGetComm((PetscObject)originalv, &comm));
769:   PetscCheck(dm, comm, PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
770:   PetscCall(PetscViewerGetFormat(viewer, &format));
771:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
772:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
773:   if (format == PETSC_VIEWER_NATIVE) {
774:     /* Natural ordering is the common case for DMDA, NATIVE means plain vector, for PLEX is the opposite */
775:     /* this need a better fix */
776:     if (dm->useNatural) {
777:       if (dm->sfNatural) {
778:         const char *vecname;
779:         PetscInt    n, nroots;

781:         PetscCall(VecGetLocalSize(originalv, &n));
782:         PetscCall(PetscSFGetGraph(dm->sfNatural, &nroots, NULL, NULL, NULL));
783:         if (n == nroots) {
784:           PetscCall(DMPlexCreateNaturalVector(dm, &v));
785:           PetscCall(DMPlexGlobalToNaturalBegin(dm, originalv, v));
786:           PetscCall(DMPlexGlobalToNaturalEnd(dm, originalv, v));
787:           PetscCall(PetscObjectGetName((PetscObject)originalv, &vecname));
788:           PetscCall(PetscObjectSetName((PetscObject)v, vecname));
789:         } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "DM global to natural SF only handles global vectors");
790:       } else SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created");
791:     } else v = originalv;
792:   } else v = originalv;

794:   if (ishdf5) {
795: #if defined(PETSC_HAVE_HDF5)
796:     PetscCall(VecView_Plex_HDF5_Native_Internal(v, viewer));
797: #else
798:     SETERRQ(comm, PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
799: #endif
800:   } else if (isvtk) {
801:     SETERRQ(comm, PETSC_ERR_SUP, "VTK format does not support viewing in natural order. Please switch to HDF5.");
802:   } else {
803:     PetscBool isseq;

805:     PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
806:     if (isseq) PetscCall(VecView_Seq(v, viewer));
807:     else PetscCall(VecView_MPI(v, viewer));
808:   }
809:   if (v != originalv) PetscCall(VecDestroy(&v));
810:   PetscFunctionReturn(PETSC_SUCCESS);
811: }

813: PetscErrorCode VecLoad_Plex_Local(Vec v, PetscViewer viewer)
814: {
815:   DM        dm;
816:   PetscBool ishdf5;

818:   PetscFunctionBegin;
819:   PetscCall(VecGetDM(v, &dm));
820:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
821:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
822:   if (ishdf5) {
823:     DM          dmBC;
824:     Vec         gv;
825:     const char *name;

827:     PetscCall(DMGetOutputDM(dm, &dmBC));
828:     PetscCall(DMGetGlobalVector(dmBC, &gv));
829:     PetscCall(PetscObjectGetName((PetscObject)v, &name));
830:     PetscCall(PetscObjectSetName((PetscObject)gv, name));
831:     PetscCall(VecLoad_Default(gv, viewer));
832:     PetscCall(DMGlobalToLocalBegin(dmBC, gv, INSERT_VALUES, v));
833:     PetscCall(DMGlobalToLocalEnd(dmBC, gv, INSERT_VALUES, v));
834:     PetscCall(DMRestoreGlobalVector(dmBC, &gv));
835:   } else PetscCall(VecLoad_Default(v, viewer));
836:   PetscFunctionReturn(PETSC_SUCCESS);
837: }

839: PetscErrorCode VecLoad_Plex(Vec v, PetscViewer viewer)
840: {
841:   DM        dm;
842:   PetscBool ishdf5, isexodusii, iscgns;

844:   PetscFunctionBegin;
845:   PetscCall(VecGetDM(v, &dm));
846:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
847:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
848:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWEREXODUSII, &isexodusii));
849:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
850:   if (ishdf5) {
851: #if defined(PETSC_HAVE_HDF5)
852:     PetscCall(VecLoad_Plex_HDF5_Internal(v, viewer));
853: #else
854:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
855: #endif
856:   } else if (isexodusii) {
857: #if defined(PETSC_HAVE_EXODUSII)
858:     PetscCall(VecLoad_PlexExodusII_Internal(v, viewer));
859: #else
860:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "ExodusII not supported in this build.\nPlease reconfigure using --download-exodusii");
861: #endif
862:   } else if (iscgns) {
863: #if defined(PETSC_HAVE_CGNS)
864:     PetscCall(VecLoad_Plex_CGNS_Internal(v, viewer));
865: #else
866:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CGNS not supported in this build.\nPlease reconfigure using --download-cgns");
867: #endif
868:   } else PetscCall(VecLoad_Default(v, viewer));
869:   PetscFunctionReturn(PETSC_SUCCESS);
870: }

872: PetscErrorCode VecLoad_Plex_Native(Vec originalv, PetscViewer viewer)
873: {
874:   DM                dm;
875:   PetscViewerFormat format;
876:   PetscBool         ishdf5;

878:   PetscFunctionBegin;
879:   PetscCall(VecGetDM(originalv, &dm));
880:   PetscCheck(dm, PetscObjectComm((PetscObject)originalv), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
881:   PetscCall(PetscViewerGetFormat(viewer, &format));
882:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
883:   if (format == PETSC_VIEWER_NATIVE) {
884:     if (dm->useNatural) {
885:       if (dm->sfNatural) {
886:         if (ishdf5) {
887: #if defined(PETSC_HAVE_HDF5)
888:           Vec         v;
889:           const char *vecname;

891:           PetscCall(DMPlexCreateNaturalVector(dm, &v));
892:           PetscCall(PetscObjectGetName((PetscObject)originalv, &vecname));
893:           PetscCall(PetscObjectSetName((PetscObject)v, vecname));
894:           PetscCall(VecLoad_Plex_HDF5_Native_Internal(v, viewer));
895:           PetscCall(DMPlexNaturalToGlobalBegin(dm, v, originalv));
896:           PetscCall(DMPlexNaturalToGlobalEnd(dm, v, originalv));
897:           PetscCall(VecDestroy(&v));
898: #else
899:           SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
900: #endif
901:         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Reading in natural order is not supported for anything but HDF5.");
902:       }
903:     } else PetscCall(VecLoad_Default(originalv, viewer));
904:   }
905:   PetscFunctionReturn(PETSC_SUCCESS);
906: }

908: PETSC_UNUSED static PetscErrorCode DMPlexView_Ascii_Geometry(DM dm, PetscViewer viewer)
909: {
910:   PetscSection       coordSection;
911:   Vec                coordinates;
912:   DMLabel            depthLabel, celltypeLabel;
913:   const char        *name[4];
914:   const PetscScalar *a;
915:   PetscInt           dim, pStart, pEnd, cStart, cEnd, c;

917:   PetscFunctionBegin;
918:   PetscCall(DMGetDimension(dm, &dim));
919:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
920:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
921:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
922:   PetscCall(DMPlexGetCellTypeLabel(dm, &celltypeLabel));
923:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
924:   PetscCall(PetscSectionGetChart(coordSection, &pStart, &pEnd));
925:   PetscCall(VecGetArrayRead(coordinates, &a));
926:   name[0]       = "vertex";
927:   name[1]       = "edge";
928:   name[dim - 1] = "face";
929:   name[dim]     = "cell";
930:   for (c = cStart; c < cEnd; ++c) {
931:     PetscInt *closure = NULL;
932:     PetscInt  closureSize, cl, ct;

934:     PetscCall(DMLabelGetValue(celltypeLabel, c, &ct));
935:     PetscCall(PetscViewerASCIIPrintf(viewer, "Geometry for cell %" PetscInt_FMT " polytope type %s:\n", c, DMPolytopeTypes[ct]));
936:     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
937:     PetscCall(PetscViewerASCIIPushTab(viewer));
938:     for (cl = 0; cl < closureSize * 2; cl += 2) {
939:       PetscInt point = closure[cl], depth, dof, off, d, p;

941:       if ((point < pStart) || (point >= pEnd)) continue;
942:       PetscCall(PetscSectionGetDof(coordSection, point, &dof));
943:       if (!dof) continue;
944:       PetscCall(DMLabelGetValue(depthLabel, point, &depth));
945:       PetscCall(PetscSectionGetOffset(coordSection, point, &off));
946:       PetscCall(PetscViewerASCIIPrintf(viewer, "%s %" PetscInt_FMT " coords:", name[depth], point));
947:       for (p = 0; p < dof / dim; ++p) {
948:         PetscCall(PetscViewerASCIIPrintf(viewer, " ("));
949:         for (d = 0; d < dim; ++d) {
950:           if (d > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
951:           PetscCall(PetscViewerASCIIPrintf(viewer, "%g", (double)PetscRealPart(a[off + p * dim + d])));
952:         }
953:         PetscCall(PetscViewerASCIIPrintf(viewer, ")"));
954:       }
955:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
956:     }
957:     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
958:     PetscCall(PetscViewerASCIIPopTab(viewer));
959:   }
960:   PetscCall(VecRestoreArrayRead(coordinates, &a));
961:   PetscFunctionReturn(PETSC_SUCCESS);
962: }

964: typedef enum {
965:   CS_CARTESIAN,
966:   CS_POLAR,
967:   CS_CYLINDRICAL,
968:   CS_SPHERICAL
969: } CoordSystem;
970: const char *CoordSystems[] = {"cartesian", "polar", "cylindrical", "spherical", "CoordSystem", "CS_", NULL};

972: static PetscErrorCode DMPlexView_Ascii_Coordinates(PetscViewer viewer, CoordSystem cs, PetscInt dim, const PetscScalar x[])
973: {
974:   PetscInt i;

976:   PetscFunctionBegin;
977:   if (dim > 3) {
978:     for (i = 0; i < dim; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %g", (double)PetscRealPart(x[i])));
979:   } else {
980:     PetscReal coords[3], trcoords[3] = {0., 0., 0.};

982:     for (i = 0; i < dim; ++i) coords[i] = PetscRealPart(x[i]);
983:     switch (cs) {
984:     case CS_CARTESIAN:
985:       for (i = 0; i < dim; ++i) trcoords[i] = coords[i];
986:       break;
987:     case CS_POLAR:
988:       PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Polar coordinates are for 2 dimension, not %" PetscInt_FMT, dim);
989:       trcoords[0] = PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1]));
990:       trcoords[1] = PetscAtan2Real(coords[1], coords[0]);
991:       break;
992:     case CS_CYLINDRICAL:
993:       PetscCheck(dim == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cylindrical coordinates are for 3 dimension, not %" PetscInt_FMT, dim);
994:       trcoords[0] = PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1]));
995:       trcoords[1] = PetscAtan2Real(coords[1], coords[0]);
996:       trcoords[2] = coords[2];
997:       break;
998:     case CS_SPHERICAL:
999:       PetscCheck(dim == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Spherical coordinates are for 3 dimension, not %" PetscInt_FMT, dim);
1000:       trcoords[0] = PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1]) + PetscSqr(coords[2]));
1001:       trcoords[1] = PetscAtan2Real(PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1])), coords[2]);
1002:       trcoords[2] = PetscAtan2Real(coords[1], coords[0]);
1003:       break;
1004:     }
1005:     for (i = 0; i < dim; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %g", (double)trcoords[i]));
1006:   }
1007:   PetscFunctionReturn(PETSC_SUCCESS);
1008: }

1010: static PetscErrorCode DMPlexView_Ascii(DM dm, PetscViewer viewer)
1011: {
1012:   DM_Plex          *mesh = (DM_Plex *)dm->data;
1013:   DM                cdm, cdmCell;
1014:   PetscSection      coordSection, coordSectionCell;
1015:   Vec               coordinates, coordinatesCell;
1016:   PetscViewerFormat format;

1018:   PetscFunctionBegin;
1019:   PetscCall(PetscViewerGetFormat(viewer, &format));
1020:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1021:     const char *name;
1022:     PetscInt    dim, cellHeight, maxConeSize, maxSupportSize;
1023:     PetscInt    pStart, pEnd, p, numLabels, l;
1024:     PetscMPIInt rank, size;

1026:     PetscCall(DMGetCoordinateDM(dm, &cdm));
1027:     PetscCall(DMGetCoordinateSection(dm, &coordSection));
1028:     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1029:     PetscCall(DMGetCellCoordinateDM(dm, &cdmCell));
1030:     PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell));
1031:     PetscCall(DMGetCellCoordinatesLocal(dm, &coordinatesCell));
1032:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1033:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1034:     PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1035:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1036:     PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
1037:     PetscCall(DMGetDimension(dm, &dim));
1038:     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1039:     if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
1040:     else PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
1041:     if (cellHeight) PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells are at height %" PetscInt_FMT "\n", cellHeight));
1042:     PetscCall(PetscViewerASCIIPrintf(viewer, "Supports:\n"));
1043:     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1044:     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max support size: %" PetscInt_FMT "\n", rank, maxSupportSize));
1045:     for (p = pStart; p < pEnd; ++p) {
1046:       PetscInt dof, off, s;

1048:       PetscCall(PetscSectionGetDof(mesh->supportSection, p, &dof));
1049:       PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
1050:       for (s = off; s < off + dof; ++s) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", rank, p, mesh->supports[s]));
1051:     }
1052:     PetscCall(PetscViewerFlush(viewer));
1053:     PetscCall(PetscViewerASCIIPrintf(viewer, "Cones:\n"));
1054:     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max cone size: %" PetscInt_FMT "\n", rank, maxConeSize));
1055:     for (p = pStart; p < pEnd; ++p) {
1056:       PetscInt dof, off, c;

1058:       PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
1059:       PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
1060:       for (c = off; c < off + dof; ++c) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %" PetscInt_FMT " <---- %" PetscInt_FMT " (%" PetscInt_FMT ")\n", rank, p, mesh->cones[c], mesh->coneOrientations[c]));
1061:     }
1062:     PetscCall(PetscViewerFlush(viewer));
1063:     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1064:     if (coordSection && coordinates) {
1065:       CoordSystem        cs = CS_CARTESIAN;
1066:       const PetscScalar *array, *arrayCell = NULL;
1067:       PetscInt           Nf, Nc, pvStart, pvEnd, pcStart = PETSC_INT_MAX, pcEnd = PETSC_INT_MIN, pStart, pEnd, p;
1068:       PetscMPIInt        rank;
1069:       const char        *name;

1071:       PetscCall(PetscOptionsGetEnum(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_coord_system", CoordSystems, (PetscEnum *)&cs, NULL));
1072:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1073:       PetscCall(PetscSectionGetNumFields(coordSection, &Nf));
1074:       PetscCheck(Nf == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coordinate section should have 1 field, not %" PetscInt_FMT, Nf);
1075:       PetscCall(PetscSectionGetFieldComponents(coordSection, 0, &Nc));
1076:       PetscCall(PetscSectionGetChart(coordSection, &pvStart, &pvEnd));
1077:       if (coordSectionCell) PetscCall(PetscSectionGetChart(coordSectionCell, &pcStart, &pcEnd));
1078:       pStart = PetscMin(pvStart, pcStart);
1079:       pEnd   = PetscMax(pvEnd, pcEnd);
1080:       PetscCall(PetscObjectGetName((PetscObject)coordinates, &name));
1081:       PetscCall(PetscViewerASCIIPrintf(viewer, "%s with %" PetscInt_FMT " fields\n", name, Nf));
1082:       PetscCall(PetscViewerASCIIPrintf(viewer, "  field 0 with %" PetscInt_FMT " components\n", Nc));
1083:       if (cs != CS_CARTESIAN) PetscCall(PetscViewerASCIIPrintf(viewer, "  output coordinate system: %s\n", CoordSystems[cs]));

1085:       PetscCall(VecGetArrayRead(coordinates, &array));
1086:       if (coordinatesCell) PetscCall(VecGetArrayRead(coordinatesCell, &arrayCell));
1087:       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1088:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank));
1089:       for (p = pStart; p < pEnd; ++p) {
1090:         PetscInt dof, off;

1092:         if (p >= pvStart && p < pvEnd) {
1093:           PetscCall(PetscSectionGetDof(coordSection, p, &dof));
1094:           PetscCall(PetscSectionGetOffset(coordSection, p, &off));
1095:           if (dof) {
1096:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  (%4" PetscInt_FMT ") dof %2" PetscInt_FMT " offset %3" PetscInt_FMT, p, dof, off));
1097:             PetscCall(DMPlexView_Ascii_Coordinates(viewer, cs, dof, &array[off]));
1098:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1099:           }
1100:         }
1101:         if (cdmCell && p >= pcStart && p < pcEnd) {
1102:           PetscCall(PetscSectionGetDof(coordSectionCell, p, &dof));
1103:           PetscCall(PetscSectionGetOffset(coordSectionCell, p, &off));
1104:           if (dof) {
1105:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  (%4" PetscInt_FMT ") dof %2" PetscInt_FMT " offset %3" PetscInt_FMT, p, dof, off));
1106:             PetscCall(DMPlexView_Ascii_Coordinates(viewer, cs, dof, &arrayCell[off]));
1107:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
1108:           }
1109:         }
1110:       }
1111:       PetscCall(PetscViewerFlush(viewer));
1112:       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1113:       PetscCall(VecRestoreArrayRead(coordinates, &array));
1114:       if (coordinatesCell) PetscCall(VecRestoreArrayRead(coordinatesCell, &arrayCell));
1115:     }
1116:     PetscCall(DMGetNumLabels(dm, &numLabels));
1117:     if (numLabels) PetscCall(PetscViewerASCIIPrintf(viewer, "Labels:\n"));
1118:     for (l = 0; l < numLabels; ++l) {
1119:       DMLabel     label;
1120:       PetscBool   isdepth;
1121:       const char *name;

1123:       PetscCall(DMGetLabelName(dm, l, &name));
1124:       PetscCall(PetscStrcmp(name, "depth", &isdepth));
1125:       if (isdepth) continue;
1126:       PetscCall(DMGetLabel(dm, name, &label));
1127:       PetscCall(DMLabelView(label, viewer));
1128:     }
1129:     if (size > 1) {
1130:       PetscSF sf;

1132:       PetscCall(DMGetPointSF(dm, &sf));
1133:       PetscCall(PetscSFView(sf, viewer));
1134:     }
1135:     if (mesh->periodic.face_sfs)
1136:       for (PetscInt i = 0; i < mesh->periodic.num_face_sfs; i++) PetscCall(PetscSFView(mesh->periodic.face_sfs[i], viewer));
1137:     PetscCall(PetscViewerFlush(viewer));
1138:   } else if (format == PETSC_VIEWER_ASCII_LATEX) {
1139:     const char  *name, *color;
1140:     const char  *defcolors[3]  = {"gray", "orange", "green"};
1141:     const char  *deflcolors[4] = {"blue", "cyan", "red", "magenta"};
1142:     char         lname[PETSC_MAX_PATH_LEN];
1143:     PetscReal    scale      = 2.0;
1144:     PetscReal    tikzscale  = 1.0;
1145:     PetscBool    useNumbers = PETSC_TRUE, drawNumbers[4], drawColors[4], useLabels, useColors, plotEdges, drawHasse = PETSC_FALSE;
1146:     double       tcoords[3];
1147:     PetscScalar *coords;
1148:     PetscInt     numLabels, l, numColors, numLColors, dim, d, depth, cStart, cEnd, c, vStart, vEnd, v, eStart = 0, eEnd = 0, fStart = 0, fEnd = 0, e, p, n;
1149:     PetscMPIInt  rank, size;
1150:     char       **names, **colors, **lcolors;
1151:     PetscBool    flg, lflg;
1152:     PetscBT      wp = NULL;
1153:     PetscInt     pEnd, pStart;

1155:     PetscCall(DMGetCoordinateDM(dm, &cdm));
1156:     PetscCall(DMGetCoordinateSection(dm, &coordSection));
1157:     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1158:     PetscCall(DMGetCellCoordinateDM(dm, &cdmCell));
1159:     PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell));
1160:     PetscCall(DMGetCellCoordinatesLocal(dm, &coordinatesCell));
1161:     PetscCall(DMGetDimension(dm, &dim));
1162:     PetscCall(DMPlexGetDepth(dm, &depth));
1163:     PetscCall(DMGetNumLabels(dm, &numLabels));
1164:     numLabels  = PetscMax(numLabels, 10);
1165:     numColors  = 10;
1166:     numLColors = 10;
1167:     PetscCall(PetscCalloc3(numLabels, &names, numColors, &colors, numLColors, &lcolors));
1168:     PetscCall(PetscOptionsGetReal(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_scale", &scale, NULL));
1169:     PetscCall(PetscOptionsGetReal(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_tikzscale", &tikzscale, NULL));
1170:     PetscCall(PetscOptionsGetBool(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_numbers", &useNumbers, NULL));
1171:     for (d = 0; d < 4; ++d) drawNumbers[d] = useNumbers;
1172:     for (d = 0; d < 4; ++d) drawColors[d] = PETSC_TRUE;
1173:     n = 4;
1174:     PetscCall(PetscOptionsGetBoolArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_numbers_depth", drawNumbers, &n, &flg));
1175:     PetscCheck(!flg || n == dim + 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of flags %" PetscInt_FMT " != %" PetscInt_FMT " dim+1", n, dim + 1);
1176:     n = 4;
1177:     PetscCall(PetscOptionsGetBoolArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_colors_depth", drawColors, &n, &flg));
1178:     PetscCheck(!flg || n == dim + 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of flags %" PetscInt_FMT " != %" PetscInt_FMT " dim+1", n, dim + 1);
1179:     PetscCall(PetscOptionsGetStringArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_labels", names, &numLabels, &useLabels));
1180:     if (!useLabels) numLabels = 0;
1181:     PetscCall(PetscOptionsGetStringArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_colors", colors, &numColors, &useColors));
1182:     if (!useColors) {
1183:       numColors = 3;
1184:       for (c = 0; c < numColors; ++c) PetscCall(PetscStrallocpy(defcolors[c], &colors[c]));
1185:     }
1186:     PetscCall(PetscOptionsGetStringArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_lcolors", lcolors, &numLColors, &useColors));
1187:     if (!useColors) {
1188:       numLColors = 4;
1189:       for (c = 0; c < numLColors; ++c) PetscCall(PetscStrallocpy(deflcolors[c], &lcolors[c]));
1190:     }
1191:     PetscCall(PetscOptionsGetString(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_label_filter", lname, sizeof(lname), &lflg));
1192:     plotEdges = (PetscBool)(depth > 1 && drawNumbers[1] && dim < 3);
1193:     PetscCall(PetscOptionsGetBool(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_edges", &plotEdges, &flg));
1194:     PetscCheck(!flg || !plotEdges || depth >= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh must be interpolated");
1195:     if (depth < dim) plotEdges = PETSC_FALSE;
1196:     PetscCall(PetscOptionsGetBool(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_hasse", &drawHasse, NULL));

1198:     /* filter points with labelvalue != labeldefaultvalue */
1199:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1200:     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1201:     PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
1202:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1203:     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
1204:     if (lflg) {
1205:       DMLabel lbl;

1207:       PetscCall(DMGetLabel(dm, lname, &lbl));
1208:       if (lbl) {
1209:         PetscInt val, defval;

1211:         PetscCall(DMLabelGetDefaultValue(lbl, &defval));
1212:         PetscCall(PetscBTCreate(pEnd - pStart, &wp));
1213:         for (c = pStart; c < pEnd; c++) {
1214:           PetscInt *closure = NULL;
1215:           PetscInt  closureSize;

1217:           PetscCall(DMLabelGetValue(lbl, c, &val));
1218:           if (val == defval) continue;

1220:           PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1221:           for (p = 0; p < closureSize * 2; p += 2) PetscCall(PetscBTSet(wp, closure[p] - pStart));
1222:           PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1223:         }
1224:       }
1225:     }

1227:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1228:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1229:     PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1230:     PetscCall(PetscViewerASCIIPrintf(viewer, "\
1231: \\documentclass[tikz]{standalone}\n\n\
1232: \\usepackage{pgflibraryshapes}\n\
1233: \\usetikzlibrary{backgrounds}\n\
1234: \\usetikzlibrary{arrows}\n\
1235: \\begin{document}\n"));
1236:     if (size > 1) {
1237:       PetscCall(PetscViewerASCIIPrintf(viewer, "%s for process ", name));
1238:       for (p = 0; p < size; ++p) {
1239:         if (p) PetscCall(PetscViewerASCIIPrintf(viewer, (p == size - 1) ? ", and " : ", "));
1240:         PetscCall(PetscViewerASCIIPrintf(viewer, "{\\textcolor{%s}%" PetscInt_FMT "}", colors[p % numColors], p));
1241:       }
1242:       PetscCall(PetscViewerASCIIPrintf(viewer, ".\n\n\n"));
1243:     }
1244:     if (drawHasse) {
1245:       PetscInt maxStratum = PetscMax(vEnd - vStart, PetscMax(eEnd - eStart, PetscMax(fEnd - fStart, cEnd - cStart)));

1247:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vStart}{%" PetscInt_FMT "}\n", vStart));
1248:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vEnd}{%" PetscInt_FMT "}\n", vEnd - 1));
1249:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numVertices}{%" PetscInt_FMT "}\n", vEnd - vStart));
1250:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vShift}{%.2f}\n", 3 + (maxStratum - (vEnd - vStart)) / 2.));
1251:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eStart}{%" PetscInt_FMT "}\n", eStart));
1252:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eEnd}{%" PetscInt_FMT "}\n", eEnd - 1));
1253:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eShift}{%.2f}\n", 3 + (maxStratum - (eEnd - eStart)) / 2.));
1254:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numEdges}{%" PetscInt_FMT "}\n", eEnd - eStart));
1255:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\fStart}{%" PetscInt_FMT "}\n", fStart));
1256:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\fEnd}{%" PetscInt_FMT "}\n", fEnd - 1));
1257:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\fShift}{%.2f}\n", 3 + (maxStratum - (fEnd - fStart)) / 2.));
1258:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numFaces}{%" PetscInt_FMT "}\n", fEnd - fStart));
1259:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cStart}{%" PetscInt_FMT "}\n", cStart));
1260:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cEnd}{%" PetscInt_FMT "}\n", cEnd - 1));
1261:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numCells}{%" PetscInt_FMT "}\n", cEnd - cStart));
1262:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cShift}{%.2f}\n", 3 + (maxStratum - (cEnd - cStart)) / 2.));
1263:     }
1264:     PetscCall(PetscViewerASCIIPrintf(viewer, "\\begin{tikzpicture}[scale = %g,font=\\fontsize{8}{8}\\selectfont]\n", (double)tikzscale));

1266:     /* Plot vertices */
1267:     PetscCall(VecGetArray(coordinates, &coords));
1268:     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1269:     for (v = vStart; v < vEnd; ++v) {
1270:       PetscInt  off, dof, d;
1271:       PetscBool isLabeled = PETSC_FALSE;

1273:       if (wp && !PetscBTLookup(wp, v - pStart)) continue;
1274:       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
1275:       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
1276:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\path ("));
1277:       PetscCheck(dof <= 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "coordSection vertex %" PetscInt_FMT " has dof %" PetscInt_FMT " > 3", v, dof);
1278:       for (d = 0; d < dof; ++d) {
1279:         tcoords[d] = (double)(scale * PetscRealPart(coords[off + d]));
1280:         tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1281:       }
1282:       /* Rotate coordinates since PGF makes z point out of the page instead of up */
1283:       if (dim == 3) {
1284:         PetscReal tmp = tcoords[1];
1285:         tcoords[1]    = tcoords[2];
1286:         tcoords[2]    = -tmp;
1287:       }
1288:       for (d = 0; d < dof; ++d) {
1289:         if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ","));
1290:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", tcoords[d]));
1291:       }
1292:       if (drawHasse) color = colors[0 % numColors];
1293:       else color = colors[rank % numColors];
1294:       for (l = 0; l < numLabels; ++l) {
1295:         PetscInt val;
1296:         PetscCall(DMGetLabelValue(dm, names[l], v, &val));
1297:         if (val >= 0) {
1298:           color     = lcolors[l % numLColors];
1299:           isLabeled = PETSC_TRUE;
1300:           break;
1301:         }
1302:       }
1303:       if (drawNumbers[0]) {
1304:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [draw,shape=circle,color=%s] {%" PetscInt_FMT "};\n", v, rank, color, v));
1305:       } else if (drawColors[0]) {
1306:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", v, rank, !isLabeled ? 1 : 2, color));
1307:       } else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [] {};\n", v, rank));
1308:     }
1309:     PetscCall(VecRestoreArray(coordinates, &coords));
1310:     PetscCall(PetscViewerFlush(viewer));
1311:     /* Plot edges */
1312:     if (plotEdges) {
1313:       PetscCall(VecGetArray(coordinates, &coords));
1314:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\path\n"));
1315:       for (e = eStart; e < eEnd; ++e) {
1316:         const PetscInt *cone;
1317:         PetscInt        coneSize, offA, offB, dof, d;

1319:         if (wp && !PetscBTLookup(wp, e - pStart)) continue;
1320:         PetscCall(DMPlexGetConeSize(dm, e, &coneSize));
1321:         PetscCheck(coneSize == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " cone should have two vertices, not %" PetscInt_FMT, e, coneSize);
1322:         PetscCall(DMPlexGetCone(dm, e, &cone));
1323:         PetscCall(PetscSectionGetDof(coordSection, cone[0], &dof));
1324:         PetscCall(PetscSectionGetOffset(coordSection, cone[0], &offA));
1325:         PetscCall(PetscSectionGetOffset(coordSection, cone[1], &offB));
1326:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "("));
1327:         for (d = 0; d < dof; ++d) {
1328:           tcoords[d] = (double)(scale * PetscRealPart(coords[offA + d] + coords[offB + d]) / 2);
1329:           tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1330:         }
1331:         /* Rotate coordinates since PGF makes z point out of the page instead of up */
1332:         if (dim == 3) {
1333:           PetscReal tmp = tcoords[1];
1334:           tcoords[1]    = tcoords[2];
1335:           tcoords[2]    = -tmp;
1336:         }
1337:         for (d = 0; d < dof; ++d) {
1338:           if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ","));
1339:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", tcoords[d]));
1340:         }
1341:         if (drawHasse) color = colors[1 % numColors];
1342:         else color = colors[rank % numColors];
1343:         for (l = 0; l < numLabels; ++l) {
1344:           PetscInt val;
1345:           PetscCall(DMGetLabelValue(dm, names[l], e, &val));
1346:           if (val >= 0) {
1347:             color = lcolors[l % numLColors];
1348:             break;
1349:           }
1350:         }
1351:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [draw,shape=circle,color=%s] {%" PetscInt_FMT "} --\n", e, rank, color, e));
1352:       }
1353:       PetscCall(VecRestoreArray(coordinates, &coords));
1354:       PetscCall(PetscViewerFlush(viewer));
1355:       PetscCall(PetscViewerASCIIPrintf(viewer, "(0,0);\n"));
1356:     }
1357:     /* Plot cells */
1358:     if (dim == 3 || !drawNumbers[1]) {
1359:       for (e = eStart; e < eEnd; ++e) {
1360:         const PetscInt *cone;

1362:         if (wp && !PetscBTLookup(wp, e - pStart)) continue;
1363:         color = colors[rank % numColors];
1364:         for (l = 0; l < numLabels; ++l) {
1365:           PetscInt val;
1366:           PetscCall(DMGetLabelValue(dm, names[l], e, &val));
1367:           if (val >= 0) {
1368:             color = lcolors[l % numLColors];
1369:             break;
1370:           }
1371:         }
1372:         PetscCall(DMPlexGetCone(dm, e, &cone));
1373:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] (%" PetscInt_FMT "_%d) -- (%" PetscInt_FMT "_%d);\n", color, cone[0], rank, cone[1], rank));
1374:       }
1375:     } else {
1376:       DMPolytopeType ct;

1378:       /* Drawing a 2D polygon */
1379:       for (c = cStart; c < cEnd; ++c) {
1380:         if (wp && !PetscBTLookup(wp, c - pStart)) continue;
1381:         PetscCall(DMPlexGetCellType(dm, c, &ct));
1382:         if (DMPolytopeTypeIsHybrid(ct)) {
1383:           const PetscInt *cone;
1384:           PetscInt        coneSize, e;

1386:           PetscCall(DMPlexGetCone(dm, c, &cone));
1387:           PetscCall(DMPlexGetConeSize(dm, c, &coneSize));
1388:           for (e = 0; e < coneSize; ++e) {
1389:             const PetscInt *econe;

1391:             PetscCall(DMPlexGetCone(dm, cone[e], &econe));
1392:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] (%" PetscInt_FMT "_%d) -- (%" PetscInt_FMT "_%d) -- (%" PetscInt_FMT "_%d);\n", colors[rank % numColors], econe[0], rank, cone[e], rank, econe[1], rank));
1393:           }
1394:         } else {
1395:           PetscInt *closure = NULL;
1396:           PetscInt  closureSize, Nv = 0, v;

1398:           PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1399:           for (p = 0; p < closureSize * 2; p += 2) {
1400:             const PetscInt point = closure[p];

1402:             if ((point >= vStart) && (point < vEnd)) closure[Nv++] = point;
1403:           }
1404:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] ", colors[rank % numColors]));
1405:           for (v = 0; v <= Nv; ++v) {
1406:             const PetscInt vertex = closure[v % Nv];

1408:             if (v > 0) {
1409:               if (plotEdges) {
1410:                 const PetscInt *edge;
1411:                 PetscInt        endpoints[2], ne;

1413:                 endpoints[0] = closure[v - 1];
1414:                 endpoints[1] = vertex;
1415:                 PetscCall(DMPlexGetJoin(dm, 2, endpoints, &ne, &edge));
1416:                 PetscCheck(ne == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find edge for vertices %" PetscInt_FMT ", %" PetscInt_FMT, endpoints[0], endpoints[1]);
1417:                 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " -- (%" PetscInt_FMT "_%d) -- ", edge[0], rank));
1418:                 PetscCall(DMPlexRestoreJoin(dm, 2, endpoints, &ne, &edge));
1419:               } else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " -- "));
1420:             }
1421:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "(%" PetscInt_FMT "_%d)", vertex, rank));
1422:           }
1423:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ";\n"));
1424:           PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1425:         }
1426:       }
1427:     }
1428:     for (c = cStart; c < cEnd; ++c) {
1429:       double             ccoords[3] = {0.0, 0.0, 0.0};
1430:       PetscBool          isLabeled  = PETSC_FALSE;
1431:       PetscScalar       *cellCoords = NULL;
1432:       const PetscScalar *array;
1433:       PetscInt           numCoords, cdim, d;
1434:       PetscBool          isDG;

1436:       if (wp && !PetscBTLookup(wp, c - pStart)) continue;
1437:       PetscCall(DMGetCoordinateDim(dm, &cdim));
1438:       PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &cellCoords));
1439:       PetscCheck(!(numCoords % cdim), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "coordinate dim %" PetscInt_FMT " does not divide numCoords %" PetscInt_FMT, cdim, numCoords);
1440:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\path ("));
1441:       for (p = 0; p < numCoords / cdim; ++p) {
1442:         for (d = 0; d < cdim; ++d) {
1443:           tcoords[d] = (double)(scale * PetscRealPart(cellCoords[p * cdim + d]));
1444:           tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1445:         }
1446:         /* Rotate coordinates since PGF makes z point out of the page instead of up */
1447:         if (cdim == 3) {
1448:           PetscReal tmp = tcoords[1];
1449:           tcoords[1]    = tcoords[2];
1450:           tcoords[2]    = -tmp;
1451:         }
1452:         for (d = 0; d < dim; ++d) ccoords[d] += tcoords[d];
1453:       }
1454:       for (d = 0; d < cdim; ++d) ccoords[d] /= (numCoords / cdim);
1455:       PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &cellCoords));
1456:       for (d = 0; d < cdim; ++d) {
1457:         if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ","));
1458:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", ccoords[d]));
1459:       }
1460:       if (drawHasse) color = colors[depth % numColors];
1461:       else color = colors[rank % numColors];
1462:       for (l = 0; l < numLabels; ++l) {
1463:         PetscInt val;
1464:         PetscCall(DMGetLabelValue(dm, names[l], c, &val));
1465:         if (val >= 0) {
1466:           color     = lcolors[l % numLColors];
1467:           isLabeled = PETSC_TRUE;
1468:           break;
1469:         }
1470:       }
1471:       if (drawNumbers[dim]) {
1472:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [draw,shape=circle,color=%s] {%" PetscInt_FMT "};\n", c, rank, color, c));
1473:       } else if (drawColors[dim]) {
1474:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", c, rank, !isLabeled ? 1 : 2, color));
1475:       } else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [] {};\n", c, rank));
1476:     }
1477:     if (drawHasse) {
1478:       int height = 0;

1480:       color = colors[depth % numColors];
1481:       PetscCall(PetscViewerASCIIPrintf(viewer, "%% Cells\n"));
1482:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\c in {\\cStart,...,\\cEnd}\n"));
1483:       PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1484:       PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\c_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\cShift+\\c-\\cStart,%d) {\\c};\n", rank, color, height++));
1485:       PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));

1487:       if (depth > 2) {
1488:         color = colors[1 % numColors];
1489:         PetscCall(PetscViewerASCIIPrintf(viewer, "%% Faces\n"));
1490:         PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\f in {\\fStart,...,\\fEnd}\n"));
1491:         PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1492:         PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\f_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\fShift+\\f-\\fStart,%d) {\\f};\n", rank, color, height++));
1493:         PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));
1494:       }

1496:       color = colors[1 % numColors];
1497:       PetscCall(PetscViewerASCIIPrintf(viewer, "%% Edges\n"));
1498:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\e in {\\eStart,...,\\eEnd}\n"));
1499:       PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1500:       PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\e_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\eShift+\\e-\\eStart,%d) {\\e};\n", rank, color, height++));
1501:       PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));

1503:       color = colors[0 % numColors];
1504:       PetscCall(PetscViewerASCIIPrintf(viewer, "%% Vertices\n"));
1505:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\v in {\\vStart,...,\\vEnd}\n"));
1506:       PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1507:       PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\v_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\vShift+\\v-\\vStart,%d) {\\v};\n", rank, color, height++));
1508:       PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));

1510:       for (p = pStart; p < pEnd; ++p) {
1511:         const PetscInt *cone;
1512:         PetscInt        coneSize, cp;

1514:         PetscCall(DMPlexGetCone(dm, p, &cone));
1515:         PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
1516:         for (cp = 0; cp < coneSize; ++cp) PetscCall(PetscViewerASCIIPrintf(viewer, "\\draw[->, shorten >=1pt] (%" PetscInt_FMT "_%d) -- (%" PetscInt_FMT "_%d);\n", cone[cp], rank, p, rank));
1517:       }
1518:     }
1519:     PetscCall(PetscViewerFlush(viewer));
1520:     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1521:     PetscCall(PetscViewerASCIIPrintf(viewer, "\\end{tikzpicture}\n"));
1522:     PetscCall(PetscViewerASCIIPrintf(viewer, "\\end{document}\n"));
1523:     for (l = 0; l < numLabels; ++l) PetscCall(PetscFree(names[l]));
1524:     for (c = 0; c < numColors; ++c) PetscCall(PetscFree(colors[c]));
1525:     for (c = 0; c < numLColors; ++c) PetscCall(PetscFree(lcolors[c]));
1526:     PetscCall(PetscFree3(names, colors, lcolors));
1527:     PetscCall(PetscBTDestroy(&wp));
1528:   } else if (format == PETSC_VIEWER_LOAD_BALANCE) {
1529:     Vec                    cown, acown;
1530:     VecScatter             sct;
1531:     ISLocalToGlobalMapping g2l;
1532:     IS                     gid, acis;
1533:     MPI_Comm               comm, ncomm = MPI_COMM_NULL;
1534:     MPI_Group              ggroup, ngroup;
1535:     PetscScalar           *array, nid;
1536:     const PetscInt        *idxs;
1537:     PetscInt              *idxs2, *start, *adjacency, *work;
1538:     PetscInt64             lm[3], gm[3];
1539:     PetscInt               i, c, cStart, cEnd, cum, numVertices, ect, ectn, cellHeight;
1540:     PetscMPIInt            d1, d2, rank;

1542:     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1543:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1544: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
1545:     PetscCallMPI(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, rank, MPI_INFO_NULL, &ncomm));
1546: #endif
1547:     if (ncomm != MPI_COMM_NULL) {
1548:       PetscCallMPI(MPI_Comm_group(comm, &ggroup));
1549:       PetscCallMPI(MPI_Comm_group(ncomm, &ngroup));
1550:       d1 = 0;
1551:       PetscCallMPI(MPI_Group_translate_ranks(ngroup, 1, &d1, ggroup, &d2));
1552:       nid = d2;
1553:       PetscCallMPI(MPI_Group_free(&ggroup));
1554:       PetscCallMPI(MPI_Group_free(&ngroup));
1555:       PetscCallMPI(MPI_Comm_free(&ncomm));
1556:     } else nid = 0.0;

1558:     /* Get connectivity */
1559:     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1560:     PetscCall(DMPlexCreatePartitionerGraph(dm, cellHeight, &numVertices, &start, &adjacency, &gid));

1562:     /* filter overlapped local cells */
1563:     PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
1564:     PetscCall(ISGetIndices(gid, &idxs));
1565:     PetscCall(ISGetLocalSize(gid, &cum));
1566:     PetscCall(PetscMalloc1(cum, &idxs2));
1567:     for (c = cStart, cum = 0; c < cEnd; c++) {
1568:       if (idxs[c - cStart] < 0) continue;
1569:       idxs2[cum++] = idxs[c - cStart];
1570:     }
1571:     PetscCall(ISRestoreIndices(gid, &idxs));
1572:     PetscCheck(numVertices == cum, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected %" PetscInt_FMT " != %" PetscInt_FMT, numVertices, cum);
1573:     PetscCall(ISDestroy(&gid));
1574:     PetscCall(ISCreateGeneral(comm, numVertices, idxs2, PETSC_OWN_POINTER, &gid));

1576:     /* support for node-aware cell locality */
1577:     PetscCall(ISCreateGeneral(comm, start[numVertices], adjacency, PETSC_USE_POINTER, &acis));
1578:     PetscCall(VecCreateSeq(PETSC_COMM_SELF, start[numVertices], &acown));
1579:     PetscCall(VecCreateMPI(comm, numVertices, PETSC_DECIDE, &cown));
1580:     PetscCall(VecGetArray(cown, &array));
1581:     for (c = 0; c < numVertices; c++) array[c] = nid;
1582:     PetscCall(VecRestoreArray(cown, &array));
1583:     PetscCall(VecScatterCreate(cown, acis, acown, NULL, &sct));
1584:     PetscCall(VecScatterBegin(sct, cown, acown, INSERT_VALUES, SCATTER_FORWARD));
1585:     PetscCall(VecScatterEnd(sct, cown, acown, INSERT_VALUES, SCATTER_FORWARD));
1586:     PetscCall(ISDestroy(&acis));
1587:     PetscCall(VecScatterDestroy(&sct));
1588:     PetscCall(VecDestroy(&cown));

1590:     /* compute edgeCut */
1591:     for (c = 0, cum = 0; c < numVertices; c++) cum = PetscMax(cum, start[c + 1] - start[c]);
1592:     PetscCall(PetscMalloc1(cum, &work));
1593:     PetscCall(ISLocalToGlobalMappingCreateIS(gid, &g2l));
1594:     PetscCall(ISLocalToGlobalMappingSetType(g2l, ISLOCALTOGLOBALMAPPINGHASH));
1595:     PetscCall(ISDestroy(&gid));
1596:     PetscCall(VecGetArray(acown, &array));
1597:     for (c = 0, ect = 0, ectn = 0; c < numVertices; c++) {
1598:       PetscInt totl;

1600:       totl = start[c + 1] - start[c];
1601:       PetscCall(ISGlobalToLocalMappingApply(g2l, IS_GTOLM_MASK, totl, adjacency + start[c], NULL, work));
1602:       for (i = 0; i < totl; i++) {
1603:         if (work[i] < 0) {
1604:           ect += 1;
1605:           ectn += (array[i + start[c]] != nid) ? 0 : 1;
1606:         }
1607:       }
1608:     }
1609:     PetscCall(PetscFree(work));
1610:     PetscCall(VecRestoreArray(acown, &array));
1611:     lm[0] = numVertices > 0 ? numVertices : PETSC_INT_MAX;
1612:     lm[1] = -numVertices;
1613:     PetscCallMPI(MPIU_Allreduce(lm, gm, 2, MPIU_INT64, MPI_MIN, comm));
1614:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cell balance: %.2f (max %" PetscInt64_FMT ", min %" PetscInt64_FMT, -((double)gm[1]) / ((double)gm[0]), -gm[1], gm[0]));
1615:     lm[0] = ect;                     /* edgeCut */
1616:     lm[1] = ectn;                    /* node-aware edgeCut */
1617:     lm[2] = numVertices > 0 ? 0 : 1; /* empty processes */
1618:     PetscCallMPI(MPIU_Allreduce(lm, gm, 3, MPIU_INT64, MPI_SUM, comm));
1619:     PetscCall(PetscViewerASCIIPrintf(viewer, ", empty %" PetscInt64_FMT ")\n", gm[2]));
1620: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
1621:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge Cut: %" PetscInt64_FMT " (on node %.3f)\n", gm[0] / 2, gm[0] ? ((double)gm[1]) / ((double)gm[0]) : 1.));
1622: #else
1623:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge Cut: %" PetscInt64_FMT " (on node %.3f)\n", gm[0] / 2, 0.0));
1624: #endif
1625:     PetscCall(ISLocalToGlobalMappingDestroy(&g2l));
1626:     PetscCall(PetscFree(start));
1627:     PetscCall(PetscFree(adjacency));
1628:     PetscCall(VecDestroy(&acown));
1629:   } else {
1630:     const char    *name;
1631:     PetscInt      *sizes, *hybsizes, *ghostsizes;
1632:     PetscInt       locDepth, depth, cellHeight, dim, d;
1633:     PetscInt       pStart, pEnd, p, gcStart, gcEnd, gcNum;
1634:     PetscInt       numLabels, l, maxSize = 17;
1635:     DMPolytopeType ct0 = DM_POLYTOPE_UNKNOWN;
1636:     MPI_Comm       comm;
1637:     PetscMPIInt    size, rank;

1639:     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1640:     PetscCallMPI(MPI_Comm_size(comm, &size));
1641:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1642:     PetscCall(DMGetDimension(dm, &dim));
1643:     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1644:     PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1645:     if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
1646:     else PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
1647:     if (cellHeight) PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells are at height %" PetscInt_FMT "\n", cellHeight));
1648:     PetscCall(DMPlexGetDepth(dm, &locDepth));
1649:     PetscCallMPI(MPIU_Allreduce(&locDepth, &depth, 1, MPIU_INT, MPI_MAX, comm));
1650:     PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &gcStart, &gcEnd));
1651:     gcNum = gcEnd - gcStart;
1652:     if (size < maxSize) PetscCall(PetscCalloc3(size, &sizes, size, &hybsizes, size, &ghostsizes));
1653:     else PetscCall(PetscCalloc3(3, &sizes, 3, &hybsizes, 3, &ghostsizes));
1654:     for (d = 0; d <= depth; d++) {
1655:       PetscInt Nc[2] = {0, 0}, ict;

1657:       PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
1658:       if (pStart < pEnd) PetscCall(DMPlexGetCellType(dm, pStart, &ct0));
1659:       ict = ct0;
1660:       PetscCallMPI(MPI_Bcast(&ict, 1, MPIU_INT, 0, comm));
1661:       ct0 = (DMPolytopeType)ict;
1662:       for (p = pStart; p < pEnd; ++p) {
1663:         DMPolytopeType ct;

1665:         PetscCall(DMPlexGetCellType(dm, p, &ct));
1666:         if (ct == ct0) ++Nc[0];
1667:         else ++Nc[1];
1668:       }
1669:       if (size < maxSize) {
1670:         PetscCallMPI(MPI_Gather(&Nc[0], 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
1671:         PetscCallMPI(MPI_Gather(&Nc[1], 1, MPIU_INT, hybsizes, 1, MPIU_INT, 0, comm));
1672:         if (d == depth) PetscCallMPI(MPI_Gather(&gcNum, 1, MPIU_INT, ghostsizes, 1, MPIU_INT, 0, comm));
1673:         PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of %" PetscInt_FMT "-cells per rank:", (depth == 1) && d ? dim : d));
1674:         for (p = 0; p < size; ++p) {
1675:           if (rank == 0) {
1676:             PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p] + hybsizes[p]));
1677:             if (hybsizes[p] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ")", hybsizes[p]));
1678:             if (ghostsizes[p] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " [%" PetscInt_FMT "]", ghostsizes[p]));
1679:           }
1680:         }
1681:       } else {
1682:         PetscInt locMinMax[2];

1684:         locMinMax[0] = Nc[0] + Nc[1];
1685:         locMinMax[1] = Nc[0] + Nc[1];
1686:         PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
1687:         locMinMax[0] = Nc[1];
1688:         locMinMax[1] = Nc[1];
1689:         PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, hybsizes));
1690:         if (d == depth) {
1691:           locMinMax[0] = gcNum;
1692:           locMinMax[1] = gcNum;
1693:           PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, ghostsizes));
1694:         }
1695:         PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of %" PetscInt_FMT "-cells per rank:", (depth == 1) && d ? dim : d));
1696:         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
1697:         if (hybsizes[0] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT "/%" PetscInt_FMT ")", hybsizes[0], hybsizes[1]));
1698:         if (ghostsizes[0] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " [%" PetscInt_FMT "/%" PetscInt_FMT "]", ghostsizes[0], ghostsizes[1]));
1699:       }
1700:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1701:     }
1702:     PetscCall(PetscFree3(sizes, hybsizes, ghostsizes));
1703:     {
1704:       const PetscReal *maxCell;
1705:       const PetscReal *L;
1706:       PetscBool        localized;

1708:       PetscCall(DMGetPeriodicity(dm, &maxCell, NULL, &L));
1709:       PetscCall(DMGetCoordinatesLocalized(dm, &localized));
1710:       if (L || localized) {
1711:         PetscCall(PetscViewerASCIIPrintf(viewer, "Periodic mesh"));
1712:         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1713:         if (L) {
1714:           PetscCall(PetscViewerASCIIPrintf(viewer, " ("));
1715:           for (d = 0; d < dim; ++d) {
1716:             if (d > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1717:             PetscCall(PetscViewerASCIIPrintf(viewer, "%s", L[d] > 0.0 ? "PERIODIC" : "NONE"));
1718:           }
1719:           PetscCall(PetscViewerASCIIPrintf(viewer, ")"));
1720:         }
1721:         PetscCall(PetscViewerASCIIPrintf(viewer, " coordinates %s\n", localized ? "localized" : "not localized"));
1722:         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1723:       }
1724:     }
1725:     PetscCall(DMGetNumLabels(dm, &numLabels));
1726:     if (numLabels) PetscCall(PetscViewerASCIIPrintf(viewer, "Labels:\n"));
1727:     for (l = 0; l < numLabels; ++l) {
1728:       DMLabel     label;
1729:       const char *name;
1730:       PetscInt   *values;
1731:       PetscInt    numValues, v;

1733:       PetscCall(DMGetLabelName(dm, l, &name));
1734:       PetscCall(DMGetLabel(dm, name, &label));
1735:       PetscCall(DMLabelGetNumValues(label, &numValues));
1736:       PetscCall(PetscViewerASCIIPrintf(viewer, "  %s: %" PetscInt_FMT " strata with value/size (", name, numValues));

1738:       { // Extract array of DMLabel values so it can be sorted
1739:         IS              is_values;
1740:         const PetscInt *is_values_local = NULL;

1742:         PetscCall(DMLabelGetValueIS(label, &is_values));
1743:         PetscCall(ISGetIndices(is_values, &is_values_local));
1744:         PetscCall(PetscMalloc1(numValues, &values));
1745:         PetscCall(PetscArraycpy(values, is_values_local, numValues));
1746:         PetscCall(PetscSortInt(numValues, values));
1747:         PetscCall(ISRestoreIndices(is_values, &is_values_local));
1748:         PetscCall(ISDestroy(&is_values));
1749:       }
1750:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1751:       for (v = 0; v < numValues; ++v) {
1752:         PetscInt size;

1754:         PetscCall(DMLabelGetStratumSize(label, values[v], &size));
1755:         if (v > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1756:         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " (%" PetscInt_FMT ")", values[v], size));
1757:       }
1758:       PetscCall(PetscViewerASCIIPrintf(viewer, ")\n"));
1759:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1760:       PetscCall(PetscFree(values));
1761:     }
1762:     {
1763:       char    **labelNames;
1764:       PetscInt  Nl = numLabels;
1765:       PetscBool flg;

1767:       PetscCall(PetscMalloc1(Nl, &labelNames));
1768:       PetscCall(PetscOptionsGetStringArray(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_view_labels", labelNames, &Nl, &flg));
1769:       for (l = 0; l < Nl; ++l) {
1770:         DMLabel label;

1772:         PetscCall(DMHasLabel(dm, labelNames[l], &flg));
1773:         if (flg) {
1774:           PetscCall(DMGetLabel(dm, labelNames[l], &label));
1775:           PetscCall(DMLabelView(label, viewer));
1776:         }
1777:         PetscCall(PetscFree(labelNames[l]));
1778:       }
1779:       PetscCall(PetscFree(labelNames));
1780:     }
1781:     /* If no fields are specified, people do not want to see adjacency */
1782:     if (dm->Nf) {
1783:       PetscInt f;

1785:       for (f = 0; f < dm->Nf; ++f) {
1786:         const char *name;

1788:         PetscCall(PetscObjectGetName(dm->fields[f].disc, &name));
1789:         if (numLabels) PetscCall(PetscViewerASCIIPrintf(viewer, "Field %s:\n", name));
1790:         PetscCall(PetscViewerASCIIPushTab(viewer));
1791:         if (dm->fields[f].label) PetscCall(DMLabelView(dm->fields[f].label, viewer));
1792:         if (dm->fields[f].adjacency[0]) {
1793:           if (dm->fields[f].adjacency[1]) PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FVM++\n"));
1794:           else PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FVM\n"));
1795:         } else {
1796:           if (dm->fields[f].adjacency[1]) PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FEM\n"));
1797:           else PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FUNKY\n"));
1798:         }
1799:         PetscCall(PetscViewerASCIIPopTab(viewer));
1800:       }
1801:     }
1802:     PetscCall(DMGetCoarseDM(dm, &cdm));
1803:     if (cdm) {
1804:       PetscCall(PetscViewerASCIIPushTab(viewer));
1805:       PetscCall(PetscViewerASCIIPrintf(viewer, "Defined by transform from:\n"));
1806:       PetscCall(DMPlexView_Ascii(cdm, viewer));
1807:       PetscCall(PetscViewerASCIIPopTab(viewer));
1808:     }
1809:   }
1810:   PetscFunctionReturn(PETSC_SUCCESS);
1811: }

1813: static PetscErrorCode DMPlexDrawCell(DM dm, PetscDraw draw, PetscInt cell, const PetscScalar coords[])
1814: {
1815:   DMPolytopeType ct;
1816:   PetscMPIInt    rank;
1817:   PetscInt       cdim;

1819:   PetscFunctionBegin;
1820:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1821:   PetscCall(DMPlexGetCellType(dm, cell, &ct));
1822:   PetscCall(DMGetCoordinateDim(dm, &cdim));
1823:   switch (ct) {
1824:   case DM_POLYTOPE_SEGMENT:
1825:   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1826:     switch (cdim) {
1827:     case 1: {
1828:       const PetscReal y  = 0.5;  /* TODO Put it in the middle of the viewport */
1829:       const PetscReal dy = 0.05; /* TODO Make it a fraction of the total length */

1831:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), y, PetscRealPart(coords[1]), y, PETSC_DRAW_BLACK));
1832:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), y + dy, PetscRealPart(coords[0]), y - dy, PETSC_DRAW_BLACK));
1833:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[1]), y + dy, PetscRealPart(coords[1]), y - dy, PETSC_DRAW_BLACK));
1834:     } break;
1835:     case 2: {
1836:       const PetscReal dx = (PetscRealPart(coords[3]) - PetscRealPart(coords[1]));
1837:       const PetscReal dy = (PetscRealPart(coords[2]) - PetscRealPart(coords[0]));
1838:       const PetscReal l  = 0.1 / PetscSqrtReal(dx * dx + dy * dy);

1840:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1841:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]) + l * dx, PetscRealPart(coords[1]) + l * dy, PetscRealPart(coords[0]) - l * dx, PetscRealPart(coords[1]) - l * dy, PETSC_DRAW_BLACK));
1842:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]) + l * dx, PetscRealPart(coords[3]) + l * dy, PetscRealPart(coords[2]) - l * dx, PetscRealPart(coords[3]) - l * dy, PETSC_DRAW_BLACK));
1843:     } break;
1844:     default:
1845:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of dimension %" PetscInt_FMT, cdim);
1846:     }
1847:     break;
1848:   case DM_POLYTOPE_TRIANGLE:
1849:     PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2));
1850:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1851:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK));
1852:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK));
1853:     break;
1854:   case DM_POLYTOPE_QUADRILATERAL:
1855:     PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2));
1856:     PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2));
1857:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1858:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK));
1859:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_BLACK));
1860:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK));
1861:     break;
1862:   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1863:     PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2));
1864:     PetscCall(PetscDrawTriangle(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2, PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2));
1865:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1866:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_BLACK));
1867:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK));
1868:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK));
1869:     break;
1870:   case DM_POLYTOPE_FV_GHOST:
1871:     break;
1872:   default:
1873:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of type %s", DMPolytopeTypes[ct]);
1874:   }
1875:   PetscFunctionReturn(PETSC_SUCCESS);
1876: }

1878: static PetscErrorCode DrawPolygon_Private(DM dm, PetscDraw draw, PetscInt cell, PetscInt Nv, const PetscReal refVertices[], const PetscScalar coords[], PetscInt edgeDiv, PetscReal refCoords[], PetscReal edgeCoords[])
1879: {
1880:   PetscReal   centroid[2] = {0., 0.};
1881:   PetscMPIInt rank;
1882:   PetscMPIInt fillColor;

1884:   PetscFunctionBegin;
1885:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1886:   fillColor = PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2;
1887:   for (PetscInt v = 0; v < Nv; ++v) {
1888:     centroid[0] += PetscRealPart(coords[v * 2 + 0]) / Nv;
1889:     centroid[1] += PetscRealPart(coords[v * 2 + 1]) / Nv;
1890:   }
1891:   for (PetscInt e = 0; e < Nv; ++e) {
1892:     refCoords[0] = refVertices[e * 2 + 0];
1893:     refCoords[1] = refVertices[e * 2 + 1];
1894:     for (PetscInt d = 1; d <= edgeDiv; ++d) {
1895:       refCoords[d * 2 + 0] = refCoords[0] + (refVertices[(e + 1) % Nv * 2 + 0] - refCoords[0]) * d / edgeDiv;
1896:       refCoords[d * 2 + 1] = refCoords[1] + (refVertices[(e + 1) % Nv * 2 + 1] - refCoords[1]) * d / edgeDiv;
1897:     }
1898:     PetscCall(DMPlexReferenceToCoordinates(dm, cell, edgeDiv + 1, refCoords, edgeCoords));
1899:     for (PetscInt d = 0; d < edgeDiv; ++d) {
1900:       PetscCall(PetscDrawTriangle(draw, centroid[0], centroid[1], edgeCoords[d * 2 + 0], edgeCoords[d * 2 + 1], edgeCoords[(d + 1) * 2 + 0], edgeCoords[(d + 1) * 2 + 1], fillColor, fillColor, fillColor));
1901:       PetscCall(PetscDrawLine(draw, edgeCoords[d * 2 + 0], edgeCoords[d * 2 + 1], edgeCoords[(d + 1) * 2 + 0], edgeCoords[(d + 1) * 2 + 1], PETSC_DRAW_BLACK));
1902:     }
1903:   }
1904:   PetscFunctionReturn(PETSC_SUCCESS);
1905: }

1907: static PetscErrorCode DMPlexDrawCellHighOrder(DM dm, PetscDraw draw, PetscInt cell, const PetscScalar coords[], PetscInt edgeDiv, PetscReal refCoords[], PetscReal edgeCoords[])
1908: {
1909:   DMPolytopeType ct;

1911:   PetscFunctionBegin;
1912:   PetscCall(DMPlexGetCellType(dm, cell, &ct));
1913:   switch (ct) {
1914:   case DM_POLYTOPE_TRIANGLE: {
1915:     PetscReal refVertices[6] = {-1., -1., 1., -1., -1., 1.};

1917:     PetscCall(DrawPolygon_Private(dm, draw, cell, 3, refVertices, coords, edgeDiv, refCoords, edgeCoords));
1918:   } break;
1919:   case DM_POLYTOPE_QUADRILATERAL: {
1920:     PetscReal refVertices[8] = {-1., -1., 1., -1., 1., 1., -1., 1.};

1922:     PetscCall(DrawPolygon_Private(dm, draw, cell, 4, refVertices, coords, edgeDiv, refCoords, edgeCoords));
1923:   } break;
1924:   default:
1925:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of type %s", DMPolytopeTypes[ct]);
1926:   }
1927:   PetscFunctionReturn(PETSC_SUCCESS);
1928: }

1930: static PetscErrorCode DMPlexView_Draw(DM dm, PetscViewer viewer)
1931: {
1932:   PetscDraw    draw;
1933:   DM           cdm;
1934:   PetscSection coordSection;
1935:   Vec          coordinates;
1936:   PetscReal    xyl[3], xyr[3];
1937:   PetscReal   *refCoords, *edgeCoords;
1938:   PetscBool    isnull, drawAffine;
1939:   PetscInt     dim, vStart, vEnd, cStart, cEnd, c, cDegree, edgeDiv;

1941:   PetscFunctionBegin;
1942:   PetscCall(DMGetCoordinateDim(dm, &dim));
1943:   PetscCheck(dim <= 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot draw meshes of dimension %" PetscInt_FMT, dim);
1944:   PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree));
1945:   drawAffine = cDegree > 1 ? PETSC_FALSE : PETSC_TRUE;
1946:   edgeDiv    = cDegree + 1;
1947:   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_view_draw_affine", &drawAffine, NULL));
1948:   if (!drawAffine) PetscCall(PetscMalloc2((edgeDiv + 1) * dim, &refCoords, (edgeDiv + 1) * dim, &edgeCoords));
1949:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1950:   PetscCall(DMGetLocalSection(cdm, &coordSection));
1951:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1952:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1953:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));

1955:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1956:   PetscCall(PetscDrawIsNull(draw, &isnull));
1957:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1958:   PetscCall(PetscDrawSetTitle(draw, "Mesh"));

1960:   PetscCall(DMGetBoundingBox(dm, xyl, xyr));
1961:   PetscCall(PetscDrawSetCoordinates(draw, xyl[0], xyl[1], xyr[0], xyr[1]));
1962:   PetscCall(PetscDrawClear(draw));

1964:   for (c = cStart; c < cEnd; ++c) {
1965:     PetscScalar       *coords = NULL;
1966:     const PetscScalar *coords_arr;
1967:     PetscInt           numCoords;
1968:     PetscBool          isDG;

1970:     PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
1971:     if (drawAffine) PetscCall(DMPlexDrawCell(dm, draw, c, coords));
1972:     else PetscCall(DMPlexDrawCellHighOrder(dm, draw, c, coords, edgeDiv, refCoords, edgeCoords));
1973:     PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
1974:   }
1975:   if (!drawAffine) PetscCall(PetscFree2(refCoords, edgeCoords));
1976:   PetscCall(PetscDrawFlush(draw));
1977:   PetscCall(PetscDrawPause(draw));
1978:   PetscCall(PetscDrawSave(draw));
1979:   PetscFunctionReturn(PETSC_SUCCESS);
1980: }

1982: static PetscErrorCode DMPlexCreateHighOrderSurrogate_Internal(DM dm, DM *hdm)
1983: {
1984:   DM           odm = dm, rdm = dm, cdm;
1985:   PetscFE      fe;
1986:   PetscSpace   sp;
1987:   PetscClassId id;
1988:   PetscInt     degree;
1989:   PetscBool    hoView = PETSC_TRUE;

1991:   PetscFunctionBegin;
1992:   PetscObjectOptionsBegin((PetscObject)dm);
1993:   PetscCall(PetscOptionsBool("-dm_plex_high_order_view", "Subsample to view meshes with high order coordinates", "DMPlexCreateHighOrderSurrogate_Internal", hoView, &hoView, NULL));
1994:   PetscOptionsEnd();
1995:   PetscCall(PetscObjectReference((PetscObject)dm));
1996:   *hdm = dm;
1997:   if (!hoView) PetscFunctionReturn(PETSC_SUCCESS);
1998:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1999:   PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
2000:   PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
2001:   if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
2002:   PetscCall(PetscFEGetBasisSpace(fe, &sp));
2003:   PetscCall(PetscSpaceGetDegree(sp, &degree, NULL));
2004:   for (PetscInt r = 0, rd = PetscCeilReal(((PetscReal)degree) / 2.); r < (PetscInt)PetscCeilReal(PetscLog2Real(degree)); ++r, rd = PetscCeilReal(((PetscReal)rd) / 2.)) {
2005:     DM  cdm, rcdm;
2006:     Mat In;
2007:     Vec cl, rcl;

2009:     PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)odm), &rdm));
2010:     PetscCall(DMPlexCreateCoordinateSpace(rdm, rd, PETSC_FALSE, NULL));
2011:     PetscCall(PetscObjectSetName((PetscObject)rdm, "Refined Mesh with Linear Coordinates"));
2012:     PetscCall(DMGetCoordinateDM(odm, &cdm));
2013:     PetscCall(DMGetCoordinateDM(rdm, &rcdm));
2014:     PetscCall(DMGetCoordinatesLocal(odm, &cl));
2015:     PetscCall(DMGetCoordinatesLocal(rdm, &rcl));
2016:     PetscCall(DMSetCoarseDM(rcdm, cdm));
2017:     PetscCall(DMCreateInterpolation(cdm, rcdm, &In, NULL));
2018:     PetscCall(MatMult(In, cl, rcl));
2019:     PetscCall(MatDestroy(&In));
2020:     PetscCall(DMSetCoordinatesLocal(rdm, rcl));
2021:     PetscCall(DMDestroy(&odm));
2022:     odm = rdm;
2023:   }
2024:   *hdm = rdm;
2025:   PetscFunctionReturn(PETSC_SUCCESS);
2026: }

2028: #if defined(PETSC_HAVE_EXODUSII)
2029:   #include <exodusII.h>
2030: #include <petscviewerexodusii.h>
2031: #endif

2033: PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer)
2034: {
2035:   PetscBool iascii, ishdf5, isvtk, isdraw, flg, isglvis, isexodus, iscgns, ispython;
2036:   char      name[PETSC_MAX_PATH_LEN];

2038:   PetscFunctionBegin;
2041:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2042:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
2043:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2044:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
2045:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
2046:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWEREXODUSII, &isexodus));
2047:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
2048:   PetscCall(PetscObjectHasFunction((PetscObject)viewer, "PetscViewerPythonViewObject_C", &ispython));
2049:   if (iascii) {
2050:     PetscViewerFormat format;
2051:     PetscCall(PetscViewerGetFormat(viewer, &format));
2052:     if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMPlexView_GLVis(dm, viewer));
2053:     else PetscCall(DMPlexView_Ascii(dm, viewer));
2054:   } else if (ishdf5) {
2055: #if defined(PETSC_HAVE_HDF5)
2056:     PetscCall(DMPlexView_HDF5_Internal(dm, viewer));
2057: #else
2058:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2059: #endif
2060:   } else if (isvtk) {
2061:     PetscCall(DMPlexVTKWriteAll((PetscObject)dm, viewer));
2062:   } else if (isdraw) {
2063:     DM hdm;

2065:     PetscCall(DMPlexCreateHighOrderSurrogate_Internal(dm, &hdm));
2066:     PetscCall(DMPlexView_Draw(hdm, viewer));
2067:     PetscCall(DMDestroy(&hdm));
2068:   } else if (isglvis) {
2069:     PetscCall(DMPlexView_GLVis(dm, viewer));
2070: #if defined(PETSC_HAVE_EXODUSII)
2071:   } else if (isexodus) {
2072:     /*
2073:       exodusII requires that all sets be part of exactly one cell set.
2074:       If the dm does not have a "Cell Sets" label defined, we create one
2075:       with ID 1, containing all cells.
2076:       Note that if the Cell Sets label is defined but does not cover all cells,
2077:       we may still have a problem. This should probably be checked here or in the viewer;
2078:     */
2079:     PetscInt numCS;
2080:     PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS));
2081:     if (!numCS) {
2082:       PetscInt cStart, cEnd, c;
2083:       PetscCall(DMCreateLabel(dm, "Cell Sets"));
2084:       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2085:       for (c = cStart; c < cEnd; ++c) PetscCall(DMSetLabelValue(dm, "Cell Sets", c, 1));
2086:     }
2087:     PetscCall(DMView_PlexExodusII(dm, viewer));
2088: #endif
2089: #if defined(PETSC_HAVE_CGNS)
2090:   } else if (iscgns) {
2091:     PetscCall(DMView_PlexCGNS(dm, viewer));
2092: #endif
2093:   } else if (ispython) {
2094:     PetscCall(PetscViewerPythonViewObject(viewer, (PetscObject)dm));
2095:   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlex writing", ((PetscObject)viewer)->type_name);
2096:   /* Optionally view the partition */
2097:   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_partition_view", &flg));
2098:   if (flg) {
2099:     Vec ranks;
2100:     PetscCall(DMPlexCreateRankField(dm, &ranks));
2101:     PetscCall(VecView(ranks, viewer));
2102:     PetscCall(VecDestroy(&ranks));
2103:   }
2104:   /* Optionally view a label */
2105:   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_label_view", name, sizeof(name), &flg));
2106:   if (flg) {
2107:     DMLabel label;
2108:     Vec     val;

2110:     PetscCall(DMGetLabel(dm, name, &label));
2111:     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Label %s provided to -dm_label_view does not exist in this DM", name);
2112:     PetscCall(DMPlexCreateLabelField(dm, label, &val));
2113:     PetscCall(VecView(val, viewer));
2114:     PetscCall(VecDestroy(&val));
2115:   }
2116:   PetscFunctionReturn(PETSC_SUCCESS);
2117: }

2119: /*@
2120:   DMPlexTopologyView - Saves a `DMPLEX` topology into a file

2122:   Collective

2124:   Input Parameters:
2125: + dm     - The `DM` whose topology is to be saved
2126: - viewer - The `PetscViewer` to save it in

2128:   Level: advanced

2130: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexCoordinatesView()`, `DMPlexLabelsView()`, `DMPlexTopologyLoad()`, `PetscViewer`
2131: @*/
2132: PetscErrorCode DMPlexTopologyView(DM dm, PetscViewer viewer)
2133: {
2134:   PetscBool ishdf5;

2136:   PetscFunctionBegin;
2139:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2140:   PetscCall(PetscLogEventBegin(DMPLEX_TopologyView, viewer, 0, 0, 0));
2141:   if (ishdf5) {
2142: #if defined(PETSC_HAVE_HDF5)
2143:     PetscViewerFormat format;
2144:     PetscCall(PetscViewerGetFormat(viewer, &format));
2145:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2146:       IS globalPointNumbering;

2148:       PetscCall(DMPlexCreatePointNumbering(dm, &globalPointNumbering));
2149:       PetscCall(DMPlexTopologyView_HDF5_Internal(dm, globalPointNumbering, viewer));
2150:       PetscCall(ISDestroy(&globalPointNumbering));
2151:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
2152: #else
2153:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2154: #endif
2155:   }
2156:   PetscCall(PetscLogEventEnd(DMPLEX_TopologyView, viewer, 0, 0, 0));
2157:   PetscFunctionReturn(PETSC_SUCCESS);
2158: }

2160: /*@
2161:   DMPlexCoordinatesView - Saves `DMPLEX` coordinates into a file

2163:   Collective

2165:   Input Parameters:
2166: + dm     - The `DM` whose coordinates are to be saved
2167: - viewer - The `PetscViewer` for saving

2169:   Level: advanced

2171: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexTopologyView()`, `DMPlexLabelsView()`, `DMPlexCoordinatesLoad()`, `PetscViewer`
2172: @*/
2173: PetscErrorCode DMPlexCoordinatesView(DM dm, PetscViewer viewer)
2174: {
2175:   PetscBool ishdf5;

2177:   PetscFunctionBegin;
2180:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2181:   PetscCall(PetscLogEventBegin(DMPLEX_CoordinatesView, viewer, 0, 0, 0));
2182:   if (ishdf5) {
2183: #if defined(PETSC_HAVE_HDF5)
2184:     PetscViewerFormat format;
2185:     PetscCall(PetscViewerGetFormat(viewer, &format));
2186:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2187:       PetscCall(DMPlexCoordinatesView_HDF5_Internal(dm, viewer));
2188:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
2189: #else
2190:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2191: #endif
2192:   }
2193:   PetscCall(PetscLogEventEnd(DMPLEX_CoordinatesView, viewer, 0, 0, 0));
2194:   PetscFunctionReturn(PETSC_SUCCESS);
2195: }

2197: /*@
2198:   DMPlexLabelsView - Saves `DMPLEX` labels into a file

2200:   Collective

2202:   Input Parameters:
2203: + dm     - The `DM` whose labels are to be saved
2204: - viewer - The `PetscViewer` for saving

2206:   Level: advanced

2208: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexTopologyView()`, `DMPlexCoordinatesView()`, `DMPlexLabelsLoad()`, `PetscViewer`
2209: @*/
2210: PetscErrorCode DMPlexLabelsView(DM dm, PetscViewer viewer)
2211: {
2212:   PetscBool ishdf5;

2214:   PetscFunctionBegin;
2217:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2218:   PetscCall(PetscLogEventBegin(DMPLEX_LabelsView, viewer, 0, 0, 0));
2219:   if (ishdf5) {
2220: #if defined(PETSC_HAVE_HDF5)
2221:     IS                globalPointNumbering;
2222:     PetscViewerFormat format;

2224:     PetscCall(PetscViewerGetFormat(viewer, &format));
2225:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2226:       PetscCall(DMPlexCreatePointNumbering(dm, &globalPointNumbering));
2227:       PetscCall(DMPlexLabelsView_HDF5_Internal(dm, globalPointNumbering, viewer));
2228:       PetscCall(ISDestroy(&globalPointNumbering));
2229:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2230: #else
2231:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2232: #endif
2233:   }
2234:   PetscCall(PetscLogEventEnd(DMPLEX_LabelsView, viewer, 0, 0, 0));
2235:   PetscFunctionReturn(PETSC_SUCCESS);
2236: }

2238: /*@
2239:   DMPlexSectionView - Saves a section associated with a `DMPLEX`

2241:   Collective

2243:   Input Parameters:
2244: + dm        - The `DM` that contains the topology on which the section to be saved is defined
2245: . viewer    - The `PetscViewer` for saving
2246: - sectiondm - The `DM` that contains the section to be saved, can be `NULL`

2248:   Level: advanced

2250:   Notes:
2251:   This function is a wrapper around `PetscSectionView()`; in addition to the raw section, it saves information that associates the section points to the topology (`dm`) points. When the topology (`dm`) and the section are later loaded with `DMPlexTopologyLoad()` and `DMPlexSectionLoad()`, respectively, this information is used to match section points with topology points.

2253:   In general `dm` and `sectiondm` are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object (or in case `sectiondm` is `NULL`) if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.

2255: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexTopologyView()`, `DMPlexCoordinatesView()`, `DMPlexLabelsView()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`, `PetscSectionView()`, `DMPlexSectionLoad()`, `PetscViewer`
2256: @*/
2257: PetscErrorCode DMPlexSectionView(DM dm, PetscViewer viewer, DM sectiondm)
2258: {
2259:   PetscBool ishdf5;

2261:   PetscFunctionBegin;
2264:   if (!sectiondm) sectiondm = dm;
2266:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2267:   PetscCall(PetscLogEventBegin(DMPLEX_SectionView, viewer, 0, 0, 0));
2268:   if (ishdf5) {
2269: #if defined(PETSC_HAVE_HDF5)
2270:     PetscCall(DMPlexSectionView_HDF5_Internal(dm, viewer, sectiondm));
2271: #else
2272:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2273: #endif
2274:   }
2275:   PetscCall(PetscLogEventEnd(DMPLEX_SectionView, viewer, 0, 0, 0));
2276:   PetscFunctionReturn(PETSC_SUCCESS);
2277: }

2279: /*@
2280:   DMPlexGlobalVectorView - Saves a global vector

2282:   Collective

2284:   Input Parameters:
2285: + dm        - The `DM` that represents the topology
2286: . viewer    - The `PetscViewer` to save data with
2287: . sectiondm - The `DM` that contains the global section on which vec is defined, can be `NULL`
2288: - vec       - The global vector to be saved

2290:   Level: advanced

2292:   Notes:
2293:   In general `dm` and `sectiondm` are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object (or in case `sectiondm` is `NULL`) if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.

2295:   Calling sequence:
2296: .vb
2297:        DMCreate(PETSC_COMM_WORLD, &dm);
2298:        DMSetType(dm, DMPLEX);
2299:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2300:        DMClone(dm, &sectiondm);
2301:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2302:        PetscSectionCreate(PETSC_COMM_WORLD, &section);
2303:        DMPlexGetChart(sectiondm, &pStart, &pEnd);
2304:        PetscSectionSetChart(section, pStart, pEnd);
2305:        PetscSectionSetUp(section);
2306:        DMSetLocalSection(sectiondm, section);
2307:        PetscSectionDestroy(&section);
2308:        DMGetGlobalVector(sectiondm, &vec);
2309:        PetscObjectSetName((PetscObject)vec, "vec_name");
2310:        DMPlexTopologyView(dm, viewer);
2311:        DMPlexSectionView(dm, viewer, sectiondm);
2312:        DMPlexGlobalVectorView(dm, viewer, sectiondm, vec);
2313:        DMRestoreGlobalVector(sectiondm, &vec);
2314:        DMDestroy(&sectiondm);
2315:        DMDestroy(&dm);
2316: .ve

2318: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyView()`, `DMPlexSectionView()`, `DMPlexLocalVectorView()`, `DMPlexGlobalVectorLoad()`, `DMPlexLocalVectorLoad()`
2319: @*/
2320: PetscErrorCode DMPlexGlobalVectorView(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
2321: {
2322:   PetscBool ishdf5;

2324:   PetscFunctionBegin;
2327:   if (!sectiondm) sectiondm = dm;
2330:   /* Check consistency */
2331:   {
2332:     PetscSection section;
2333:     PetscBool    includesConstraints;
2334:     PetscInt     m, m1;

2336:     PetscCall(VecGetLocalSize(vec, &m1));
2337:     PetscCall(DMGetGlobalSection(sectiondm, &section));
2338:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2339:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2340:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2341:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global vector size (%" PetscInt_FMT ") != global section storage size (%" PetscInt_FMT ")", m1, m);
2342:   }
2343:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2344:   PetscCall(PetscLogEventBegin(DMPLEX_GlobalVectorView, viewer, 0, 0, 0));
2345:   if (ishdf5) {
2346: #if defined(PETSC_HAVE_HDF5)
2347:     PetscCall(DMPlexGlobalVectorView_HDF5_Internal(dm, viewer, sectiondm, vec));
2348: #else
2349:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2350: #endif
2351:   }
2352:   PetscCall(PetscLogEventEnd(DMPLEX_GlobalVectorView, viewer, 0, 0, 0));
2353:   PetscFunctionReturn(PETSC_SUCCESS);
2354: }

2356: /*@
2357:   DMPlexLocalVectorView - Saves a local vector

2359:   Collective

2361:   Input Parameters:
2362: + dm        - The `DM` that represents the topology
2363: . viewer    - The `PetscViewer` to save data with
2364: . sectiondm - The `DM` that contains the local section on which `vec` is defined, can be `NULL`
2365: - vec       - The local vector to be saved

2367:   Level: advanced

2369:   Note:
2370:   In general `dm` and `sectiondm` are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object (or in case `sectiondm` is `NULL`) if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.

2372:   Calling sequence:
2373: .vb
2374:        DMCreate(PETSC_COMM_WORLD, &dm);
2375:        DMSetType(dm, DMPLEX);
2376:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2377:        DMClone(dm, &sectiondm);
2378:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2379:        PetscSectionCreate(PETSC_COMM_WORLD, &section);
2380:        DMPlexGetChart(sectiondm, &pStart, &pEnd);
2381:        PetscSectionSetChart(section, pStart, pEnd);
2382:        PetscSectionSetUp(section);
2383:        DMSetLocalSection(sectiondm, section);
2384:        DMGetLocalVector(sectiondm, &vec);
2385:        PetscObjectSetName((PetscObject)vec, "vec_name");
2386:        DMPlexTopologyView(dm, viewer);
2387:        DMPlexSectionView(dm, viewer, sectiondm);
2388:        DMPlexLocalVectorView(dm, viewer, sectiondm, vec);
2389:        DMRestoreLocalVector(sectiondm, &vec);
2390:        DMDestroy(&sectiondm);
2391:        DMDestroy(&dm);
2392: .ve

2394: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyView()`, `DMPlexSectionView()`, `DMPlexGlobalVectorView()`, `DMPlexGlobalVectorLoad()`, `DMPlexLocalVectorLoad()`
2395: @*/
2396: PetscErrorCode DMPlexLocalVectorView(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
2397: {
2398:   PetscBool ishdf5;

2400:   PetscFunctionBegin;
2403:   if (!sectiondm) sectiondm = dm;
2406:   /* Check consistency */
2407:   {
2408:     PetscSection section;
2409:     PetscBool    includesConstraints;
2410:     PetscInt     m, m1;

2412:     PetscCall(VecGetLocalSize(vec, &m1));
2413:     PetscCall(DMGetLocalSection(sectiondm, &section));
2414:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2415:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2416:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2417:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local vector size (%" PetscInt_FMT ") != local section storage size (%" PetscInt_FMT ")", m1, m);
2418:   }
2419:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2420:   PetscCall(PetscLogEventBegin(DMPLEX_LocalVectorView, viewer, 0, 0, 0));
2421:   if (ishdf5) {
2422: #if defined(PETSC_HAVE_HDF5)
2423:     PetscCall(DMPlexLocalVectorView_HDF5_Internal(dm, viewer, sectiondm, vec));
2424: #else
2425:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2426: #endif
2427:   }
2428:   PetscCall(PetscLogEventEnd(DMPLEX_LocalVectorView, viewer, 0, 0, 0));
2429:   PetscFunctionReturn(PETSC_SUCCESS);
2430: }

2432: PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer)
2433: {
2434:   PetscBool ishdf5;

2436:   PetscFunctionBegin;
2439:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2440:   if (ishdf5) {
2441: #if defined(PETSC_HAVE_HDF5)
2442:     PetscViewerFormat format;
2443:     PetscCall(PetscViewerGetFormat(viewer, &format));
2444:     if (format == PETSC_VIEWER_HDF5_XDMF || format == PETSC_VIEWER_HDF5_VIZ) {
2445:       PetscCall(DMPlexLoad_HDF5_Xdmf_Internal(dm, viewer));
2446:     } else if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2447:       PetscCall(DMPlexLoad_HDF5_Internal(dm, viewer));
2448:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2449:     PetscFunctionReturn(PETSC_SUCCESS);
2450: #else
2451:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2452: #endif
2453:   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlex loading", ((PetscObject)viewer)->type_name);
2454: }

2456: /*@
2457:   DMPlexTopologyLoad - Loads a topology into a `DMPLEX`

2459:   Collective

2461:   Input Parameters:
2462: + dm     - The `DM` into which the topology is loaded
2463: - viewer - The `PetscViewer` for the saved topology

2465:   Output Parameter:
2466: . globalToLocalPointSF - The `PetscSF` that pushes points in [0, N) to the associated points in the loaded `DMPLEX`, where N is the global number of points;
2467:   `NULL` if unneeded

2469:   Level: advanced

2471: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexCoordinatesLoad()`, `DMPlexLabelsLoad()`, `DMView()`, `PetscViewerHDF5Open()`, `PetscViewerPushFormat()`,
2472:           `PetscViewer`, `PetscSF`
2473: @*/
2474: PetscErrorCode DMPlexTopologyLoad(DM dm, PetscViewer viewer, PetscSF *globalToLocalPointSF)
2475: {
2476:   PetscBool ishdf5;

2478:   PetscFunctionBegin;
2481:   if (globalToLocalPointSF) PetscAssertPointer(globalToLocalPointSF, 3);
2482:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2483:   PetscCall(PetscLogEventBegin(DMPLEX_TopologyLoad, viewer, 0, 0, 0));
2484:   if (ishdf5) {
2485: #if defined(PETSC_HAVE_HDF5)
2486:     PetscViewerFormat format;
2487:     PetscCall(PetscViewerGetFormat(viewer, &format));
2488:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2489:       PetscCall(DMPlexTopologyLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF));
2490:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2491: #else
2492:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2493: #endif
2494:   }
2495:   PetscCall(PetscLogEventEnd(DMPLEX_TopologyLoad, viewer, 0, 0, 0));
2496:   PetscFunctionReturn(PETSC_SUCCESS);
2497: }

2499: /*@
2500:   DMPlexCoordinatesLoad - Loads coordinates into a `DMPLEX`

2502:   Collective

2504:   Input Parameters:
2505: + dm                   - The `DM` into which the coordinates are loaded
2506: . viewer               - The `PetscViewer` for the saved coordinates
2507: - globalToLocalPointSF - The `PetscSF` returned by `DMPlexTopologyLoad()` when loading dm from viewer

2509:   Level: advanced

2511: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexTopologyLoad()`, `DMPlexLabelsLoad()`, `DMView()`, `PetscViewerHDF5Open()`, `PetscViewerPushFormat()`,
2512:           `PetscSF`, `PetscViewer`
2513: @*/
2514: PetscErrorCode DMPlexCoordinatesLoad(DM dm, PetscViewer viewer, PetscSF globalToLocalPointSF)
2515: {
2516:   PetscBool ishdf5;

2518:   PetscFunctionBegin;
2522:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2523:   PetscCall(PetscLogEventBegin(DMPLEX_CoordinatesLoad, viewer, 0, 0, 0));
2524:   if (ishdf5) {
2525: #if defined(PETSC_HAVE_HDF5)
2526:     PetscViewerFormat format;
2527:     PetscCall(PetscViewerGetFormat(viewer, &format));
2528:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2529:       PetscCall(DMPlexCoordinatesLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF));
2530:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2531: #else
2532:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2533: #endif
2534:   }
2535:   PetscCall(PetscLogEventEnd(DMPLEX_CoordinatesLoad, viewer, 0, 0, 0));
2536:   PetscFunctionReturn(PETSC_SUCCESS);
2537: }

2539: /*@
2540:   DMPlexLabelsLoad - Loads labels into a `DMPLEX`

2542:   Collective

2544:   Input Parameters:
2545: + dm                   - The `DM` into which the labels are loaded
2546: . viewer               - The `PetscViewer` for the saved labels
2547: - globalToLocalPointSF - The `PetscSF` returned by `DMPlexTopologyLoad()` when loading `dm` from viewer

2549:   Level: advanced

2551:   Note:
2552:   The `PetscSF` argument must not be `NULL` if the `DM` is distributed, otherwise an error occurs.

2554: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexTopologyLoad()`, `DMPlexCoordinatesLoad()`, `DMView()`, `PetscViewerHDF5Open()`, `PetscViewerPushFormat()`,
2555:           `PetscSF`, `PetscViewer`
2556: @*/
2557: PetscErrorCode DMPlexLabelsLoad(DM dm, PetscViewer viewer, PetscSF globalToLocalPointSF)
2558: {
2559:   PetscBool ishdf5;

2561:   PetscFunctionBegin;
2565:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2566:   PetscCall(PetscLogEventBegin(DMPLEX_LabelsLoad, viewer, 0, 0, 0));
2567:   if (ishdf5) {
2568: #if defined(PETSC_HAVE_HDF5)
2569:     PetscViewerFormat format;

2571:     PetscCall(PetscViewerGetFormat(viewer, &format));
2572:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2573:       PetscCall(DMPlexLabelsLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF));
2574:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2575: #else
2576:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2577: #endif
2578:   }
2579:   PetscCall(PetscLogEventEnd(DMPLEX_LabelsLoad, viewer, 0, 0, 0));
2580:   PetscFunctionReturn(PETSC_SUCCESS);
2581: }

2583: /*@
2584:   DMPlexSectionLoad - Loads section into a `DMPLEX`

2586:   Collective

2588:   Input Parameters:
2589: + dm                   - The `DM` that represents the topology
2590: . viewer               - The `PetscViewer` that represents the on-disk section (sectionA)
2591: . sectiondm            - The `DM` into which the on-disk section (sectionA) is migrated, can be `NULL`
2592: - globalToLocalPointSF - The `PetscSF` returned by `DMPlexTopologyLoad(`) when loading dm from viewer

2594:   Output Parameters:
2595: + globalDofSF - The `PetscSF` that migrates any on-disk `Vec` data associated with sectionA into a global `Vec` associated with the `sectiondm`'s global section (`NULL` if not needed)
2596: - localDofSF  - The `PetscSF` that migrates any on-disk `Vec` data associated with sectionA into a local `Vec` associated with the `sectiondm`'s local section (`NULL` if not needed)

2598:   Level: advanced

2600:   Notes:
2601:   This function is a wrapper around `PetscSectionLoad()`; it loads, in addition to the raw section, a list of global point numbers that associates each on-disk section point with a global point number in [0, NX), where NX is the number of topology points in `dm`. Noting that globalToLocalPointSF associates each topology point in dm with a global number in [0, NX), one can readily establish an association of the on-disk section points with the topology points.

2603:   In general `dm` and `sectiondm` are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object (or in case `sectiondm` is `NULL`) if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.

2605:   The output parameter, `globalDofSF` (`localDofSF`), can later be used with `DMPlexGlobalVectorLoad()` (`DMPlexLocalVectorLoad()`) to load on-disk vectors into global (local) vectors associated with sectiondm's global (local) section.

2607:   Example using 2 processes:
2608: .vb
2609:   NX (number of points on dm): 4
2610:   sectionA                   : the on-disk section
2611:   vecA                       : a vector associated with sectionA
2612:   sectionB                   : sectiondm's local section constructed in this function
2613:   vecB (local)               : a vector associated with sectiondm's local section
2614:   vecB (global)              : a vector associated with sectiondm's global section

2616:                                      rank 0    rank 1
2617:   vecA (global)                  : [.0 .4 .1 | .2 .3]        <- to be loaded in DMPlexGlobalVectorLoad() or DMPlexLocalVectorLoad()
2618:   sectionA->atlasOff             :       0 2 | 1             <- loaded in PetscSectionLoad()
2619:   sectionA->atlasDof             :       1 3 | 1             <- loaded in PetscSectionLoad()
2620:   sectionA's global point numbers:       0 2 | 3             <- loaded in DMPlexSectionLoad()
2621:   [0, NX)                        :       0 1 | 2 3           <- conceptual partition used in globalToLocalPointSF
2622:   sectionB's global point numbers:     0 1 3 | 3 2           <- associated with [0, NX) by globalToLocalPointSF
2623:   sectionB->atlasDof             :     1 0 1 | 1 3
2624:   sectionB->atlasOff (no perm)   :     0 1 1 | 0 1
2625:   vecB (local)                   :   [.0 .4] | [.4 .1 .2 .3] <- to be constructed by calling DMPlexLocalVectorLoad() with localDofSF
2626:   vecB (global)                  :    [.0 .4 | .1 .2 .3]     <- to be constructed by calling DMPlexGlobalVectorLoad() with globalDofSF
2627: .ve
2628:   where "|" represents a partition of loaded data, and global point 3 is assumed to be owned by rank 0.

2630: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexTopologyLoad()`, `DMPlexCoordinatesLoad()`, `DMPlexLabelsLoad()`, `DMPlexGlobalVectorLoad()`, `DMPlexLocalVectorLoad()`, `PetscSectionLoad()`, `DMPlexSectionView()`, `PetscSF`, `PetscViewer`
2631: @*/
2632: PetscErrorCode DMPlexSectionLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF globalToLocalPointSF, PetscSF *globalDofSF, PetscSF *localDofSF)
2633: {
2634:   PetscBool ishdf5;

2636:   PetscFunctionBegin;
2639:   if (!sectiondm) sectiondm = dm;
2642:   if (globalDofSF) PetscAssertPointer(globalDofSF, 5);
2643:   if (localDofSF) PetscAssertPointer(localDofSF, 6);
2644:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2645:   PetscCall(PetscLogEventBegin(DMPLEX_SectionLoad, viewer, 0, 0, 0));
2646:   if (ishdf5) {
2647: #if defined(PETSC_HAVE_HDF5)
2648:     PetscCall(DMPlexSectionLoad_HDF5_Internal(dm, viewer, sectiondm, globalToLocalPointSF, globalDofSF, localDofSF));
2649: #else
2650:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2651: #endif
2652:   }
2653:   PetscCall(PetscLogEventEnd(DMPLEX_SectionLoad, viewer, 0, 0, 0));
2654:   PetscFunctionReturn(PETSC_SUCCESS);
2655: }

2657: /*@
2658:   DMPlexGlobalVectorLoad - Loads on-disk vector data into a global vector

2660:   Collective

2662:   Input Parameters:
2663: + dm        - The `DM` that represents the topology
2664: . viewer    - The `PetscViewer` that represents the on-disk vector data
2665: . sectiondm - The `DM` that contains the global section on which vec is defined, can be `NULL`
2666: . sf        - The `PetscSF` that migrates the on-disk vector data into vec
2667: - vec       - The global vector to set values of

2669:   Level: advanced

2671:   Notes:
2672:   In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object (or in case `sectiondm` is `NULL`) if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.

2674:   Calling sequence:
2675: .vb
2676:        DMCreate(PETSC_COMM_WORLD, &dm);
2677:        DMSetType(dm, DMPLEX);
2678:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2679:        DMPlexTopologyLoad(dm, viewer, &sfX);
2680:        DMClone(dm, &sectiondm);
2681:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2682:        DMPlexSectionLoad(dm, viewer, sectiondm, sfX, &gsf, NULL);
2683:        DMGetGlobalVector(sectiondm, &vec);
2684:        PetscObjectSetName((PetscObject)vec, "vec_name");
2685:        DMPlexGlobalVectorLoad(dm, viewer, sectiondm, gsf, vec);
2686:        DMRestoreGlobalVector(sectiondm, &vec);
2687:        PetscSFDestroy(&gsf);
2688:        PetscSFDestroy(&sfX);
2689:        DMDestroy(&sectiondm);
2690:        DMDestroy(&dm);
2691: .ve

2693: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyLoad()`, `DMPlexSectionLoad()`, `DMPlexLocalVectorLoad()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`,
2694:           `PetscSF`, `PetscViewer`
2695: @*/
2696: PetscErrorCode DMPlexGlobalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2697: {
2698:   PetscBool ishdf5;

2700:   PetscFunctionBegin;
2703:   if (!sectiondm) sectiondm = dm;
2707:   /* Check consistency */
2708:   {
2709:     PetscSection section;
2710:     PetscBool    includesConstraints;
2711:     PetscInt     m, m1;

2713:     PetscCall(VecGetLocalSize(vec, &m1));
2714:     PetscCall(DMGetGlobalSection(sectiondm, &section));
2715:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2716:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2717:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2718:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global vector size (%" PetscInt_FMT ") != global section storage size (%" PetscInt_FMT ")", m1, m);
2719:   }
2720:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2721:   PetscCall(PetscLogEventBegin(DMPLEX_GlobalVectorLoad, viewer, 0, 0, 0));
2722:   if (ishdf5) {
2723: #if defined(PETSC_HAVE_HDF5)
2724:     PetscCall(DMPlexVecLoad_HDF5_Internal(dm, viewer, sectiondm, sf, vec));
2725: #else
2726:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2727: #endif
2728:   }
2729:   PetscCall(PetscLogEventEnd(DMPLEX_GlobalVectorLoad, viewer, 0, 0, 0));
2730:   PetscFunctionReturn(PETSC_SUCCESS);
2731: }

2733: /*@
2734:   DMPlexLocalVectorLoad - Loads on-disk vector data into a local vector

2736:   Collective

2738:   Input Parameters:
2739: + dm        - The `DM` that represents the topology
2740: . viewer    - The `PetscViewer` that represents the on-disk vector data
2741: . sectiondm - The `DM` that contains the local section on which vec is defined, can be `NULL`
2742: . sf        - The `PetscSF` that migrates the on-disk vector data into vec
2743: - vec       - The local vector to set values of

2745:   Level: advanced

2747:   Notes:
2748:   In general `dm` and `sectiondm` are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with `PetscObjectSetName()`. In practice, however, they can be the same object (or in case `sectiondm` is `NULL`) if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.

2750:   Calling sequence:
2751: .vb
2752:        DMCreate(PETSC_COMM_WORLD, &dm);
2753:        DMSetType(dm, DMPLEX);
2754:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2755:        DMPlexTopologyLoad(dm, viewer, &sfX);
2756:        DMClone(dm, &sectiondm);
2757:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2758:        DMPlexSectionLoad(dm, viewer, sectiondm, sfX, NULL, &lsf);
2759:        DMGetLocalVector(sectiondm, &vec);
2760:        PetscObjectSetName((PetscObject)vec, "vec_name");
2761:        DMPlexLocalVectorLoad(dm, viewer, sectiondm, lsf, vec);
2762:        DMRestoreLocalVector(sectiondm, &vec);
2763:        PetscSFDestroy(&lsf);
2764:        PetscSFDestroy(&sfX);
2765:        DMDestroy(&sectiondm);
2766:        DMDestroy(&dm);
2767: .ve

2769: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyLoad()`, `DMPlexSectionLoad()`, `DMPlexGlobalVectorLoad()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`,
2770:           `PetscSF`, `PetscViewer`
2771: @*/
2772: PetscErrorCode DMPlexLocalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2773: {
2774:   PetscBool ishdf5;

2776:   PetscFunctionBegin;
2779:   if (!sectiondm) sectiondm = dm;
2783:   /* Check consistency */
2784:   {
2785:     PetscSection section;
2786:     PetscBool    includesConstraints;
2787:     PetscInt     m, m1;

2789:     PetscCall(VecGetLocalSize(vec, &m1));
2790:     PetscCall(DMGetLocalSection(sectiondm, &section));
2791:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2792:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2793:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2794:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local vector size (%" PetscInt_FMT ") != local section storage size (%" PetscInt_FMT ")", m1, m);
2795:   }
2796:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2797:   PetscCall(PetscLogEventBegin(DMPLEX_LocalVectorLoad, viewer, 0, 0, 0));
2798:   if (ishdf5) {
2799: #if defined(PETSC_HAVE_HDF5)
2800:     PetscCall(DMPlexVecLoad_HDF5_Internal(dm, viewer, sectiondm, sf, vec));
2801: #else
2802:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2803: #endif
2804:   }
2805:   PetscCall(PetscLogEventEnd(DMPLEX_LocalVectorLoad, viewer, 0, 0, 0));
2806:   PetscFunctionReturn(PETSC_SUCCESS);
2807: }

2809: PetscErrorCode DMDestroy_Plex(DM dm)
2810: {
2811:   DM_Plex *mesh = (DM_Plex *)dm->data;

2813:   PetscFunctionBegin;
2814:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
2815:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", NULL));
2816:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", NULL));
2817:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMInterpolateSolution_C", NULL));
2818:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertTimeDerivativeBoundaryValues_C", NULL));
2819:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", NULL));
2820:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeGetDefault_C", NULL));
2821:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeSetDefault_C", NULL));
2822:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", NULL));
2823:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderGetDefault_C", NULL));
2824:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderSetDefault_C", NULL));
2825:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionGetDefault_C", NULL));
2826:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionSetDefault_C", NULL));
2827:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionGetType_C", NULL));
2828:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionSetType_C", NULL));
2829:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", NULL));
2830:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetOverlap_C", NULL));
2831:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetUseCeed_C", NULL));
2832:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetUseCeed_C", NULL));
2833:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", NULL));
2834:   if (--mesh->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
2835:   PetscCall(PetscSectionDestroy(&mesh->coneSection));
2836:   PetscCall(PetscFree(mesh->cones));
2837:   PetscCall(PetscFree(mesh->coneOrientations));
2838:   PetscCall(PetscSectionDestroy(&mesh->supportSection));
2839:   PetscCall(PetscSectionDestroy(&mesh->subdomainSection));
2840:   PetscCall(PetscFree(mesh->supports));
2841:   PetscCall(PetscFree(mesh->cellTypes));
2842:   PetscCall(DMPlexTransformDestroy(&mesh->tr));
2843:   PetscCall(PetscFree(mesh->tetgenOpts));
2844:   PetscCall(PetscFree(mesh->triangleOpts));
2845:   PetscCall(PetscFree(mesh->transformType));
2846:   PetscCall(PetscFree(mesh->distributionName));
2847:   PetscCall(PetscPartitionerDestroy(&mesh->partitioner));
2848:   PetscCall(DMLabelDestroy(&mesh->subpointMap));
2849:   PetscCall(ISDestroy(&mesh->subpointIS));
2850:   PetscCall(ISDestroy(&mesh->globalVertexNumbers));
2851:   PetscCall(ISDestroy(&mesh->globalCellNumbers));
2852:   if (mesh->periodic.face_sfs) {
2853:     for (PetscInt i = 0; i < mesh->periodic.num_face_sfs; i++) PetscCall(PetscSFDestroy(&mesh->periodic.face_sfs[i]));
2854:     PetscCall(PetscFree(mesh->periodic.face_sfs));
2855:   }
2856:   PetscCall(PetscSFDestroy(&mesh->periodic.composed_sf));
2857:   if (mesh->periodic.periodic_points) {
2858:     for (PetscInt i = 0; i < mesh->periodic.num_face_sfs; i++) PetscCall(ISDestroy(&mesh->periodic.periodic_points[i]));
2859:     PetscCall(PetscFree(mesh->periodic.periodic_points));
2860:   }
2861:   if (mesh->periodic.transform) PetscCall(PetscFree(mesh->periodic.transform));
2862:   PetscCall(PetscSectionDestroy(&mesh->anchorSection));
2863:   PetscCall(ISDestroy(&mesh->anchorIS));
2864:   PetscCall(PetscSectionDestroy(&mesh->parentSection));
2865:   PetscCall(PetscFree(mesh->parents));
2866:   PetscCall(PetscFree(mesh->childIDs));
2867:   PetscCall(PetscSectionDestroy(&mesh->childSection));
2868:   PetscCall(PetscFree(mesh->children));
2869:   PetscCall(DMDestroy(&mesh->referenceTree));
2870:   PetscCall(PetscGridHashDestroy(&mesh->lbox));
2871:   PetscCall(PetscFree(mesh->neighbors));
2872:   if (mesh->metricCtx) PetscCall(PetscFree(mesh->metricCtx));
2873:   if (mesh->nonempty_comm != MPI_COMM_NULL && mesh->nonempty_comm != MPI_COMM_SELF) PetscCallMPI(MPI_Comm_free(&mesh->nonempty_comm));
2874:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
2875:   PetscCall(PetscFree(mesh));
2876:   PetscFunctionReturn(PETSC_SUCCESS);
2877: }

2879: PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J)
2880: {
2881:   PetscSection           sectionGlobal, sectionLocal;
2882:   PetscInt               bs = -1, mbs;
2883:   PetscInt               localSize, localStart = 0;
2884:   PetscBool              isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isMatIS;
2885:   MatType                mtype;
2886:   ISLocalToGlobalMapping ltog;

2888:   PetscFunctionBegin;
2889:   PetscCall(MatInitializePackage());
2890:   mtype = dm->mattype;
2891:   PetscCall(DMGetLocalSection(dm, &sectionLocal));
2892:   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
2893:   /* PetscCall(PetscSectionGetStorageSize(sectionGlobal, &localSize)); */
2894:   PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize));
2895:   PetscCallMPI(MPI_Exscan(&localSize, &localStart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
2896:   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
2897:   PetscCall(MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE));
2898:   PetscCall(MatSetType(*J, mtype));
2899:   PetscCall(MatSetFromOptions(*J));
2900:   PetscCall(MatGetBlockSize(*J, &mbs));
2901:   if (mbs > 1) bs = mbs;
2902:   PetscCall(PetscStrcmp(mtype, MATSHELL, &isShell));
2903:   PetscCall(PetscStrcmp(mtype, MATBAIJ, &isBlock));
2904:   PetscCall(PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock));
2905:   PetscCall(PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock));
2906:   PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
2907:   PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
2908:   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
2909:   PetscCall(PetscStrcmp(mtype, MATIS, &isMatIS));
2910:   if (!isShell) {
2911:     // There are three states with pblocks, since block starts can have no dofs:
2912:     // UNKNOWN) New Block:   An open block has been signalled by pblocks[p] == 1
2913:     // TRUE)    Block Start: The first entry in a block has been added
2914:     // FALSE)   Block Add:   An additional block entry has been added, since pblocks[p] == 0
2915:     PetscBT         blst;
2916:     PetscBool3      bstate     = PETSC_BOOL3_UNKNOWN;
2917:     PetscBool       fillMatrix = (PetscBool)(!dm->prealloc_only && !isMatIS);
2918:     const PetscInt *perm       = NULL;
2919:     PetscInt       *dnz, *onz, *dnzu, *onzu, bsLocal[2], bsMinMax[2], *pblocks;
2920:     PetscInt        pStart, pEnd, dof, cdof, num_fields;

2922:     PetscCall(DMGetLocalToGlobalMapping(dm, &ltog));
2923:     PetscCall(PetscSectionGetBlockStarts(sectionLocal, &blst));
2924:     if (sectionLocal->perm) PetscCall(ISGetIndices(sectionLocal->perm, &perm));

2926:     PetscCall(PetscCalloc1(localSize, &pblocks));
2927:     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
2928:     PetscCall(PetscSectionGetNumFields(sectionGlobal, &num_fields));
2929:     // We need to process in the permuted order to get block sizes right
2930:     for (PetscInt point = pStart; point < pEnd; ++point) {
2931:       const PetscInt p = perm ? perm[point] : point;

2933:       switch (dm->blocking_type) {
2934:       case DM_BLOCKING_TOPOLOGICAL_POINT: { // One block per topological point
2935:         PetscInt bdof, offset;

2937:         PetscCall(PetscSectionGetDof(sectionGlobal, p, &dof));
2938:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &offset));
2939:         PetscCall(PetscSectionGetConstraintDof(sectionGlobal, p, &cdof));
2940:         if (blst && PetscBTLookup(blst, p)) bstate = PETSC_BOOL3_UNKNOWN;
2941:         if (dof > 0) {
2942:           // State change
2943:           if (bstate == PETSC_BOOL3_UNKNOWN) bstate = PETSC_BOOL3_TRUE;
2944:           else if (bstate == PETSC_BOOL3_TRUE && blst && !PetscBTLookup(blst, p)) bstate = PETSC_BOOL3_FALSE;

2946:           for (PetscInt i = 0; i < dof - cdof; ++i) pblocks[offset - localStart + i] = dof - cdof;
2947:           // Signal block concatenation
2948:           if (bstate == PETSC_BOOL3_FALSE && dof - cdof) pblocks[offset - localStart] = -(dof - cdof);
2949:         }
2950:         dof  = dof < 0 ? -(dof + 1) : dof;
2951:         bdof = cdof && (dof - cdof) ? 1 : dof;
2952:         if (dof) {
2953:           if (bs < 0) {
2954:             bs = bdof;
2955:           } else if (bs != bdof) {
2956:             bs = 1;
2957:           }
2958:         }
2959:       } break;
2960:       case DM_BLOCKING_FIELD_NODE: {
2961:         for (PetscInt field = 0; field < num_fields; field++) {
2962:           PetscInt num_comp, bdof, offset;
2963:           PetscCall(PetscSectionGetFieldComponents(sectionGlobal, field, &num_comp));
2964:           PetscCall(PetscSectionGetFieldDof(sectionGlobal, p, field, &dof));
2965:           if (dof < 0) continue;
2966:           PetscCall(PetscSectionGetFieldOffset(sectionGlobal, p, field, &offset));
2967:           PetscCall(PetscSectionGetFieldConstraintDof(sectionGlobal, p, field, &cdof));
2968:           PetscAssert(dof % num_comp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point %" PetscInt_FMT " field %" PetscInt_FMT " has %" PetscInt_FMT " dof, not divisible by %" PetscInt_FMT " component ", p, field, dof, num_comp);
2969:           PetscInt num_nodes = dof / num_comp;
2970:           for (PetscInt i = 0; i < dof - cdof; i++) pblocks[offset - localStart + i] = (dof - cdof) / num_nodes;
2971:           // Handle possibly constant block size (unlikely)
2972:           bdof = cdof && (dof - cdof) ? 1 : dof;
2973:           if (dof) {
2974:             if (bs < 0) {
2975:               bs = bdof;
2976:             } else if (bs != bdof) {
2977:               bs = 1;
2978:             }
2979:           }
2980:         }
2981:       } break;
2982:       }
2983:     }
2984:     if (sectionLocal->perm) PetscCall(ISRestoreIndices(sectionLocal->perm, &perm));
2985:     /* Must have same blocksize on all procs (some might have no points) */
2986:     bsLocal[0] = bs < 0 ? PETSC_INT_MAX : bs;
2987:     bsLocal[1] = bs;
2988:     PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax));
2989:     if (bsMinMax[0] != bsMinMax[1]) bs = 1;
2990:     else bs = bsMinMax[0];
2991:     bs = PetscMax(1, bs);
2992:     PetscCall(MatSetLocalToGlobalMapping(*J, ltog, ltog));
2993:     if (dm->prealloc_skip) { // User will likely use MatSetPreallocationCOO(), but still set structural parameters
2994:       PetscCall(MatSetBlockSize(*J, bs));
2995:       PetscCall(MatSetUp(*J));
2996:     } else {
2997:       PetscCall(PetscCalloc4(localSize / bs, &dnz, localSize / bs, &onz, localSize / bs, &dnzu, localSize / bs, &onzu));
2998:       PetscCall(DMPlexPreallocateOperator(dm, bs, dnz, onz, dnzu, onzu, *J, fillMatrix));
2999:       PetscCall(PetscFree4(dnz, onz, dnzu, onzu));
3000:     }
3001:     if (pblocks) { // Consolidate blocks
3002:       PetscInt nblocks = 0;
3003:       pblocks[0]       = PetscAbs(pblocks[0]);
3004:       for (PetscInt i = 0; i < localSize; i += PetscMax(1, pblocks[i])) {
3005:         if (pblocks[i] == 0) continue;
3006:         // Negative block size indicates the blocks should be concatenated
3007:         if (pblocks[i] < 0) {
3008:           pblocks[i] = -pblocks[i];
3009:           pblocks[nblocks - 1] += pblocks[i];
3010:         } else {
3011:           pblocks[nblocks++] = pblocks[i]; // nblocks always <= i
3012:         }
3013:         for (PetscInt j = 1; j < pblocks[i]; j++)
3014:           PetscCheck(pblocks[i + j] == pblocks[i], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Block of size %" PetscInt_FMT " at %" PetscInt_FMT " mismatches entry %" PetscInt_FMT " at %" PetscInt_FMT, pblocks[i], i, pblocks[i + j], i + j);
3015:       }
3016:       PetscCall(MatSetVariableBlockSizes(*J, nblocks, pblocks));
3017:     }
3018:     PetscCall(PetscFree(pblocks));
3019:   }
3020:   PetscCall(MatSetDM(*J, dm));
3021:   PetscFunctionReturn(PETSC_SUCCESS);
3022: }

3024: /*@
3025:   DMPlexGetSubdomainSection - Returns the section associated with the subdomain

3027:   Not Collective

3029:   Input Parameter:
3030: . dm - The `DMPLEX`

3032:   Output Parameter:
3033: . subsection - The subdomain section

3035:   Level: developer

3037: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSection`
3038: @*/
3039: PetscErrorCode DMPlexGetSubdomainSection(DM dm, PetscSection *subsection)
3040: {
3041:   DM_Plex *mesh = (DM_Plex *)dm->data;

3043:   PetscFunctionBegin;
3045:   if (!mesh->subdomainSection) {
3046:     PetscSection section;
3047:     PetscSF      sf;

3049:     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sf));
3050:     PetscCall(DMGetLocalSection(dm, &section));
3051:     PetscCall(PetscSectionCreateGlobalSection(section, sf, PETSC_TRUE, PETSC_FALSE, PETSC_TRUE, &mesh->subdomainSection));
3052:     PetscCall(PetscSFDestroy(&sf));
3053:   }
3054:   *subsection = mesh->subdomainSection;
3055:   PetscFunctionReturn(PETSC_SUCCESS);
3056: }

3058: /*@
3059:   DMPlexGetChart - Return the interval for all mesh points [`pStart`, `pEnd`)

3061:   Not Collective

3063:   Input Parameter:
3064: . dm - The `DMPLEX`

3066:   Output Parameters:
3067: + pStart - The first mesh point
3068: - pEnd   - The upper bound for mesh points

3070:   Level: beginner

3072: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetChart()`
3073: @*/
3074: PetscErrorCode DMPlexGetChart(DM dm, PetscInt *pStart, PetscInt *pEnd)
3075: {
3076:   DM_Plex *mesh = (DM_Plex *)dm->data;

3078:   PetscFunctionBegin;
3080:   if (mesh->tr) PetscCall(DMPlexTransformGetChart(mesh->tr, pStart, pEnd));
3081:   else PetscCall(PetscSectionGetChart(mesh->coneSection, pStart, pEnd));
3082:   PetscFunctionReturn(PETSC_SUCCESS);
3083: }

3085: /*@
3086:   DMPlexSetChart - Set the interval for all mesh points [`pStart`, `pEnd`)

3088:   Not Collective

3090:   Input Parameters:
3091: + dm     - The `DMPLEX`
3092: . pStart - The first mesh point
3093: - pEnd   - The upper bound for mesh points

3095:   Level: beginner

3097: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetChart()`
3098: @*/
3099: PetscErrorCode DMPlexSetChart(DM dm, PetscInt pStart, PetscInt pEnd)
3100: {
3101:   DM_Plex *mesh = (DM_Plex *)dm->data;

3103:   PetscFunctionBegin;
3105:   PetscCall(PetscSectionSetChart(mesh->coneSection, pStart, pEnd));
3106:   PetscCall(PetscSectionSetChart(mesh->supportSection, pStart, pEnd));
3107:   PetscCall(PetscFree(mesh->cellTypes));
3108:   PetscFunctionReturn(PETSC_SUCCESS);
3109: }

3111: /*@
3112:   DMPlexGetConeSize - Return the number of in-edges for this point in the DAG

3114:   Not Collective

3116:   Input Parameters:
3117: + dm - The `DMPLEX`
3118: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

3120:   Output Parameter:
3121: . size - The cone size for point `p`

3123:   Level: beginner

3125: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`
3126: @*/
3127: PetscErrorCode DMPlexGetConeSize(DM dm, PetscInt p, PetscInt *size)
3128: {
3129:   DM_Plex *mesh = (DM_Plex *)dm->data;

3131:   PetscFunctionBegin;
3133:   PetscAssertPointer(size, 3);
3134:   if (mesh->tr) PetscCall(DMPlexTransformGetConeSize(mesh->tr, p, size));
3135:   else PetscCall(PetscSectionGetDof(mesh->coneSection, p, size));
3136:   PetscFunctionReturn(PETSC_SUCCESS);
3137: }

3139: /*@
3140:   DMPlexSetConeSize - Set the number of in-edges for this point in the DAG

3142:   Not Collective

3144:   Input Parameters:
3145: + dm   - The `DMPLEX`
3146: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
3147: - size - The cone size for point `p`

3149:   Level: beginner

3151:   Note:
3152:   This should be called after `DMPlexSetChart()`.

3154: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetCone()`, `DMPlexCreate()`, `DMPlexGetConeSize()`, `DMPlexSetChart()`
3155: @*/
3156: PetscErrorCode DMPlexSetConeSize(DM dm, PetscInt p, PetscInt size)
3157: {
3158:   DM_Plex *mesh = (DM_Plex *)dm->data;

3160:   PetscFunctionBegin;
3162:   PetscCheck(!mesh->tr, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot call DMPlexSetConeSize() on a mesh with a transform defined.");
3163:   PetscCall(PetscSectionSetDof(mesh->coneSection, p, size));
3164:   PetscFunctionReturn(PETSC_SUCCESS);
3165: }

3167: /*@C
3168:   DMPlexGetCone - Return the points on the in-edges for this point in the DAG

3170:   Not Collective

3172:   Input Parameters:
3173: + dm - The `DMPLEX`
3174: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

3176:   Output Parameter:
3177: . cone - An array of points which are on the in-edges for point `p`, the length of `cone` is the result of `DMPlexGetConeSize()`

3179:   Level: beginner

3181:   Fortran Notes:
3182:   `cone` must be declared with
3183: .vb
3184:   PetscInt, pointer :: cone(:)
3185: .ve

3187:   You must also call `DMPlexRestoreCone()` after you finish using the array.
3188:   `DMPlexRestoreCone()` is not needed/available in C.

3190: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSize()`, `DMPlexSetCone()`, `DMPlexGetConeTuple()`, `DMPlexSetChart()`, `DMPlexRestoreCone()`
3191: @*/
3192: PetscErrorCode DMPlexGetCone(DM dm, PetscInt p, const PetscInt *cone[])
3193: {
3194:   DM_Plex *mesh = (DM_Plex *)dm->data;
3195:   PetscInt off;

3197:   PetscFunctionBegin;
3199:   PetscAssertPointer(cone, 3);
3200:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3201:   *cone = PetscSafePointerPlusOffset(mesh->cones, off);
3202:   PetscFunctionReturn(PETSC_SUCCESS);
3203: }

3205: /*@
3206:   DMPlexGetConeTuple - Return the points on the in-edges of several points in the DAG

3208:   Not Collective

3210:   Input Parameters:
3211: + dm - The `DMPLEX`
3212: - p  - The `IS` of points, which must lie in the chart set with `DMPlexSetChart()`

3214:   Output Parameters:
3215: + pConesSection - `PetscSection` describing the layout of `pCones`
3216: - pCones        - An `IS` containing the points which are on the in-edges for the point set `p`

3218:   Level: intermediate

3220: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeRecursive()`, `DMPlexSetChart()`, `PetscSection`, `IS`
3221: @*/
3222: PetscErrorCode DMPlexGetConeTuple(DM dm, IS p, PetscSection *pConesSection, IS *pCones)
3223: {
3224:   PetscSection cs, newcs;
3225:   PetscInt    *cones;
3226:   PetscInt    *newarr = NULL;
3227:   PetscInt     n;

3229:   PetscFunctionBegin;
3230:   PetscCall(DMPlexGetCones(dm, &cones));
3231:   PetscCall(DMPlexGetConeSection(dm, &cs));
3232:   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_INT, cones, p, &newcs, pCones ? ((void **)&newarr) : NULL));
3233:   if (pConesSection) *pConesSection = newcs;
3234:   if (pCones) {
3235:     PetscCall(PetscSectionGetStorageSize(newcs, &n));
3236:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)p), n, newarr, PETSC_OWN_POINTER, pCones));
3237:   }
3238:   PetscFunctionReturn(PETSC_SUCCESS);
3239: }

3241: /*@
3242:   DMPlexGetConeRecursiveVertices - Expand each given point into its cone points and do that recursively until we end up just with vertices.

3244:   Not Collective

3246:   Input Parameters:
3247: + dm     - The `DMPLEX`
3248: - points - The `IS` of points, which must lie in the chart set with `DMPlexSetChart()`

3250:   Output Parameter:
3251: . expandedPoints - An `IS` containing the of vertices recursively expanded from input points

3253:   Level: advanced

3255:   Notes:
3256:   Like `DMPlexGetConeRecursive()` but returns only the 0-depth `IS` (i.e. vertices only) and no sections.

3258:   There is no corresponding Restore function, just call `ISDestroy()` on the returned `IS` to deallocate.

3260: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeTuple()`, `DMPlexGetConeRecursive()`, `DMPlexRestoreConeRecursive()`,
3261:           `DMPlexGetDepth()`, `IS`
3262: @*/
3263: PetscErrorCode DMPlexGetConeRecursiveVertices(DM dm, IS points, IS *expandedPoints)
3264: {
3265:   IS      *expandedPointsAll;
3266:   PetscInt depth;

3268:   PetscFunctionBegin;
3271:   PetscAssertPointer(expandedPoints, 3);
3272:   PetscCall(DMPlexGetConeRecursive(dm, points, &depth, &expandedPointsAll, NULL));
3273:   *expandedPoints = expandedPointsAll[0];
3274:   PetscCall(PetscObjectReference((PetscObject)expandedPointsAll[0]));
3275:   PetscCall(DMPlexRestoreConeRecursive(dm, points, &depth, &expandedPointsAll, NULL));
3276:   PetscFunctionReturn(PETSC_SUCCESS);
3277: }

3279: /*@
3280:   DMPlexGetConeRecursive - Expand each given point into its cone points and do that recursively until we end up just with vertices
3281:   (DAG points of depth 0, i.e., without cones).

3283:   Not Collective

3285:   Input Parameters:
3286: + dm     - The `DMPLEX`
3287: - points - The `IS` of points, which must lie in the chart set with `DMPlexSetChart()`

3289:   Output Parameters:
3290: + depth          - (optional) Size of the output arrays, equal to `DMPLEX` depth, returned by `DMPlexGetDepth()`
3291: . expandedPoints - (optional) An array of index sets with recursively expanded cones
3292: - sections       - (optional) An array of sections which describe mappings from points to their cone points

3294:   Level: advanced

3296:   Notes:
3297:   Like `DMPlexGetConeTuple()` but recursive.

3299:   Array `expandedPoints` has size equal to `depth`. Each `expandedPoints`[d] contains DAG points with maximum depth d, recursively cone-wise expanded from the input points.
3300:   For example, for d=0 it contains only vertices, for d=1 it can contain vertices and edges, etc.

3302:   Array section has size equal to `depth`.  Each `PetscSection` `sections`[d] realizes mapping from `expandedPoints`[d+1] (section points) to `expandedPoints`[d] (section dofs) as follows\:
3303:   (1) DAG points in `expandedPoints`[d+1] with `depth` d+1 to their cone points in `expandedPoints`[d];
3304:   (2) DAG points in `expandedPoints`[d+1] with `depth` in [0,d] to the same points in `expandedPoints`[d].

3306: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeTuple()`, `DMPlexRestoreConeRecursive()`, `DMPlexGetConeRecursiveVertices()`,
3307:           `DMPlexGetDepth()`, `PetscSection`, `IS`
3308: @*/
3309: PetscErrorCode DMPlexGetConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
3310: {
3311:   const PetscInt *arr0 = NULL, *cone = NULL;
3312:   PetscInt       *arr = NULL, *newarr = NULL;
3313:   PetscInt        d, depth_, i, n, newn, cn, co, start, end;
3314:   IS             *expandedPoints_;
3315:   PetscSection   *sections_;

3317:   PetscFunctionBegin;
3320:   if (depth) PetscAssertPointer(depth, 3);
3321:   if (expandedPoints) PetscAssertPointer(expandedPoints, 4);
3322:   if (sections) PetscAssertPointer(sections, 5);
3323:   PetscCall(ISGetLocalSize(points, &n));
3324:   PetscCall(ISGetIndices(points, &arr0));
3325:   PetscCall(DMPlexGetDepth(dm, &depth_));
3326:   PetscCall(PetscCalloc1(depth_, &expandedPoints_));
3327:   PetscCall(PetscCalloc1(depth_, &sections_));
3328:   arr = (PetscInt *)arr0; /* this is ok because first generation of arr is not modified */
3329:   for (d = depth_ - 1; d >= 0; d--) {
3330:     PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &sections_[d]));
3331:     PetscCall(PetscSectionSetChart(sections_[d], 0, n));
3332:     for (i = 0; i < n; i++) {
3333:       PetscCall(DMPlexGetDepthStratum(dm, d + 1, &start, &end));
3334:       if (arr[i] >= start && arr[i] < end) {
3335:         PetscCall(DMPlexGetConeSize(dm, arr[i], &cn));
3336:         PetscCall(PetscSectionSetDof(sections_[d], i, cn));
3337:       } else {
3338:         PetscCall(PetscSectionSetDof(sections_[d], i, 1));
3339:       }
3340:     }
3341:     PetscCall(PetscSectionSetUp(sections_[d]));
3342:     PetscCall(PetscSectionGetStorageSize(sections_[d], &newn));
3343:     PetscCall(PetscMalloc1(newn, &newarr));
3344:     for (i = 0; i < n; i++) {
3345:       PetscCall(PetscSectionGetDof(sections_[d], i, &cn));
3346:       PetscCall(PetscSectionGetOffset(sections_[d], i, &co));
3347:       if (cn > 1) {
3348:         PetscCall(DMPlexGetCone(dm, arr[i], &cone));
3349:         PetscCall(PetscMemcpy(&newarr[co], cone, cn * sizeof(PetscInt)));
3350:       } else {
3351:         newarr[co] = arr[i];
3352:       }
3353:     }
3354:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newarr, PETSC_OWN_POINTER, &expandedPoints_[d]));
3355:     arr = newarr;
3356:     n   = newn;
3357:   }
3358:   PetscCall(ISRestoreIndices(points, &arr0));
3359:   *depth = depth_;
3360:   if (expandedPoints) *expandedPoints = expandedPoints_;
3361:   else {
3362:     for (d = 0; d < depth_; d++) PetscCall(ISDestroy(&expandedPoints_[d]));
3363:     PetscCall(PetscFree(expandedPoints_));
3364:   }
3365:   if (sections) *sections = sections_;
3366:   else {
3367:     for (d = 0; d < depth_; d++) PetscCall(PetscSectionDestroy(&sections_[d]));
3368:     PetscCall(PetscFree(sections_));
3369:   }
3370:   PetscFunctionReturn(PETSC_SUCCESS);
3371: }

3373: /*@
3374:   DMPlexRestoreConeRecursive - Deallocates arrays created by `DMPlexGetConeRecursive()`

3376:   Not Collective

3378:   Input Parameters:
3379: + dm     - The `DMPLEX`
3380: - points - The `IS` of points, which must lie in the chart set with `DMPlexSetChart()`

3382:   Output Parameters:
3383: + depth          - (optional) Size of the output arrays, equal to `DMPLEX` depth, returned by `DMPlexGetDepth()`
3384: . expandedPoints - (optional) An array of recursively expanded cones
3385: - sections       - (optional) An array of sections which describe mappings from points to their cone points

3387:   Level: advanced

3389:   Note:
3390:   See `DMPlexGetConeRecursive()`

3392: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeTuple()`, `DMPlexGetConeRecursive()`, `DMPlexGetConeRecursiveVertices()`,
3393:           `DMPlexGetDepth()`, `IS`, `PetscSection`
3394: @*/
3395: PetscErrorCode DMPlexRestoreConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
3396: {
3397:   PetscInt d, depth_;

3399:   PetscFunctionBegin;
3400:   PetscCall(DMPlexGetDepth(dm, &depth_));
3401:   PetscCheck(!depth || *depth == depth_, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "depth changed since last call to DMPlexGetConeRecursive");
3402:   if (depth) *depth = 0;
3403:   if (expandedPoints) {
3404:     for (d = 0; d < depth_; d++) PetscCall(ISDestroy(&(*expandedPoints)[d]));
3405:     PetscCall(PetscFree(*expandedPoints));
3406:   }
3407:   if (sections) {
3408:     for (d = 0; d < depth_; d++) PetscCall(PetscSectionDestroy(&(*sections)[d]));
3409:     PetscCall(PetscFree(*sections));
3410:   }
3411:   PetscFunctionReturn(PETSC_SUCCESS);
3412: }

3414: /*@
3415:   DMPlexSetCone - Set the points on the in-edges for this point in the DAG; that is these are the points that cover the specific point

3417:   Not Collective

3419:   Input Parameters:
3420: + dm   - The `DMPLEX`
3421: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
3422: - cone - An array of points which are on the in-edges for point `p`, its length must have been previously provided with `DMPlexSetConeSize()`

3424:   Level: beginner

3426:   Note:
3427:   This should be called after all calls to `DMPlexSetConeSize()` and `DMSetUp()`.

3429: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`, `DMPlexSetSupport()`, `DMPlexSetSupportSize()`
3430: @*/
3431: PetscErrorCode DMPlexSetCone(DM dm, PetscInt p, const PetscInt cone[])
3432: {
3433:   DM_Plex *mesh = (DM_Plex *)dm->data;
3434:   PetscInt dof, off, c;

3436:   PetscFunctionBegin;
3438:   PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3439:   if (dof) PetscAssertPointer(cone, 3);
3440:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3441:   if (PetscDefined(USE_DEBUG)) {
3442:     PetscInt pStart, pEnd;
3443:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3444:     PetscCheck(!(p < pStart) && !(p >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
3445:     for (c = 0; c < dof; ++c) {
3446:       PetscCheck(!(cone[c] < pStart) && !(cone[c] >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", cone[c], pStart, pEnd);
3447:       mesh->cones[off + c] = cone[c];
3448:     }
3449:   } else {
3450:     for (c = 0; c < dof; ++c) mesh->cones[off + c] = cone[c];
3451:   }
3452:   PetscFunctionReturn(PETSC_SUCCESS);
3453: }

3455: /*@C
3456:   DMPlexGetConeOrientation - Return the orientations on the in-edges for this point in the DAG

3458:   Not Collective

3460:   Input Parameters:
3461: + dm - The `DMPLEX`
3462: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

3464:   Output Parameter:
3465: . coneOrientation - An array of orientations which are on the in-edges for point `p`. An orientation is an
3466:                     integer giving the prescription for cone traversal. Its length is given by the result of `DMPlexSetConeSize()`

3468:   Level: beginner

3470:   Note:
3471:   The number indexes the symmetry transformations for the cell type (see manual). Orientation 0 is always
3472:   the identity transformation. Negative orientation indicates reflection so that -(o+1) is the reflection
3473:   of o, however it is not necessarily the inverse. To get the inverse, use `DMPolytopeTypeComposeOrientationInv()`
3474:   with the identity.

3476:   Fortran Notes:
3477:   You must call `DMPlexRestoreConeOrientation()` after you finish using the returned array.
3478:   `DMPlexRestoreConeOrientation()` is not needed/available in C.

3480: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetConeSize()`, `DMPolytopeTypeComposeOrientation()`, `DMPolytopeTypeComposeOrientationInv()`,
3481:           `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetCone()`, `DMPlexSetChart()`
3482: @*/
3483: PetscErrorCode DMPlexGetConeOrientation(DM dm, PetscInt p, const PetscInt *coneOrientation[])
3484: {
3485:   DM_Plex *mesh = (DM_Plex *)dm->data;
3486:   PetscInt off;

3488:   PetscFunctionBegin;
3490:   if (PetscDefined(USE_DEBUG)) {
3491:     PetscInt dof;
3492:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3493:     if (dof) PetscAssertPointer(coneOrientation, 3);
3494:   }
3495:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));

3497:   *coneOrientation = &mesh->coneOrientations[off];
3498:   PetscFunctionReturn(PETSC_SUCCESS);
3499: }

3501: /*@
3502:   DMPlexSetConeOrientation - Set the orientations on the in-edges for this point in the DAG

3504:   Not Collective

3506:   Input Parameters:
3507: + dm              - The `DMPLEX`
3508: . p               - The point, which must lie in the chart set with `DMPlexSetChart()`
3509: - coneOrientation - An array of orientations. Its length is given by the result of `DMPlexSetConeSize()`

3511:   Level: beginner

3513:   Notes:
3514:   This should be called after all calls to `DMPlexSetConeSize()` and `DMSetUp()`.

3516:   The meaning of coneOrientation is detailed in `DMPlexGetConeOrientation()`.

3518: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetConeOrientation()`, `DMPlexSetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3519: @*/
3520: PetscErrorCode DMPlexSetConeOrientation(DM dm, PetscInt p, const PetscInt coneOrientation[])
3521: {
3522:   DM_Plex *mesh = (DM_Plex *)dm->data;
3523:   PetscInt pStart, pEnd;
3524:   PetscInt dof, off, c;

3526:   PetscFunctionBegin;
3528:   PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3529:   if (dof) PetscAssertPointer(coneOrientation, 3);
3530:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3531:   if (PetscDefined(USE_DEBUG)) {
3532:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3533:     PetscCheck(!(p < pStart) && !(p >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
3534:     for (c = 0; c < dof; ++c) {
3535:       PetscInt cdof, o = coneOrientation[c];

3537:       PetscCall(PetscSectionGetDof(mesh->coneSection, mesh->cones[off + c], &cdof));
3538:       PetscCheck(!o || (o >= -(cdof + 1) && o < cdof), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone orientation %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ". %" PetscInt_FMT ")", o, -(cdof + 1), cdof);
3539:       mesh->coneOrientations[off + c] = o;
3540:     }
3541:   } else {
3542:     for (c = 0; c < dof; ++c) mesh->coneOrientations[off + c] = coneOrientation[c];
3543:   }
3544:   PetscFunctionReturn(PETSC_SUCCESS);
3545: }

3547: /*@
3548:   DMPlexInsertCone - Insert a point into the in-edges for the point p in the DAG

3550:   Not Collective

3552:   Input Parameters:
3553: + dm        - The `DMPLEX`
3554: . p         - The point, which must lie in the chart set with `DMPlexSetChart()`
3555: . conePos   - The local index in the cone where the point should be put
3556: - conePoint - The mesh point to insert

3558:   Level: beginner

3560: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3561: @*/
3562: PetscErrorCode DMPlexInsertCone(DM dm, PetscInt p, PetscInt conePos, PetscInt conePoint)
3563: {
3564:   DM_Plex *mesh = (DM_Plex *)dm->data;
3565:   PetscInt pStart, pEnd;
3566:   PetscInt dof, off;

3568:   PetscFunctionBegin;
3570:   if (PetscDefined(USE_DEBUG)) {
3571:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3572:     PetscCheck(!(p < pStart) && !(p >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
3573:     PetscCheck(!(conePoint < pStart) && !(conePoint >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", conePoint, pStart, pEnd);
3574:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3575:     PetscCheck(!(conePos < 0) && !(conePos >= dof), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone position %" PetscInt_FMT " of point %" PetscInt_FMT " is not in the valid range [0, %" PetscInt_FMT ")", conePos, p, dof);
3576:   }
3577:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3578:   mesh->cones[off + conePos] = conePoint;
3579:   PetscFunctionReturn(PETSC_SUCCESS);
3580: }

3582: /*@
3583:   DMPlexInsertConeOrientation - Insert a point orientation for the in-edge for the point p in the DAG

3585:   Not Collective

3587:   Input Parameters:
3588: + dm              - The `DMPLEX`
3589: . p               - The point, which must lie in the chart set with `DMPlexSetChart()`
3590: . conePos         - The local index in the cone where the point should be put
3591: - coneOrientation - The point orientation to insert

3593:   Level: beginner

3595:   Note:
3596:   The meaning of coneOrientation values is detailed in `DMPlexGetConeOrientation()`.

3598: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3599: @*/
3600: PetscErrorCode DMPlexInsertConeOrientation(DM dm, PetscInt p, PetscInt conePos, PetscInt coneOrientation)
3601: {
3602:   DM_Plex *mesh = (DM_Plex *)dm->data;
3603:   PetscInt pStart, pEnd;
3604:   PetscInt dof, off;

3606:   PetscFunctionBegin;
3608:   if (PetscDefined(USE_DEBUG)) {
3609:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3610:     PetscCheck(!(p < pStart) && !(p >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
3611:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3612:     PetscCheck(!(conePos < 0) && !(conePos >= dof), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone position %" PetscInt_FMT " of point %" PetscInt_FMT " is not in the valid range [0, %" PetscInt_FMT ")", conePos, p, dof);
3613:   }
3614:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3615:   mesh->coneOrientations[off + conePos] = coneOrientation;
3616:   PetscFunctionReturn(PETSC_SUCCESS);
3617: }

3619: /*@C
3620:   DMPlexGetOrientedCone - Return the points and orientations on the in-edges for this point in the DAG

3622:   Not collective

3624:   Input Parameters:
3625: + dm - The DMPlex
3626: - p  - The point, which must lie in the chart set with DMPlexSetChart()

3628:   Output Parameters:
3629: + cone - An array of points which are on the in-edges for point `p`
3630: - ornt - An array of orientations which are on the in-edges for point `p`. An orientation is an
3631:          integer giving the prescription for cone traversal.

3633:   Level: beginner

3635:   Notes:
3636:   The number indexes the symmetry transformations for the cell type (see manual). Orientation 0 is always
3637:   the identity transformation. Negative orientation indicates reflection so that -(o+1) is the reflection
3638:   of o, however it is not necessarily the inverse. To get the inverse, use `DMPolytopeTypeComposeOrientationInv()`
3639:   with the identity.

3641:   You must also call `DMPlexRestoreOrientedCone()` after you finish using the returned array.

3643:   Fortran Notes:
3644:   `cone` and `ornt` must be declared with
3645: .vb
3646:   PetscInt, pointer :: cone(:)
3647:   PetscInt, pointer :: ornt(:)
3648: .ve

3650: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreOrientedCone()`, `DMPlexGetConeSize()`, `DMPlexGetCone()`, `DMPlexGetChart()`
3651: @*/
3652: PetscErrorCode DMPlexGetOrientedCone(DM dm, PetscInt p, const PetscInt *cone[], const PetscInt *ornt[])
3653: {
3654:   DM_Plex *mesh = (DM_Plex *)dm->data;

3656:   PetscFunctionBegin;
3658:   if (mesh->tr) {
3659:     PetscCall(DMPlexTransformGetCone(mesh->tr, p, cone, ornt));
3660:   } else {
3661:     PetscInt off;
3662:     if (PetscDefined(USE_DEBUG)) {
3663:       PetscInt dof;
3664:       PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3665:       if (dof) {
3666:         if (cone) PetscAssertPointer(cone, 3);
3667:         if (ornt) PetscAssertPointer(ornt, 4);
3668:       }
3669:     }
3670:     PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3671:     if (cone) *cone = PetscSafePointerPlusOffset(mesh->cones, off);
3672:     if (ornt) *ornt = PetscSafePointerPlusOffset(mesh->coneOrientations, off);
3673:   }
3674:   PetscFunctionReturn(PETSC_SUCCESS);
3675: }

3677: /*@C
3678:   DMPlexRestoreOrientedCone - Restore the points and orientations on the in-edges for this point in the DAG obtained with `DMPlexGetOrientedCone()`

3680:   Not Collective

3682:   Input Parameters:
3683: + dm   - The DMPlex
3684: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
3685: . cone - An array of points which are on the in-edges for point p
3686: - ornt - An array of orientations which are on the in-edges for point `p`. An orientation is an
3687:          integer giving the prescription for cone traversal.

3689:   Level: beginner

3691: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetOrientedCone()`, `DMPlexGetConeSize()`, `DMPlexGetCone()`, `DMPlexGetChart()`
3692: @*/
3693: PetscErrorCode DMPlexRestoreOrientedCone(DM dm, PetscInt p, const PetscInt *cone[], const PetscInt *ornt[])
3694: {
3695:   DM_Plex *mesh = (DM_Plex *)dm->data;

3697:   PetscFunctionBegin;
3699:   if (mesh->tr) PetscCall(DMPlexTransformRestoreCone(mesh->tr, p, cone, ornt));
3700:   PetscFunctionReturn(PETSC_SUCCESS);
3701: }

3703: /*@
3704:   DMPlexGetSupportSize - Return the number of out-edges for this point in the DAG

3706:   Not Collective

3708:   Input Parameters:
3709: + dm - The `DMPLEX`
3710: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

3712:   Output Parameter:
3713: . size - The support size for point `p`

3715:   Level: beginner

3717: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`, `DMPlexGetConeSize()`
3718: @*/
3719: PetscErrorCode DMPlexGetSupportSize(DM dm, PetscInt p, PetscInt *size)
3720: {
3721:   DM_Plex *mesh = (DM_Plex *)dm->data;

3723:   PetscFunctionBegin;
3725:   PetscAssertPointer(size, 3);
3726:   PetscCall(PetscSectionGetDof(mesh->supportSection, p, size));
3727:   PetscFunctionReturn(PETSC_SUCCESS);
3728: }

3730: /*@
3731:   DMPlexSetSupportSize - Set the number of out-edges for this point in the DAG

3733:   Not Collective

3735:   Input Parameters:
3736: + dm   - The `DMPLEX`
3737: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
3738: - size - The support size for point `p`

3740:   Level: beginner

3742:   Note:
3743:   This should be called after `DMPlexSetChart()`.

3745: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetSupportSize()`, `DMPlexSetChart()`
3746: @*/
3747: PetscErrorCode DMPlexSetSupportSize(DM dm, PetscInt p, PetscInt size)
3748: {
3749:   DM_Plex *mesh = (DM_Plex *)dm->data;

3751:   PetscFunctionBegin;
3753:   PetscCall(PetscSectionSetDof(mesh->supportSection, p, size));
3754:   PetscFunctionReturn(PETSC_SUCCESS);
3755: }

3757: /*@C
3758:   DMPlexGetSupport - Return the points on the out-edges for this point in the DAG

3760:   Not Collective

3762:   Input Parameters:
3763: + dm - The `DMPLEX`
3764: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

3766:   Output Parameter:
3767: . support - An array of points which are on the out-edges for point `p`, its length is that obtained from `DMPlexGetSupportSize()`

3769:   Level: beginner

3771:   Fortran Notes:
3772:   `support` must be declared with
3773: .vb
3774:   PetscInt, pointer :: support(:)
3775: .ve

3777:   You must also call `DMPlexRestoreSupport()` after you finish using the returned array.
3778:   `DMPlexRestoreSupport()` is not needed/available in C.

3780: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetSupportSize()`, `DMPlexSetSupport()`, `DMPlexGetCone()`, `DMPlexSetChart()`
3781: @*/
3782: PetscErrorCode DMPlexGetSupport(DM dm, PetscInt p, const PetscInt *support[])
3783: {
3784:   DM_Plex *mesh = (DM_Plex *)dm->data;
3785:   PetscInt off;

3787:   PetscFunctionBegin;
3789:   PetscAssertPointer(support, 3);
3790:   PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
3791:   *support = PetscSafePointerPlusOffset(mesh->supports, off);
3792:   PetscFunctionReturn(PETSC_SUCCESS);
3793: }

3795: /*@
3796:   DMPlexSetSupport - Set the points on the out-edges for this point in the DAG, that is the list of points that this point covers

3798:   Not Collective

3800:   Input Parameters:
3801: + dm      - The `DMPLEX`
3802: . p       - The point, which must lie in the chart set with `DMPlexSetChart()`
3803: - support - An array of points which are on the out-edges for point `p`, its length is that obtained from `DMPlexGetSupportSize()`

3805:   Level: beginner

3807:   Note:
3808:   This should be called after all calls to `DMPlexSetSupportSize()` and `DMSetUp()`.

3810: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetCone()`, `DMPlexSetConeSize()`, `DMPlexCreate()`, `DMPlexGetSupport()`, `DMPlexSetChart()`, `DMPlexSetSupportSize()`, `DMSetUp()`
3811: @*/
3812: PetscErrorCode DMPlexSetSupport(DM dm, PetscInt p, const PetscInt support[])
3813: {
3814:   DM_Plex *mesh = (DM_Plex *)dm->data;
3815:   PetscInt pStart, pEnd;
3816:   PetscInt dof, off, c;

3818:   PetscFunctionBegin;
3820:   PetscCall(PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd));
3821:   PetscCall(PetscSectionGetDof(mesh->supportSection, p, &dof));
3822:   if (dof) PetscAssertPointer(support, 3);
3823:   PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
3824:   PetscCheck(!(p < pStart) && !(p >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
3825:   for (c = 0; c < dof; ++c) {
3826:     PetscCheck(!(support[c] < pStart) && !(support[c] >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Support point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", support[c], pStart, pEnd);
3827:     mesh->supports[off + c] = support[c];
3828:   }
3829:   PetscFunctionReturn(PETSC_SUCCESS);
3830: }

3832: /*@
3833:   DMPlexInsertSupport - Insert a point into the out-edges for the point p in the DAG

3835:   Not Collective

3837:   Input Parameters:
3838: + dm           - The `DMPLEX`
3839: . p            - The point, which must lie in the chart set with `DMPlexSetChart()`
3840: . supportPos   - The local index in the cone where the point should be put
3841: - supportPoint - The mesh point to insert

3843:   Level: beginner

3845: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3846: @*/
3847: PetscErrorCode DMPlexInsertSupport(DM dm, PetscInt p, PetscInt supportPos, PetscInt supportPoint)
3848: {
3849:   DM_Plex *mesh = (DM_Plex *)dm->data;
3850:   PetscInt pStart, pEnd;
3851:   PetscInt dof, off;

3853:   PetscFunctionBegin;
3855:   PetscCall(PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd));
3856:   PetscCall(PetscSectionGetDof(mesh->supportSection, p, &dof));
3857:   PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
3858:   PetscCheck(!(p < pStart) && !(p >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
3859:   PetscCheck(!(supportPoint < pStart) && !(supportPoint >= pEnd), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Support point %" PetscInt_FMT " is not in the valid range [%" PetscInt_FMT ", %" PetscInt_FMT ")", supportPoint, pStart, pEnd);
3860:   PetscCheck(supportPos < dof, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Support position %" PetscInt_FMT " of point %" PetscInt_FMT " is not in the valid range [0, %" PetscInt_FMT ")", supportPos, p, dof);
3861:   mesh->supports[off + supportPos] = supportPoint;
3862:   PetscFunctionReturn(PETSC_SUCCESS);
3863: }

3865: /* Converts an orientation o in the current numbering to the previous scheme used in Plex */
3866: PetscInt DMPolytopeConvertNewOrientation_Internal(DMPolytopeType ct, PetscInt o)
3867: {
3868:   switch (ct) {
3869:   case DM_POLYTOPE_SEGMENT:
3870:     if (o == -1) return -2;
3871:     break;
3872:   case DM_POLYTOPE_TRIANGLE:
3873:     if (o == -3) return -1;
3874:     if (o == -2) return -3;
3875:     if (o == -1) return -2;
3876:     break;
3877:   case DM_POLYTOPE_QUADRILATERAL:
3878:     if (o == -4) return -2;
3879:     if (o == -3) return -1;
3880:     if (o == -2) return -4;
3881:     if (o == -1) return -3;
3882:     break;
3883:   default:
3884:     return o;
3885:   }
3886:   return o;
3887: }

3889: /* Converts an orientation o in the previous scheme used in Plex to the current numbering */
3890: PetscInt DMPolytopeConvertOldOrientation_Internal(DMPolytopeType ct, PetscInt o)
3891: {
3892:   switch (ct) {
3893:   case DM_POLYTOPE_SEGMENT:
3894:     if ((o == -2) || (o == 1)) return -1;
3895:     if (o == -1) return 0;
3896:     break;
3897:   case DM_POLYTOPE_TRIANGLE:
3898:     if (o == -3) return -2;
3899:     if (o == -2) return -1;
3900:     if (o == -1) return -3;
3901:     break;
3902:   case DM_POLYTOPE_QUADRILATERAL:
3903:     if (o == -4) return -2;
3904:     if (o == -3) return -1;
3905:     if (o == -2) return -4;
3906:     if (o == -1) return -3;
3907:     break;
3908:   default:
3909:     return o;
3910:   }
3911:   return o;
3912: }

3914: /* Takes in a mesh whose orientations are in the previous scheme and converts them all to the current numbering */
3915: PetscErrorCode DMPlexConvertOldOrientations_Internal(DM dm)
3916: {
3917:   PetscInt pStart, pEnd, p;

3919:   PetscFunctionBegin;
3920:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
3921:   for (p = pStart; p < pEnd; ++p) {
3922:     const PetscInt *cone, *ornt;
3923:     PetscInt        coneSize, c;

3925:     PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
3926:     PetscCall(DMPlexGetCone(dm, p, &cone));
3927:     PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
3928:     for (c = 0; c < coneSize; ++c) {
3929:       DMPolytopeType ct;
3930:       const PetscInt o = ornt[c];

3932:       PetscCall(DMPlexGetCellType(dm, cone[c], &ct));
3933:       switch (ct) {
3934:       case DM_POLYTOPE_SEGMENT:
3935:         if ((o == -2) || (o == 1)) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -1));
3936:         if (o == -1) PetscCall(DMPlexInsertConeOrientation(dm, p, c, 0));
3937:         break;
3938:       case DM_POLYTOPE_TRIANGLE:
3939:         if (o == -3) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -2));
3940:         if (o == -2) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -1));
3941:         if (o == -1) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -3));
3942:         break;
3943:       case DM_POLYTOPE_QUADRILATERAL:
3944:         if (o == -4) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -2));
3945:         if (o == -3) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -1));
3946:         if (o == -2) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -4));
3947:         if (o == -1) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -3));
3948:         break;
3949:       default:
3950:         break;
3951:       }
3952:     }
3953:   }
3954:   PetscFunctionReturn(PETSC_SUCCESS);
3955: }

3957: static inline PetscErrorCode DMPlexGetTransitiveClosure_Hot_Private(DM dm, PetscInt p, PetscBool useCone, PetscInt *size, const PetscInt *arr[], const PetscInt *ornt[])
3958: {
3959:   DM_Plex *mesh = (DM_Plex *)dm->data;

3961:   PetscFunctionBeginHot;
3962:   if (PetscDefined(USE_DEBUG) || mesh->tr) {
3963:     if (useCone) {
3964:       PetscCall(DMPlexGetConeSize(dm, p, size));
3965:       PetscCall(DMPlexGetOrientedCone(dm, p, arr, ornt));
3966:     } else {
3967:       PetscCall(DMPlexGetSupportSize(dm, p, size));
3968:       PetscCall(DMPlexGetSupport(dm, p, arr));
3969:     }
3970:   } else {
3971:     if (useCone) {
3972:       const PetscSection s   = mesh->coneSection;
3973:       const PetscInt     ps  = p - s->pStart;
3974:       const PetscInt     off = s->atlasOff[ps];

3976:       *size = s->atlasDof[ps];
3977:       *arr  = mesh->cones + off;
3978:       *ornt = mesh->coneOrientations + off;
3979:     } else {
3980:       const PetscSection s   = mesh->supportSection;
3981:       const PetscInt     ps  = p - s->pStart;
3982:       const PetscInt     off = s->atlasOff[ps];

3984:       *size = s->atlasDof[ps];
3985:       *arr  = mesh->supports + off;
3986:     }
3987:   }
3988:   PetscFunctionReturn(PETSC_SUCCESS);
3989: }

3991: static inline PetscErrorCode DMPlexRestoreTransitiveClosure_Hot_Private(DM dm, PetscInt p, PetscBool useCone, PetscInt *size, const PetscInt *arr[], const PetscInt *ornt[])
3992: {
3993:   DM_Plex *mesh = (DM_Plex *)dm->data;

3995:   PetscFunctionBeginHot;
3996:   if (PetscDefined(USE_DEBUG) || mesh->tr) {
3997:     if (useCone) PetscCall(DMPlexRestoreOrientedCone(dm, p, arr, ornt));
3998:   }
3999:   PetscFunctionReturn(PETSC_SUCCESS);
4000: }

4002: static PetscErrorCode DMPlexGetTransitiveClosure_Depth1_Private(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
4003: {
4004:   DMPolytopeType  ct = DM_POLYTOPE_UNKNOWN;
4005:   PetscInt       *closure;
4006:   const PetscInt *tmp = NULL, *tmpO = NULL;
4007:   PetscInt        off = 0, tmpSize, t;

4009:   PetscFunctionBeginHot;
4010:   if (ornt) {
4011:     PetscCall(DMPlexGetCellType(dm, p, &ct));
4012:     if (ct == DM_POLYTOPE_FV_GHOST || ct == DM_POLYTOPE_INTERIOR_GHOST || ct == DM_POLYTOPE_UNKNOWN || ct == DM_POLYTOPE_UNKNOWN_CELL || ct == DM_POLYTOPE_UNKNOWN_FACE) ct = DM_POLYTOPE_UNKNOWN;
4013:   }
4014:   if (*points) {
4015:     closure = *points;
4016:   } else {
4017:     PetscInt maxConeSize, maxSupportSize;
4018:     PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
4019:     PetscCall(DMGetWorkArray(dm, 2 * (PetscMax(maxConeSize, maxSupportSize) + 1), MPIU_INT, &closure));
4020:   }
4021:   PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, p, useCone, &tmpSize, &tmp, &tmpO));
4022:   if (ct == DM_POLYTOPE_UNKNOWN) {
4023:     closure[off++] = p;
4024:     closure[off++] = 0;
4025:     for (t = 0; t < tmpSize; ++t) {
4026:       closure[off++] = tmp[t];
4027:       closure[off++] = tmpO ? tmpO[t] : 0;
4028:     }
4029:   } else {
4030:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, ornt);

4032:     /* We assume that cells with a valid type have faces with a valid type */
4033:     closure[off++] = p;
4034:     closure[off++] = ornt;
4035:     for (t = 0; t < tmpSize; ++t) {
4036:       DMPolytopeType ft;

4038:       PetscCall(DMPlexGetCellType(dm, tmp[t], &ft));
4039:       closure[off++] = tmp[arr[t]];
4040:       closure[off++] = tmpO ? DMPolytopeTypeComposeOrientation(ft, ornt, tmpO[t]) : 0;
4041:     }
4042:   }
4043:   PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, p, useCone, &tmpSize, &tmp, &tmpO));
4044:   if (numPoints) *numPoints = tmpSize + 1;
4045:   if (points) *points = closure;
4046:   PetscFunctionReturn(PETSC_SUCCESS);
4047: }

4049: /* We need a special tensor version because we want to allow duplicate points in the endcaps for hybrid cells */
4050: static PetscErrorCode DMPlexTransitiveClosure_Tensor_Internal(DM dm, PetscInt point, DMPolytopeType ct, PetscInt o, PetscBool useCone, PetscInt *numPoints, PetscInt **points)
4051: {
4052:   const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);
4053:   const PetscInt *cone, *ornt;
4054:   PetscInt       *pts, *closure = NULL;
4055:   DMPolytopeType  ft;
4056:   PetscInt        maxConeSize, maxSupportSize, coneSeries, supportSeries, maxSize;
4057:   PetscInt        dim, coneSize, c, d, clSize, cl;

4059:   PetscFunctionBeginHot;
4060:   PetscCall(DMGetDimension(dm, &dim));
4061:   PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, point, PETSC_TRUE, &coneSize, &cone, &ornt));
4062:   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
4063:   coneSeries    = (maxConeSize > 1) ? ((PetscPowInt(maxConeSize, dim + 1) - 1) / (maxConeSize - 1)) : dim + 1;
4064:   supportSeries = (maxSupportSize > 1) ? ((PetscPowInt(maxSupportSize, dim + 1) - 1) / (maxSupportSize - 1)) : dim + 1;
4065:   maxSize       = PetscMax(coneSeries, supportSeries);
4066:   if (*points) {
4067:     pts = *points;
4068:   } else PetscCall(DMGetWorkArray(dm, 2 * maxSize, MPIU_INT, &pts));
4069:   c        = 0;
4070:   pts[c++] = point;
4071:   pts[c++] = o;
4072:   PetscCall(DMPlexGetCellType(dm, cone[arr[0 * 2 + 0]], &ft));
4073:   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, cone[arr[0 * 2 + 0]], DMPolytopeTypeComposeOrientation(ft, arr[0 * 2 + 1], ornt[0]), useCone, &clSize, &closure));
4074:   for (cl = 0; cl < clSize * 2; cl += 2) {
4075:     pts[c++] = closure[cl];
4076:     pts[c++] = closure[cl + 1];
4077:   }
4078:   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, cone[arr[1 * 2 + 0]], DMPolytopeTypeComposeOrientation(ft, arr[1 * 2 + 1], ornt[1]), useCone, &clSize, &closure));
4079:   for (cl = 0; cl < clSize * 2; cl += 2) {
4080:     pts[c++] = closure[cl];
4081:     pts[c++] = closure[cl + 1];
4082:   }
4083:   PetscCall(DMPlexRestoreTransitiveClosure(dm, cone[0], useCone, &clSize, &closure));
4084:   for (d = 2; d < coneSize; ++d) {
4085:     PetscCall(DMPlexGetCellType(dm, cone[arr[d * 2 + 0]], &ft));
4086:     pts[c++] = cone[arr[d * 2 + 0]];
4087:     pts[c++] = DMPolytopeTypeComposeOrientation(ft, arr[d * 2 + 1], ornt[d]);
4088:   }
4089:   PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, point, PETSC_TRUE, &coneSize, &cone, &ornt));
4090:   if (dim >= 3) {
4091:     for (d = 2; d < coneSize; ++d) {
4092:       const PetscInt  fpoint = cone[arr[d * 2 + 0]];
4093:       const PetscInt *fcone, *fornt;
4094:       PetscInt        fconeSize, fc, i;

4096:       PetscCall(DMPlexGetCellType(dm, fpoint, &ft));
4097:       const PetscInt *farr = DMPolytopeTypeGetArrangement(ft, DMPolytopeTypeComposeOrientation(ft, arr[d * 2 + 1], ornt[d]));
4098:       PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, fpoint, PETSC_TRUE, &fconeSize, &fcone, &fornt));
4099:       for (fc = 0; fc < fconeSize; ++fc) {
4100:         const PetscInt cp = fcone[farr[fc * 2 + 0]];
4101:         const PetscInt co = farr[fc * 2 + 1];

4103:         for (i = 0; i < c; i += 2)
4104:           if (pts[i] == cp) break;
4105:         if (i == c) {
4106:           PetscCall(DMPlexGetCellType(dm, cp, &ft));
4107:           pts[c++] = cp;
4108:           pts[c++] = DMPolytopeTypeComposeOrientation(ft, co, fornt[farr[fc * 2 + 0]]);
4109:         }
4110:       }
4111:       PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, fpoint, PETSC_TRUE, &fconeSize, &fcone, &fornt));
4112:     }
4113:   }
4114:   *numPoints = c / 2;
4115:   *points    = pts;
4116:   PetscFunctionReturn(PETSC_SUCCESS);
4117: }

4119: PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
4120: {
4121:   DMPolytopeType ct;
4122:   PetscInt      *closure, *fifo;
4123:   PetscInt       closureSize = 0, fifoStart = 0, fifoSize = 0;
4124:   PetscInt       maxConeSize, maxSupportSize, coneSeries, supportSeries;
4125:   PetscInt       depth, maxSize;

4127:   PetscFunctionBeginHot;
4128:   PetscCall(DMPlexGetDepth(dm, &depth));
4129:   if (depth == 1) {
4130:     PetscCall(DMPlexGetTransitiveClosure_Depth1_Private(dm, p, ornt, useCone, numPoints, points));
4131:     PetscFunctionReturn(PETSC_SUCCESS);
4132:   }
4133:   PetscCall(DMPlexGetCellType(dm, p, &ct));
4134:   if (ct == DM_POLYTOPE_FV_GHOST || ct == DM_POLYTOPE_INTERIOR_GHOST || ct == DM_POLYTOPE_UNKNOWN || ct == DM_POLYTOPE_UNKNOWN_CELL || ct == DM_POLYTOPE_UNKNOWN_FACE) ct = DM_POLYTOPE_UNKNOWN;
4135:   if (DMPolytopeTypeIsHybrid(ct) && ct != DM_POLYTOPE_POINT_PRISM_TENSOR) {
4136:     PetscCall(DMPlexTransitiveClosure_Tensor_Internal(dm, p, ct, ornt, useCone, numPoints, points));
4137:     PetscFunctionReturn(PETSC_SUCCESS);
4138:   }
4139:   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
4140:   coneSeries    = (maxConeSize > 1) ? ((PetscPowInt(maxConeSize, depth + 1) - 1) / (maxConeSize - 1)) : depth + 1;
4141:   supportSeries = (maxSupportSize > 1) ? ((PetscPowInt(maxSupportSize, depth + 1) - 1) / (maxSupportSize - 1)) : depth + 1;
4142:   maxSize       = PetscMax(coneSeries, supportSeries);
4143:   PetscCall(DMGetWorkArray(dm, 3 * maxSize, MPIU_INT, &fifo));
4144:   if (*points) {
4145:     closure = *points;
4146:   } else PetscCall(DMGetWorkArray(dm, 2 * maxSize, MPIU_INT, &closure));
4147:   closure[closureSize++] = p;
4148:   closure[closureSize++] = ornt;
4149:   fifo[fifoSize++]       = p;
4150:   fifo[fifoSize++]       = ornt;
4151:   fifo[fifoSize++]       = ct;
4152:   /* Should kick out early when depth is reached, rather than checking all vertices for empty cones */
4153:   while (fifoSize - fifoStart) {
4154:     const PetscInt       q    = fifo[fifoStart++];
4155:     const PetscInt       o    = fifo[fifoStart++];
4156:     const DMPolytopeType qt   = (DMPolytopeType)fifo[fifoStart++];
4157:     const PetscInt      *qarr = DMPolytopeTypeGetArrangement(qt, o);
4158:     const PetscInt      *tmp, *tmpO = NULL;
4159:     PetscInt             tmpSize, t;

4161:     if (PetscDefined(USE_DEBUG)) {
4162:       PetscInt nO = DMPolytopeTypeGetNumArrangements(qt) / 2;
4163:       PetscCheck(!o || !(o >= nO || o < -nO), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ") for %s %" PetscInt_FMT, o, -nO, nO, DMPolytopeTypes[qt], q);
4164:     }
4165:     PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, q, useCone, &tmpSize, &tmp, &tmpO));
4166:     for (t = 0; t < tmpSize; ++t) {
4167:       const PetscInt ip = useCone && qarr ? qarr[t * 2] : t;
4168:       const PetscInt io = useCone && qarr ? qarr[t * 2 + 1] : 0;
4169:       const PetscInt cp = tmp[ip];
4170:       PetscCall(DMPlexGetCellType(dm, cp, &ct));
4171:       const PetscInt co = tmpO ? DMPolytopeTypeComposeOrientation(ct, io, tmpO[ip]) : 0;
4172:       PetscInt       c;

4174:       /* Check for duplicate */
4175:       for (c = 0; c < closureSize; c += 2) {
4176:         if (closure[c] == cp) break;
4177:       }
4178:       if (c == closureSize) {
4179:         closure[closureSize++] = cp;
4180:         closure[closureSize++] = co;
4181:         fifo[fifoSize++]       = cp;
4182:         fifo[fifoSize++]       = co;
4183:         fifo[fifoSize++]       = ct;
4184:       }
4185:     }
4186:     PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, q, useCone, &tmpSize, &tmp, &tmpO));
4187:   }
4188:   PetscCall(DMRestoreWorkArray(dm, 3 * maxSize, MPIU_INT, &fifo));
4189:   if (numPoints) *numPoints = closureSize / 2;
4190:   if (points) *points = closure;
4191:   PetscFunctionReturn(PETSC_SUCCESS);
4192: }

4194: /*@C
4195:   DMPlexGetTransitiveClosure - Return the points on the transitive closure of the in-edges or out-edges for this point in the DAG

4197:   Not Collective

4199:   Input Parameters:
4200: + dm      - The `DMPLEX`
4201: . p       - The mesh point
4202: - useCone - `PETSC_TRUE` for the closure, otherwise return the star

4204:   Input/Output Parameter:
4205: . points - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...];
4206:            if *points is `NULL` on input, internal storage will be returned, use `DMPlexRestoreTransitiveClosure()`,
4207:            otherwise the provided array is used to hold the values

4209:   Output Parameter:
4210: . numPoints - The number of points in the closure, so `points` is of size 2*`numPoints`

4212:   Level: beginner

4214:   Note:
4215:   If using internal storage (points is `NULL` on input), each call overwrites the last output.

4217:   Fortran Notes:
4218:   `points` must be declared with
4219: .vb
4220:   PetscInt, pointer :: points(:)
4221: .ve
4222:   and is always allocated by the function.

4224:   The `numPoints` argument is not present in the Fortran binding.

4226: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreTransitiveClosure()`, `DMPlexCreate()`, `DMPlexSetCone()`, `DMPlexSetChart()`, `DMPlexGetCone()`
4227: @*/
4228: PetscErrorCode DMPlexGetTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
4229: {
4230:   PetscFunctionBeginHot;
4232:   if (numPoints) PetscAssertPointer(numPoints, 4);
4233:   if (points) PetscAssertPointer(points, 5);
4234:   if (PetscDefined(USE_DEBUG)) {
4235:     PetscInt pStart, pEnd;
4236:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4237:     PetscCheck(p >= pStart && p < pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " is not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", p, pStart, pEnd);
4238:   }
4239:   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, p, 0, useCone, numPoints, points));
4240:   PetscFunctionReturn(PETSC_SUCCESS);
4241: }

4243: /*@C
4244:   DMPlexRestoreTransitiveClosure - Restore the array of points on the transitive closure of the in-edges or out-edges for this point in the DAG

4246:   Not Collective

4248:   Input Parameters:
4249: + dm        - The `DMPLEX`
4250: . p         - The mesh point
4251: . useCone   - `PETSC_TRUE` for the closure, otherwise return the star
4252: . numPoints - The number of points in the closure, so points[] is of size 2*`numPoints`
4253: - points    - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...]

4255:   Level: beginner

4257:   Note:
4258:   If not using internal storage (points is not `NULL` on input), this call is unnecessary

4260: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetTransitiveClosure()`, `DMPlexCreate()`, `DMPlexSetCone()`, `DMPlexSetChart()`, `DMPlexGetCone()`
4261: @*/
4262: PetscErrorCode DMPlexRestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
4263: {
4264:   PetscFunctionBeginHot;
4266:   if (numPoints) *numPoints = 0;
4267:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, points));
4268:   PetscFunctionReturn(PETSC_SUCCESS);
4269: }

4271: /*@
4272:   DMPlexGetMaxSizes - Return the maximum number of in-edges (cone) and out-edges (support) for any point in the DAG

4274:   Not Collective

4276:   Input Parameter:
4277: . dm - The `DMPLEX`

4279:   Output Parameters:
4280: + maxConeSize    - The maximum number of in-edges
4281: - maxSupportSize - The maximum number of out-edges

4283:   Level: beginner

4285: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`
4286: @*/
4287: PetscErrorCode DMPlexGetMaxSizes(DM dm, PetscInt *maxConeSize, PetscInt *maxSupportSize)
4288: {
4289:   DM_Plex *mesh = (DM_Plex *)dm->data;

4291:   PetscFunctionBegin;
4293:   if (maxConeSize) PetscCall(PetscSectionGetMaxDof(mesh->coneSection, maxConeSize));
4294:   if (maxSupportSize) PetscCall(PetscSectionGetMaxDof(mesh->supportSection, maxSupportSize));
4295:   PetscFunctionReturn(PETSC_SUCCESS);
4296: }

4298: PetscErrorCode DMSetUp_Plex(DM dm)
4299: {
4300:   DM_Plex *mesh = (DM_Plex *)dm->data;
4301:   PetscInt size, maxSupportSize;

4303:   PetscFunctionBegin;
4305:   PetscCall(PetscSectionSetUp(mesh->coneSection));
4306:   PetscCall(PetscSectionGetStorageSize(mesh->coneSection, &size));
4307:   PetscCall(PetscMalloc1(size, &mesh->cones));
4308:   PetscCall(PetscCalloc1(size, &mesh->coneOrientations));
4309:   PetscCall(PetscSectionGetMaxDof(mesh->supportSection, &maxSupportSize));
4310:   if (maxSupportSize) {
4311:     PetscCall(PetscSectionSetUp(mesh->supportSection));
4312:     PetscCall(PetscSectionGetStorageSize(mesh->supportSection, &size));
4313:     PetscCall(PetscMalloc1(size, &mesh->supports));
4314:   }
4315:   PetscFunctionReturn(PETSC_SUCCESS);
4316: }

4318: PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
4319: {
4320:   PetscFunctionBegin;
4321:   if (subdm) PetscCall(DMClone(dm, subdm));
4322:   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, NULL, NULL, is, subdm));
4323:   if (subdm) (*subdm)->useNatural = dm->useNatural;
4324:   if (dm->useNatural && dm->sfMigration) {
4325:     PetscSF sfNatural;

4327:     (*subdm)->sfMigration = dm->sfMigration;
4328:     PetscCall(PetscObjectReference((PetscObject)dm->sfMigration));
4329:     PetscCall(DMPlexCreateGlobalToNaturalSF(*subdm, NULL, (*subdm)->sfMigration, &sfNatural));
4330:     (*subdm)->sfNatural = sfNatural;
4331:   }
4332:   PetscFunctionReturn(PETSC_SUCCESS);
4333: }

4335: PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm)
4336: {
4337:   PetscInt i = 0;

4339:   PetscFunctionBegin;
4340:   PetscCall(DMClone(dms[0], superdm));
4341:   PetscCall(DMCreateSectionSuperDM(dms, len, is, superdm));
4342:   (*superdm)->useNatural = PETSC_FALSE;
4343:   for (i = 0; i < len; i++) {
4344:     if (dms[i]->useNatural && dms[i]->sfMigration) {
4345:       PetscSF sfNatural;

4347:       (*superdm)->sfMigration = dms[i]->sfMigration;
4348:       PetscCall(PetscObjectReference((PetscObject)dms[i]->sfMigration));
4349:       (*superdm)->useNatural = PETSC_TRUE;
4350:       PetscCall(DMPlexCreateGlobalToNaturalSF(*superdm, NULL, (*superdm)->sfMigration, &sfNatural));
4351:       (*superdm)->sfNatural = sfNatural;
4352:       break;
4353:     }
4354:   }
4355:   PetscFunctionReturn(PETSC_SUCCESS);
4356: }

4358: /*@
4359:   DMPlexSymmetrize - Create support (out-edge) information from cone (in-edge) information

4361:   Not Collective

4363:   Input Parameter:
4364: . dm - The `DMPLEX`

4366:   Level: beginner

4368:   Note:
4369:   This should be called after all calls to `DMPlexSetCone()`

4371: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMPlexSetCone()`
4372: @*/
4373: PetscErrorCode DMPlexSymmetrize(DM dm)
4374: {
4375:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4376:   PetscInt *offsets;
4377:   PetscInt  supportSize;
4378:   PetscInt  pStart, pEnd, p;

4380:   PetscFunctionBegin;
4382:   PetscCheck(!mesh->supports, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Supports were already setup in this DMPlex");
4383:   PetscCall(PetscLogEventBegin(DMPLEX_Symmetrize, dm, 0, 0, 0));
4384:   /* Calculate support sizes */
4385:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4386:   for (p = pStart; p < pEnd; ++p) {
4387:     PetscInt dof, off, c;

4389:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
4390:     PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
4391:     for (c = off; c < off + dof; ++c) PetscCall(PetscSectionAddDof(mesh->supportSection, mesh->cones[c], 1));
4392:   }
4393:   PetscCall(PetscSectionSetUp(mesh->supportSection));
4394:   /* Calculate supports */
4395:   PetscCall(PetscSectionGetStorageSize(mesh->supportSection, &supportSize));
4396:   PetscCall(PetscMalloc1(supportSize, &mesh->supports));
4397:   PetscCall(PetscCalloc1(pEnd - pStart, &offsets));
4398:   for (p = pStart; p < pEnd; ++p) {
4399:     PetscInt dof, off, c;

4401:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
4402:     PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
4403:     for (c = off; c < off + dof; ++c) {
4404:       const PetscInt q = mesh->cones[c];
4405:       PetscInt       offS;

4407:       PetscCall(PetscSectionGetOffset(mesh->supportSection, q, &offS));

4409:       mesh->supports[offS + offsets[q]] = p;
4410:       ++offsets[q];
4411:     }
4412:   }
4413:   PetscCall(PetscFree(offsets));
4414:   PetscCall(PetscLogEventEnd(DMPLEX_Symmetrize, dm, 0, 0, 0));
4415:   PetscFunctionReturn(PETSC_SUCCESS);
4416: }

4418: static PetscErrorCode DMPlexCreateDepthStratum(DM dm, DMLabel label, PetscInt depth, PetscInt pStart, PetscInt pEnd)
4419: {
4420:   IS stratumIS;

4422:   PetscFunctionBegin;
4423:   if (pStart >= pEnd) PetscFunctionReturn(PETSC_SUCCESS);
4424:   if (PetscDefined(USE_DEBUG)) {
4425:     PetscInt  qStart, qEnd, numLevels, level;
4426:     PetscBool overlap = PETSC_FALSE;
4427:     PetscCall(DMLabelGetNumValues(label, &numLevels));
4428:     for (level = 0; level < numLevels; level++) {
4429:       PetscCall(DMLabelGetStratumBounds(label, level, &qStart, &qEnd));
4430:       if ((pStart >= qStart && pStart < qEnd) || (pEnd > qStart && pEnd <= qEnd)) {
4431:         overlap = PETSC_TRUE;
4432:         break;
4433:       }
4434:     }
4435:     PetscCheck(!overlap, PETSC_COMM_SELF, PETSC_ERR_PLIB, "New depth %" PetscInt_FMT " range [%" PetscInt_FMT ",%" PetscInt_FMT ") overlaps with depth %" PetscInt_FMT " range [%" PetscInt_FMT ",%" PetscInt_FMT ")", depth, pStart, pEnd, level, qStart, qEnd);
4436:   }
4437:   PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &stratumIS));
4438:   PetscCall(DMLabelSetStratumIS(label, depth, stratumIS));
4439:   PetscCall(ISDestroy(&stratumIS));
4440:   PetscFunctionReturn(PETSC_SUCCESS);
4441: }

4443: static PetscErrorCode DMPlexStratify_CellType_Private(DM dm, DMLabel label)
4444: {
4445:   PetscInt *pMin, *pMax;
4446:   PetscInt  pStart, pEnd;
4447:   PetscInt  dmin = PETSC_INT_MAX, dmax = PETSC_INT_MIN;

4449:   PetscFunctionBegin;
4450:   {
4451:     DMLabel label2;

4453:     PetscCall(DMPlexGetCellTypeLabel(dm, &label2));
4454:     PetscCall(PetscObjectViewFromOptions((PetscObject)label2, NULL, "-ct_view"));
4455:   }
4456:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4457:   for (PetscInt p = pStart; p < pEnd; ++p) {
4458:     DMPolytopeType ct;

4460:     PetscCall(DMPlexGetCellType(dm, p, &ct));
4461:     dmin = PetscMin(DMPolytopeTypeGetDim(ct), dmin);
4462:     dmax = PetscMax(DMPolytopeTypeGetDim(ct), dmax);
4463:   }
4464:   PetscCall(PetscMalloc2(dmax + 1, &pMin, dmax + 1, &pMax));
4465:   for (PetscInt d = dmin; d <= dmax; ++d) {
4466:     pMin[d] = PETSC_INT_MAX;
4467:     pMax[d] = PETSC_INT_MIN;
4468:   }
4469:   for (PetscInt p = pStart; p < pEnd; ++p) {
4470:     DMPolytopeType ct;
4471:     PetscInt       d;

4473:     PetscCall(DMPlexGetCellType(dm, p, &ct));
4474:     d       = DMPolytopeTypeGetDim(ct);
4475:     pMin[d] = PetscMin(p, pMin[d]);
4476:     pMax[d] = PetscMax(p, pMax[d]);
4477:   }
4478:   for (PetscInt d = dmin; d <= dmax; ++d) {
4479:     if (pMin[d] > pMax[d]) continue;
4480:     PetscCall(DMPlexCreateDepthStratum(dm, label, d, pMin[d], pMax[d] + 1));
4481:   }
4482:   PetscCall(PetscFree2(pMin, pMax));
4483:   PetscFunctionReturn(PETSC_SUCCESS);
4484: }

4486: static PetscErrorCode DMPlexStratify_Topological_Private(DM dm, DMLabel label)
4487: {
4488:   PetscInt pStart, pEnd;
4489:   PetscInt numRoots = 0, numLeaves = 0;

4491:   PetscFunctionBegin;
4492:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4493:   {
4494:     /* Initialize roots and count leaves */
4495:     PetscInt sMin = PETSC_INT_MAX;
4496:     PetscInt sMax = PETSC_INT_MIN;
4497:     PetscInt coneSize, supportSize;

4499:     for (PetscInt p = pStart; p < pEnd; ++p) {
4500:       PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
4501:       PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
4502:       if (!coneSize && supportSize) {
4503:         sMin = PetscMin(p, sMin);
4504:         sMax = PetscMax(p, sMax);
4505:         ++numRoots;
4506:       } else if (!supportSize && coneSize) {
4507:         ++numLeaves;
4508:       } else if (!supportSize && !coneSize) {
4509:         /* Isolated points */
4510:         sMin = PetscMin(p, sMin);
4511:         sMax = PetscMax(p, sMax);
4512:       }
4513:     }
4514:     PetscCall(DMPlexCreateDepthStratum(dm, label, 0, sMin, sMax + 1));
4515:   }

4517:   if (numRoots + numLeaves == (pEnd - pStart)) {
4518:     PetscInt sMin = PETSC_INT_MAX;
4519:     PetscInt sMax = PETSC_INT_MIN;
4520:     PetscInt coneSize, supportSize;

4522:     for (PetscInt p = pStart; p < pEnd; ++p) {
4523:       PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
4524:       PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
4525:       if (!supportSize && coneSize) {
4526:         sMin = PetscMin(p, sMin);
4527:         sMax = PetscMax(p, sMax);
4528:       }
4529:     }
4530:     PetscCall(DMPlexCreateDepthStratum(dm, label, 1, sMin, sMax + 1));
4531:   } else {
4532:     PetscInt level = 0;
4533:     PetscInt qStart, qEnd;

4535:     PetscCall(DMLabelGetStratumBounds(label, level, &qStart, &qEnd));
4536:     while (qEnd > qStart) {
4537:       PetscInt sMin = PETSC_INT_MAX;
4538:       PetscInt sMax = PETSC_INT_MIN;

4540:       for (PetscInt q = qStart; q < qEnd; ++q) {
4541:         const PetscInt *support;
4542:         PetscInt        supportSize;

4544:         PetscCall(DMPlexGetSupportSize(dm, q, &supportSize));
4545:         PetscCall(DMPlexGetSupport(dm, q, &support));
4546:         for (PetscInt s = 0; s < supportSize; ++s) {
4547:           sMin = PetscMin(support[s], sMin);
4548:           sMax = PetscMax(support[s], sMax);
4549:         }
4550:       }
4551:       PetscCall(DMLabelGetNumValues(label, &level));
4552:       PetscCall(DMPlexCreateDepthStratum(dm, label, level, sMin, sMax + 1));
4553:       PetscCall(DMLabelGetStratumBounds(label, level, &qStart, &qEnd));
4554:     }
4555:   }
4556:   PetscFunctionReturn(PETSC_SUCCESS);
4557: }

4559: /*@
4560:   DMPlexStratify - Computes the strata for all points in the `DMPLEX`

4562:   Collective

4564:   Input Parameter:
4565: . dm - The `DMPLEX`

4567:   Level: beginner

4569:   Notes:
4570:   The strata group all points of the same grade, and this function calculates the strata. This
4571:   grade can be seen as the height (or depth) of the point in the DAG.

4573:   The DAG for most topologies is a graded poset (https://en.wikipedia.org/wiki/Graded_poset), and
4574:   can be illustrated by a Hasse Diagram (https://en.wikipedia.org/wiki/Hasse_diagram).
4575:   Concretely, `DMPlexStratify()` creates a new label named "depth" containing the depth in the DAG of each point. For cell-vertex
4576:   meshes, vertices are depth 0 and cells are depth 1. For fully interpolated meshes, depth 0 for vertices, 1 for edges, and so on
4577:   until cells have depth equal to the dimension of the mesh. The depth label can be accessed through `DMPlexGetDepthLabel()` or `DMPlexGetDepthStratum()`, or
4578:   manually via `DMGetLabel()`.  The height is defined implicitly by height = maxDimension - depth, and can be accessed
4579:   via `DMPlexGetHeightStratum()`.  For example, cells have height 0 and faces have height 1.

4581:   The depth of a point is calculated by executing a breadth-first search (BFS) on the DAG. This could produce surprising results
4582:   if run on a partially interpolated mesh, meaning one that had some edges and faces, but not others. For example, suppose that
4583:   we had a mesh consisting of one triangle (c0) and three vertices (v0, v1, v2), and only one edge is on the boundary so we choose
4584:   to interpolate only that one (e0), so that
4585: .vb
4586:   cone(c0) = {e0, v2}
4587:   cone(e0) = {v0, v1}
4588: .ve
4589:   If `DMPlexStratify()` is run on this mesh, it will give depths
4590: .vb
4591:    depth 0 = {v0, v1, v2}
4592:    depth 1 = {e0, c0}
4593: .ve
4594:   where the triangle has been given depth 1, instead of 2, because it is reachable from vertex v2.

4596:   `DMPlexStratify()` should be called after all calls to `DMPlexSymmetrize()`

4598: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSymmetrize()`, `DMPlexComputeCellTypes()`
4599: @*/
4600: PetscErrorCode DMPlexStratify(DM dm)
4601: {
4602:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4603:   DMLabel   label;
4604:   PetscBool flg = PETSC_FALSE;

4606:   PetscFunctionBegin;
4608:   PetscCall(PetscLogEventBegin(DMPLEX_Stratify, dm, 0, 0, 0));

4610:   // Create depth label
4611:   PetscCall(DMRemoveLabel(dm, "depth", NULL));
4612:   PetscCall(DMCreateLabel(dm, "depth"));
4613:   PetscCall(DMPlexGetDepthLabel(dm, &label));

4615:   PetscCall(PetscOptionsGetBool(NULL, dm->hdr.prefix, "-dm_plex_stratify_celltype", &flg, NULL));
4616:   if (flg) PetscCall(DMPlexStratify_CellType_Private(dm, label));
4617:   else PetscCall(DMPlexStratify_Topological_Private(dm, label));

4619:   { /* just in case there is an empty process */
4620:     PetscInt numValues, maxValues = 0, v;

4622:     PetscCall(DMLabelGetNumValues(label, &numValues));
4623:     PetscCallMPI(MPIU_Allreduce(&numValues, &maxValues, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
4624:     for (v = numValues; v < maxValues; v++) PetscCall(DMLabelAddStratum(label, v));
4625:   }
4626:   PetscCall(PetscObjectStateGet((PetscObject)label, &mesh->depthState));
4627:   PetscCall(PetscLogEventEnd(DMPLEX_Stratify, dm, 0, 0, 0));
4628:   PetscFunctionReturn(PETSC_SUCCESS);
4629: }

4631: PetscErrorCode DMPlexComputeCellType_Internal(DM dm, PetscInt p, PetscInt pdepth, DMPolytopeType *pt)
4632: {
4633:   DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
4634:   PetscInt       dim, depth, pheight, coneSize;
4635:   PetscBool      preferTensor;

4637:   PetscFunctionBeginHot;
4638:   PetscCall(DMGetDimension(dm, &dim));
4639:   PetscCall(DMPlexGetDepth(dm, &depth));
4640:   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
4641:   PetscCall(DMPlexGetInterpolatePreferTensor(dm, &preferTensor));
4642:   pheight = depth - pdepth;
4643:   if (depth <= 1) {
4644:     switch (pdepth) {
4645:     case 0:
4646:       ct = DM_POLYTOPE_POINT;
4647:       break;
4648:     case 1:
4649:       switch (coneSize) {
4650:       case 2:
4651:         ct = DM_POLYTOPE_SEGMENT;
4652:         break;
4653:       case 3:
4654:         ct = DM_POLYTOPE_TRIANGLE;
4655:         break;
4656:       case 4:
4657:         switch (dim) {
4658:         case 2:
4659:           ct = DM_POLYTOPE_QUADRILATERAL;
4660:           break;
4661:         case 3:
4662:           ct = DM_POLYTOPE_TETRAHEDRON;
4663:           break;
4664:         default:
4665:           break;
4666:         }
4667:         break;
4668:       case 5:
4669:         ct = DM_POLYTOPE_PYRAMID;
4670:         break;
4671:       case 6:
4672:         ct = preferTensor ? DM_POLYTOPE_TRI_PRISM_TENSOR : DM_POLYTOPE_TRI_PRISM;
4673:         break;
4674:       case 8:
4675:         ct = DM_POLYTOPE_HEXAHEDRON;
4676:         break;
4677:       default:
4678:         break;
4679:       }
4680:     }
4681:   } else {
4682:     if (pdepth == 0) {
4683:       ct = DM_POLYTOPE_POINT;
4684:     } else if (pheight == 0) {
4685:       switch (dim) {
4686:       case 1:
4687:         switch (coneSize) {
4688:         case 2:
4689:           ct = DM_POLYTOPE_SEGMENT;
4690:           break;
4691:         default:
4692:           break;
4693:         }
4694:         break;
4695:       case 2:
4696:         switch (coneSize) {
4697:         case 3:
4698:           ct = DM_POLYTOPE_TRIANGLE;
4699:           break;
4700:         case 4:
4701:           ct = DM_POLYTOPE_QUADRILATERAL;
4702:           break;
4703:         default:
4704:           break;
4705:         }
4706:         break;
4707:       case 3:
4708:         switch (coneSize) {
4709:         case 4:
4710:           ct = DM_POLYTOPE_TETRAHEDRON;
4711:           break;
4712:         case 5: {
4713:           const PetscInt *cone;
4714:           PetscInt        faceConeSize;

4716:           PetscCall(DMPlexGetCone(dm, p, &cone));
4717:           PetscCall(DMPlexGetConeSize(dm, cone[0], &faceConeSize));
4718:           switch (faceConeSize) {
4719:           case 3:
4720:             ct = preferTensor ? DM_POLYTOPE_TRI_PRISM_TENSOR : DM_POLYTOPE_TRI_PRISM;
4721:             break;
4722:           case 4:
4723:             ct = DM_POLYTOPE_PYRAMID;
4724:             break;
4725:           }
4726:         } break;
4727:         case 6:
4728:           ct = DM_POLYTOPE_HEXAHEDRON;
4729:           break;
4730:         default:
4731:           break;
4732:         }
4733:         break;
4734:       default:
4735:         break;
4736:       }
4737:     } else if (pheight > 0) {
4738:       switch (coneSize) {
4739:       case 2:
4740:         ct = DM_POLYTOPE_SEGMENT;
4741:         break;
4742:       case 3:
4743:         ct = DM_POLYTOPE_TRIANGLE;
4744:         break;
4745:       case 4:
4746:         ct = DM_POLYTOPE_QUADRILATERAL;
4747:         break;
4748:       default:
4749:         break;
4750:       }
4751:     }
4752:   }
4753:   *pt = ct;
4754:   PetscFunctionReturn(PETSC_SUCCESS);
4755: }

4757: /*@
4758:   DMPlexComputeCellTypes - Infer the polytope type of every cell using its dimension and cone size.

4760:   Collective

4762:   Input Parameter:
4763: . dm - The `DMPLEX`

4765:   Level: developer

4767:   Note:
4768:   This function is normally called automatically when a cell type is requested. It creates an
4769:   internal `DMLabel` named "celltype" which can be directly accessed using `DMGetLabel()`. A user may disable
4770:   automatic creation by creating the label manually, using `DMCreateLabel`(dm, "celltype").

4772:   `DMPlexComputeCellTypes()` should be called after all calls to `DMPlexSymmetrize()` and `DMPlexStratify()`

4774: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSymmetrize()`, `DMPlexStratify()`, `DMGetLabel()`, `DMCreateLabel()`
4775: @*/
4776: PetscErrorCode DMPlexComputeCellTypes(DM dm)
4777: {
4778:   DM_Plex *mesh;
4779:   DMLabel  ctLabel;
4780:   PetscInt pStart, pEnd, p;

4782:   PetscFunctionBegin;
4784:   mesh = (DM_Plex *)dm->data;
4785:   PetscCall(DMCreateLabel(dm, "celltype"));
4786:   PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel));
4787:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4788:   PetscCall(PetscFree(mesh->cellTypes));
4789:   PetscCall(PetscMalloc1(pEnd - pStart, &mesh->cellTypes));
4790:   for (p = pStart; p < pEnd; ++p) {
4791:     DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
4792:     PetscInt       pdepth;

4794:     PetscCall(DMPlexGetPointDepth(dm, p, &pdepth));
4795:     PetscCall(DMPlexComputeCellType_Internal(dm, p, pdepth, &ct));
4796:     PetscCheck(ct != DM_POLYTOPE_UNKNOWN && ct != DM_POLYTOPE_UNKNOWN_CELL && ct != DM_POLYTOPE_UNKNOWN_FACE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Point %" PetscInt_FMT " has invalid celltype (%s)", p, DMPolytopeTypes[ct]);
4797:     PetscCall(DMLabelSetValue(ctLabel, p, ct));
4798:     mesh->cellTypes[p - pStart].value_as_uint8 = (uint8_t)ct;
4799:   }
4800:   PetscCall(PetscObjectStateGet((PetscObject)ctLabel, &mesh->celltypeState));
4801:   PetscCall(PetscObjectViewFromOptions((PetscObject)ctLabel, NULL, "-dm_plex_celltypes_view"));
4802:   PetscFunctionReturn(PETSC_SUCCESS);
4803: }

4805: /*@C
4806:   DMPlexGetJoin - Get an array for the join of the set of points

4808:   Not Collective

4810:   Input Parameters:
4811: + dm        - The `DMPLEX` object
4812: . numPoints - The number of input points for the join
4813: - points    - The input points

4815:   Output Parameters:
4816: + numCoveredPoints - The number of points in the join
4817: - coveredPoints    - The points in the join

4819:   Level: intermediate

4821:   Note:
4822:   Currently, this is restricted to a single level join

4824:   Fortran Notes:
4825:   `converedPoints` must be declared with
4826: .vb
4827:   PetscInt, pointer :: coveredPints(:)
4828: .ve

4830:   The `numCoveredPoints` argument is not present in the Fortran binding.

4832: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreJoin()`, `DMPlexGetMeet()`
4833: @*/
4834: PetscErrorCode DMPlexGetJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4835: {
4836:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4837:   PetscInt *join[2];
4838:   PetscInt  joinSize, i = 0;
4839:   PetscInt  dof, off, p, c, m;
4840:   PetscInt  maxSupportSize;

4842:   PetscFunctionBegin;
4844:   PetscAssertPointer(points, 3);
4845:   PetscAssertPointer(numCoveredPoints, 4);
4846:   PetscAssertPointer(coveredPoints, 5);
4847:   PetscCall(PetscSectionGetMaxDof(mesh->supportSection, &maxSupportSize));
4848:   PetscCall(DMGetWorkArray(dm, maxSupportSize, MPIU_INT, &join[0]));
4849:   PetscCall(DMGetWorkArray(dm, maxSupportSize, MPIU_INT, &join[1]));
4850:   /* Copy in support of first point */
4851:   PetscCall(PetscSectionGetDof(mesh->supportSection, points[0], &dof));
4852:   PetscCall(PetscSectionGetOffset(mesh->supportSection, points[0], &off));
4853:   for (joinSize = 0; joinSize < dof; ++joinSize) join[i][joinSize] = mesh->supports[off + joinSize];
4854:   /* Check each successive support */
4855:   for (p = 1; p < numPoints; ++p) {
4856:     PetscInt newJoinSize = 0;

4858:     PetscCall(PetscSectionGetDof(mesh->supportSection, points[p], &dof));
4859:     PetscCall(PetscSectionGetOffset(mesh->supportSection, points[p], &off));
4860:     for (c = 0; c < dof; ++c) {
4861:       const PetscInt point = mesh->supports[off + c];

4863:       for (m = 0; m < joinSize; ++m) {
4864:         if (point == join[i][m]) {
4865:           join[1 - i][newJoinSize++] = point;
4866:           break;
4867:         }
4868:       }
4869:     }
4870:     joinSize = newJoinSize;
4871:     i        = 1 - i;
4872:   }
4873:   *numCoveredPoints = joinSize;
4874:   *coveredPoints    = join[i];
4875:   PetscCall(DMRestoreWorkArray(dm, maxSupportSize, MPIU_INT, &join[1 - i]));
4876:   PetscFunctionReturn(PETSC_SUCCESS);
4877: }

4879: /*@C
4880:   DMPlexRestoreJoin - Restore an array for the join of the set of points obtained with `DMPlexGetJoin()`

4882:   Not Collective

4884:   Input Parameters:
4885: + dm        - The `DMPLEX` object
4886: . numPoints - The number of input points for the join
4887: - points    - The input points

4889:   Output Parameters:
4890: + numCoveredPoints - The number of points in the join
4891: - coveredPoints    - The points in the join

4893:   Level: intermediate

4895:   Fortran Notes:
4896:   The `numCoveredPoints` argument is not present in the Fortran binding since it is internal to the array.

4898: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetJoin()`, `DMPlexGetFullJoin()`, `DMPlexGetMeet()`
4899: @*/
4900: PetscErrorCode DMPlexRestoreJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4901: {
4902:   PetscFunctionBegin;
4904:   if (points) PetscAssertPointer(points, 3);
4905:   if (numCoveredPoints) PetscAssertPointer(numCoveredPoints, 4);
4906:   PetscAssertPointer(coveredPoints, 5);
4907:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, (void *)coveredPoints));
4908:   if (numCoveredPoints) *numCoveredPoints = 0;
4909:   PetscFunctionReturn(PETSC_SUCCESS);
4910: }

4912: /*@C
4913:   DMPlexGetFullJoin - Get an array for the join of the set of points

4915:   Not Collective

4917:   Input Parameters:
4918: + dm        - The `DMPLEX` object
4919: . numPoints - The number of input points for the join
4920: - points    - The input points, its length is `numPoints`

4922:   Output Parameters:
4923: + numCoveredPoints - The number of points in the join
4924: - coveredPoints    - The points in the join, its length is `numCoveredPoints`

4926:   Level: intermediate

4928:   Fortran Notes:
4929:   `points` and `converedPoints` must be declared with
4930: .vb
4931:   PetscInt, pointer :: points(:)
4932:   PetscInt, pointer :: coveredPints(:)
4933: .ve

4935:   The `numCoveredPoints` argument is not present in the Fortran binding.

4937: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetJoin()`, `DMPlexRestoreJoin()`, `DMPlexGetMeet()`
4938: @*/
4939: PetscErrorCode DMPlexGetFullJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4940: {
4941:   PetscInt *offsets, **closures;
4942:   PetscInt *join[2];
4943:   PetscInt  depth = 0, maxSize, joinSize = 0, i = 0;
4944:   PetscInt  p, d, c, m, ms;

4946:   PetscFunctionBegin;
4948:   PetscAssertPointer(points, 3);
4949:   PetscAssertPointer(numCoveredPoints, 4);
4950:   PetscAssertPointer(coveredPoints, 5);

4952:   PetscCall(DMPlexGetDepth(dm, &depth));
4953:   PetscCall(PetscCalloc1(numPoints, &closures));
4954:   PetscCall(DMGetWorkArray(dm, numPoints * (depth + 2), MPIU_INT, &offsets));
4955:   PetscCall(DMPlexGetMaxSizes(dm, NULL, &ms));
4956:   maxSize = (ms > 1) ? ((PetscPowInt(ms, depth + 1) - 1) / (ms - 1)) : depth + 1;
4957:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &join[0]));
4958:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &join[1]));

4960:   for (p = 0; p < numPoints; ++p) {
4961:     PetscInt closureSize;

4963:     PetscCall(DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &closureSize, &closures[p]));

4965:     offsets[p * (depth + 2) + 0] = 0;
4966:     for (d = 0; d < depth + 1; ++d) {
4967:       PetscInt pStart, pEnd, i;

4969:       PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
4970:       for (i = offsets[p * (depth + 2) + d]; i < closureSize; ++i) {
4971:         if ((pStart > closures[p][i * 2]) || (pEnd <= closures[p][i * 2])) {
4972:           offsets[p * (depth + 2) + d + 1] = i;
4973:           break;
4974:         }
4975:       }
4976:       if (i == closureSize) offsets[p * (depth + 2) + d + 1] = i;
4977:     }
4978:     PetscCheck(offsets[p * (depth + 2) + depth + 1] == closureSize, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Total size of closure %" PetscInt_FMT " should be %" PetscInt_FMT, offsets[p * (depth + 2) + depth + 1], closureSize);
4979:   }
4980:   for (d = 0; d < depth + 1; ++d) {
4981:     PetscInt dof;

4983:     /* Copy in support of first point */
4984:     dof = offsets[d + 1] - offsets[d];
4985:     for (joinSize = 0; joinSize < dof; ++joinSize) join[i][joinSize] = closures[0][(offsets[d] + joinSize) * 2];
4986:     /* Check each successive cone */
4987:     for (p = 1; p < numPoints && joinSize; ++p) {
4988:       PetscInt newJoinSize = 0;

4990:       dof = offsets[p * (depth + 2) + d + 1] - offsets[p * (depth + 2) + d];
4991:       for (c = 0; c < dof; ++c) {
4992:         const PetscInt point = closures[p][(offsets[p * (depth + 2) + d] + c) * 2];

4994:         for (m = 0; m < joinSize; ++m) {
4995:           if (point == join[i][m]) {
4996:             join[1 - i][newJoinSize++] = point;
4997:             break;
4998:           }
4999:         }
5000:       }
5001:       joinSize = newJoinSize;
5002:       i        = 1 - i;
5003:     }
5004:     if (joinSize) break;
5005:   }
5006:   *numCoveredPoints = joinSize;
5007:   *coveredPoints    = join[i];
5008:   for (p = 0; p < numPoints; ++p) PetscCall(DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, NULL, &closures[p]));
5009:   PetscCall(PetscFree(closures));
5010:   PetscCall(DMRestoreWorkArray(dm, numPoints * (depth + 2), MPIU_INT, &offsets));
5011:   PetscCall(DMRestoreWorkArray(dm, ms, MPIU_INT, &join[1 - i]));
5012:   PetscFunctionReturn(PETSC_SUCCESS);
5013: }

5015: /*@C
5016:   DMPlexGetMeet - Get an array for the meet of the set of points

5018:   Not Collective

5020:   Input Parameters:
5021: + dm        - The `DMPLEX` object
5022: . numPoints - The number of input points for the meet
5023: - points    - The input points, of length `numPoints`

5025:   Output Parameters:
5026: + numCoveringPoints - The number of points in the meet
5027: - coveringPoints    - The points in the meet, of length `numCoveringPoints`

5029:   Level: intermediate

5031:   Note:
5032:   Currently, this is restricted to a single level meet

5034:   Fortran Notes:
5035:   `coveringPoints` must be declared with
5036: .vb
5037:   PetscInt, pointer :: coveringPoints(:)
5038: .ve

5040:   The `numCoveredPoints` argument is not present in the Fortran binding since it is internal to the array.

5042: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreMeet()`, `DMPlexGetJoin()`
5043: @*/
5044: PetscErrorCode DMPlexGetMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveringPoints, const PetscInt *coveringPoints[])
5045: {
5046:   DM_Plex  *mesh = (DM_Plex *)dm->data;
5047:   PetscInt *meet[2];
5048:   PetscInt  meetSize, i = 0;
5049:   PetscInt  dof, off, p, c, m;
5050:   PetscInt  maxConeSize;

5052:   PetscFunctionBegin;
5054:   PetscAssertPointer(points, 3);
5055:   PetscAssertPointer(numCoveringPoints, 4);
5056:   PetscAssertPointer(coveringPoints, 5);
5057:   PetscCall(PetscSectionGetMaxDof(mesh->coneSection, &maxConeSize));
5058:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &meet[0]));
5059:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &meet[1]));
5060:   /* Copy in cone of first point */
5061:   PetscCall(PetscSectionGetDof(mesh->coneSection, points[0], &dof));
5062:   PetscCall(PetscSectionGetOffset(mesh->coneSection, points[0], &off));
5063:   for (meetSize = 0; meetSize < dof; ++meetSize) meet[i][meetSize] = mesh->cones[off + meetSize];
5064:   /* Check each successive cone */
5065:   for (p = 1; p < numPoints; ++p) {
5066:     PetscInt newMeetSize = 0;

5068:     PetscCall(PetscSectionGetDof(mesh->coneSection, points[p], &dof));
5069:     PetscCall(PetscSectionGetOffset(mesh->coneSection, points[p], &off));
5070:     for (c = 0; c < dof; ++c) {
5071:       const PetscInt point = mesh->cones[off + c];

5073:       for (m = 0; m < meetSize; ++m) {
5074:         if (point == meet[i][m]) {
5075:           meet[1 - i][newMeetSize++] = point;
5076:           break;
5077:         }
5078:       }
5079:     }
5080:     meetSize = newMeetSize;
5081:     i        = 1 - i;
5082:   }
5083:   *numCoveringPoints = meetSize;
5084:   *coveringPoints    = meet[i];
5085:   PetscCall(DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &meet[1 - i]));
5086:   PetscFunctionReturn(PETSC_SUCCESS);
5087: }

5089: /*@C
5090:   DMPlexRestoreMeet - Restore an array for the meet of the set of points obtained with `DMPlexGetMeet()`

5092:   Not Collective

5094:   Input Parameters:
5095: + dm        - The `DMPLEX` object
5096: . numPoints - The number of input points for the meet
5097: - points    - The input points

5099:   Output Parameters:
5100: + numCoveredPoints - The number of points in the meet
5101: - coveredPoints    - The points in the meet

5103:   Level: intermediate

5105:   Fortran Notes:
5106:   The `numCoveredPoints` argument is not present in the Fortran binding since it is internal to the array.

5108: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetMeet()`, `DMPlexGetFullMeet()`, `DMPlexGetJoin()`
5109: @*/
5110: PetscErrorCode DMPlexRestoreMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
5111: {
5112:   PetscFunctionBegin;
5114:   if (points) PetscAssertPointer(points, 3);
5115:   if (numCoveredPoints) PetscAssertPointer(numCoveredPoints, 4);
5116:   PetscAssertPointer(coveredPoints, 5);
5117:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, (void *)coveredPoints));
5118:   if (numCoveredPoints) *numCoveredPoints = 0;
5119:   PetscFunctionReturn(PETSC_SUCCESS);
5120: }

5122: /*@C
5123:   DMPlexGetFullMeet - Get an array for the meet of the set of points

5125:   Not Collective

5127:   Input Parameters:
5128: + dm        - The `DMPLEX` object
5129: . numPoints - The number of input points for the meet
5130: - points    - The input points, of length  `numPoints`

5132:   Output Parameters:
5133: + numCoveredPoints - The number of points in the meet
5134: - coveredPoints    - The points in the meet, of length  `numCoveredPoints`

5136:   Level: intermediate

5138:   Fortran Notes:
5139:   `points` and `coveredPoints` must be declared with
5140: .vb
5141:   PetscInt, pointer :: points(:)
5142:   PetscInt, pointer :: coveredPoints(:)
5143: .ve

5145:   The `numCoveredPoints` argument is not present in the Fortran binding since it is internal to the array.

5147: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetMeet()`, `DMPlexRestoreMeet()`, `DMPlexGetJoin()`
5148: @*/
5149: PetscErrorCode DMPlexGetFullMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
5150: {
5151:   PetscInt *offsets, **closures;
5152:   PetscInt *meet[2];
5153:   PetscInt  height = 0, maxSize, meetSize = 0, i = 0;
5154:   PetscInt  p, h, c, m, mc;

5156:   PetscFunctionBegin;
5158:   PetscAssertPointer(points, 3);
5159:   PetscAssertPointer(numCoveredPoints, 4);
5160:   PetscAssertPointer(coveredPoints, 5);

5162:   PetscCall(DMPlexGetDepth(dm, &height));
5163:   PetscCall(PetscMalloc1(numPoints, &closures));
5164:   PetscCall(DMGetWorkArray(dm, numPoints * (height + 2), MPIU_INT, &offsets));
5165:   PetscCall(DMPlexGetMaxSizes(dm, &mc, NULL));
5166:   maxSize = (mc > 1) ? ((PetscPowInt(mc, height + 1) - 1) / (mc - 1)) : height + 1;
5167:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &meet[0]));
5168:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &meet[1]));

5170:   for (p = 0; p < numPoints; ++p) {
5171:     PetscInt closureSize;

5173:     PetscCall(DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closures[p]));

5175:     offsets[p * (height + 2) + 0] = 0;
5176:     for (h = 0; h < height + 1; ++h) {
5177:       PetscInt pStart, pEnd, i;

5179:       PetscCall(DMPlexGetHeightStratum(dm, h, &pStart, &pEnd));
5180:       for (i = offsets[p * (height + 2) + h]; i < closureSize; ++i) {
5181:         if ((pStart > closures[p][i * 2]) || (pEnd <= closures[p][i * 2])) {
5182:           offsets[p * (height + 2) + h + 1] = i;
5183:           break;
5184:         }
5185:       }
5186:       if (i == closureSize) offsets[p * (height + 2) + h + 1] = i;
5187:     }
5188:     PetscCheck(offsets[p * (height + 2) + height + 1] == closureSize, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Total size of closure %" PetscInt_FMT " should be %" PetscInt_FMT, offsets[p * (height + 2) + height + 1], closureSize);
5189:   }
5190:   for (h = 0; h < height + 1; ++h) {
5191:     PetscInt dof;

5193:     /* Copy in cone of first point */
5194:     dof = offsets[h + 1] - offsets[h];
5195:     for (meetSize = 0; meetSize < dof; ++meetSize) meet[i][meetSize] = closures[0][(offsets[h] + meetSize) * 2];
5196:     /* Check each successive cone */
5197:     for (p = 1; p < numPoints && meetSize; ++p) {
5198:       PetscInt newMeetSize = 0;

5200:       dof = offsets[p * (height + 2) + h + 1] - offsets[p * (height + 2) + h];
5201:       for (c = 0; c < dof; ++c) {
5202:         const PetscInt point = closures[p][(offsets[p * (height + 2) + h] + c) * 2];

5204:         for (m = 0; m < meetSize; ++m) {
5205:           if (point == meet[i][m]) {
5206:             meet[1 - i][newMeetSize++] = point;
5207:             break;
5208:           }
5209:         }
5210:       }
5211:       meetSize = newMeetSize;
5212:       i        = 1 - i;
5213:     }
5214:     if (meetSize) break;
5215:   }
5216:   *numCoveredPoints = meetSize;
5217:   *coveredPoints    = meet[i];
5218:   for (p = 0; p < numPoints; ++p) PetscCall(DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, NULL, &closures[p]));
5219:   PetscCall(PetscFree(closures));
5220:   PetscCall(DMRestoreWorkArray(dm, numPoints * (height + 2), MPIU_INT, &offsets));
5221:   PetscCall(DMRestoreWorkArray(dm, mc, MPIU_INT, &meet[1 - i]));
5222:   PetscFunctionReturn(PETSC_SUCCESS);
5223: }

5225: /*@
5226:   DMPlexEqual - Determine if two `DM` have the same topology

5228:   Not Collective

5230:   Input Parameters:
5231: + dmA - A `DMPLEX` object
5232: - dmB - A `DMPLEX` object

5234:   Output Parameter:
5235: . equal - `PETSC_TRUE` if the topologies are identical

5237:   Level: intermediate

5239:   Note:
5240:   We are not solving graph isomorphism, so we do not permute.

5242: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCone()`
5243: @*/
5244: PetscErrorCode DMPlexEqual(DM dmA, DM dmB, PetscBool *equal)
5245: {
5246:   PetscInt depth, depthB, pStart, pEnd, pStartB, pEndB, p;

5248:   PetscFunctionBegin;
5251:   PetscAssertPointer(equal, 3);

5253:   *equal = PETSC_FALSE;
5254:   PetscCall(DMPlexGetDepth(dmA, &depth));
5255:   PetscCall(DMPlexGetDepth(dmB, &depthB));
5256:   if (depth != depthB) PetscFunctionReturn(PETSC_SUCCESS);
5257:   PetscCall(DMPlexGetChart(dmA, &pStart, &pEnd));
5258:   PetscCall(DMPlexGetChart(dmB, &pStartB, &pEndB));
5259:   if ((pStart != pStartB) || (pEnd != pEndB)) PetscFunctionReturn(PETSC_SUCCESS);
5260:   for (p = pStart; p < pEnd; ++p) {
5261:     const PetscInt *cone, *coneB, *ornt, *orntB, *support, *supportB;
5262:     PetscInt        coneSize, coneSizeB, c, supportSize, supportSizeB, s;

5264:     PetscCall(DMPlexGetConeSize(dmA, p, &coneSize));
5265:     PetscCall(DMPlexGetCone(dmA, p, &cone));
5266:     PetscCall(DMPlexGetConeOrientation(dmA, p, &ornt));
5267:     PetscCall(DMPlexGetConeSize(dmB, p, &coneSizeB));
5268:     PetscCall(DMPlexGetCone(dmB, p, &coneB));
5269:     PetscCall(DMPlexGetConeOrientation(dmB, p, &orntB));
5270:     if (coneSize != coneSizeB) PetscFunctionReturn(PETSC_SUCCESS);
5271:     for (c = 0; c < coneSize; ++c) {
5272:       if (cone[c] != coneB[c]) PetscFunctionReturn(PETSC_SUCCESS);
5273:       if (ornt[c] != orntB[c]) PetscFunctionReturn(PETSC_SUCCESS);
5274:     }
5275:     PetscCall(DMPlexGetSupportSize(dmA, p, &supportSize));
5276:     PetscCall(DMPlexGetSupport(dmA, p, &support));
5277:     PetscCall(DMPlexGetSupportSize(dmB, p, &supportSizeB));
5278:     PetscCall(DMPlexGetSupport(dmB, p, &supportB));
5279:     if (supportSize != supportSizeB) PetscFunctionReturn(PETSC_SUCCESS);
5280:     for (s = 0; s < supportSize; ++s) {
5281:       if (support[s] != supportB[s]) PetscFunctionReturn(PETSC_SUCCESS);
5282:     }
5283:   }
5284:   *equal = PETSC_TRUE;
5285:   PetscFunctionReturn(PETSC_SUCCESS);
5286: }

5288: /*@
5289:   DMPlexGetNumFaceVertices - Returns the number of vertices on a face

5291:   Not Collective

5293:   Input Parameters:
5294: + dm         - The `DMPLEX`
5295: . cellDim    - The cell dimension
5296: - numCorners - The number of vertices on a cell

5298:   Output Parameter:
5299: . numFaceVertices - The number of vertices on a face

5301:   Level: developer

5303:   Note:
5304:   Of course this can only work for a restricted set of symmetric shapes

5306: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCone()`
5307: @*/
5308: PetscErrorCode DMPlexGetNumFaceVertices(DM dm, PetscInt cellDim, PetscInt numCorners, PetscInt *numFaceVertices)
5309: {
5310:   MPI_Comm comm;

5312:   PetscFunctionBegin;
5313:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5314:   PetscAssertPointer(numFaceVertices, 4);
5315:   switch (cellDim) {
5316:   case 0:
5317:     *numFaceVertices = 0;
5318:     break;
5319:   case 1:
5320:     *numFaceVertices = 1;
5321:     break;
5322:   case 2:
5323:     switch (numCorners) {
5324:     case 3:                 /* triangle */
5325:       *numFaceVertices = 2; /* Edge has 2 vertices */
5326:       break;
5327:     case 4:                 /* quadrilateral */
5328:       *numFaceVertices = 2; /* Edge has 2 vertices */
5329:       break;
5330:     case 6:                 /* quadratic triangle, tri and quad cohesive Lagrange cells */
5331:       *numFaceVertices = 3; /* Edge has 3 vertices */
5332:       break;
5333:     case 9:                 /* quadratic quadrilateral, quadratic quad cohesive Lagrange cells */
5334:       *numFaceVertices = 3; /* Edge has 3 vertices */
5335:       break;
5336:     default:
5337:       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %" PetscInt_FMT " for dimension %" PetscInt_FMT, numCorners, cellDim);
5338:     }
5339:     break;
5340:   case 3:
5341:     switch (numCorners) {
5342:     case 4:                 /* tetradehdron */
5343:       *numFaceVertices = 3; /* Face has 3 vertices */
5344:       break;
5345:     case 6:                 /* tet cohesive cells */
5346:       *numFaceVertices = 4; /* Face has 4 vertices */
5347:       break;
5348:     case 8:                 /* hexahedron */
5349:       *numFaceVertices = 4; /* Face has 4 vertices */
5350:       break;
5351:     case 9:                 /* tet cohesive Lagrange cells */
5352:       *numFaceVertices = 6; /* Face has 6 vertices */
5353:       break;
5354:     case 10:                /* quadratic tetrahedron */
5355:       *numFaceVertices = 6; /* Face has 6 vertices */
5356:       break;
5357:     case 12:                /* hex cohesive Lagrange cells */
5358:       *numFaceVertices = 6; /* Face has 6 vertices */
5359:       break;
5360:     case 18:                /* quadratic tet cohesive Lagrange cells */
5361:       *numFaceVertices = 6; /* Face has 6 vertices */
5362:       break;
5363:     case 27:                /* quadratic hexahedron, quadratic hex cohesive Lagrange cells */
5364:       *numFaceVertices = 9; /* Face has 9 vertices */
5365:       break;
5366:     default:
5367:       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %" PetscInt_FMT " for dimension %" PetscInt_FMT, numCorners, cellDim);
5368:     }
5369:     break;
5370:   default:
5371:     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid cell dimension %" PetscInt_FMT, cellDim);
5372:   }
5373:   PetscFunctionReturn(PETSC_SUCCESS);
5374: }

5376: /*@
5377:   DMPlexGetDepthLabel - Get the `DMLabel` recording the depth of each point

5379:   Not Collective

5381:   Input Parameter:
5382: . dm - The `DMPLEX` object

5384:   Output Parameter:
5385: . depthLabel - The `DMLabel` recording point depth

5387:   Level: developer

5389: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetDepth()`, `DMPlexGetHeightStratum()`, `DMPlexGetDepthStratum()`, `DMPlexGetPointDepth()`,
5390: @*/
5391: PetscErrorCode DMPlexGetDepthLabel(DM dm, DMLabel *depthLabel)
5392: {
5393:   PetscFunctionBegin;
5395:   PetscAssertPointer(depthLabel, 2);
5396:   *depthLabel = dm->depthLabel;
5397:   PetscFunctionReturn(PETSC_SUCCESS);
5398: }

5400: /*@
5401:   DMPlexGetDepth - Get the depth of the DAG representing this mesh

5403:   Not Collective

5405:   Input Parameter:
5406: . dm - The `DMPLEX` object

5408:   Output Parameter:
5409: . depth - The number of strata (breadth first levels) in the DAG

5411:   Level: developer

5413:   Notes:
5414:   This returns maximum of point depths over all points, i.e. maximum value of the label returned by `DMPlexGetDepthLabel()`.

5416:   The point depth is described more in detail in `DMPlexGetDepthStratum()`.

5418:   An empty mesh gives -1.

5420: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetDepthLabel()`, `DMPlexGetDepthStratum()`, `DMPlexGetPointDepth()`, `DMPlexSymmetrize()`
5421: @*/
5422: PetscErrorCode DMPlexGetDepth(DM dm, PetscInt *depth)
5423: {
5424:   DM_Plex *mesh = (DM_Plex *)dm->data;
5425:   DMLabel  label;
5426:   PetscInt d = -1;

5428:   PetscFunctionBegin;
5430:   PetscAssertPointer(depth, 2);
5431:   if (mesh->tr) {
5432:     PetscCall(DMPlexTransformGetDepth(mesh->tr, depth));
5433:   } else {
5434:     PetscCall(DMPlexGetDepthLabel(dm, &label));
5435:     // Allow missing depths
5436:     if (label) PetscCall(DMLabelGetValueBounds(label, NULL, &d));
5437:     *depth = d;
5438:   }
5439:   PetscFunctionReturn(PETSC_SUCCESS);
5440: }

5442: /*@
5443:   DMPlexGetDepthStratum - Get the bounds [`start`, `end`) for all points at a certain depth.

5445:   Not Collective

5447:   Input Parameters:
5448: + dm    - The `DMPLEX` object
5449: - depth - The requested depth

5451:   Output Parameters:
5452: + start - The first point at this `depth`
5453: - end   - One beyond the last point at this `depth`

5455:   Level: developer

5457:   Notes:
5458:   Depth indexing is related to topological dimension.  Depth stratum 0 contains the lowest topological dimension points,
5459:   often "vertices".  If the mesh is "interpolated" (see `DMPlexInterpolate()`), then depth stratum 1 contains the next
5460:   higher dimension, e.g., "edges".

5462: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetHeightStratum()`, `DMPlexGetCellTypeStratum()`, `DMPlexGetDepth()`, `DMPlexGetDepthLabel()`, `DMPlexGetPointDepth()`, `DMPlexSymmetrize()`, `DMPlexInterpolate()`
5463: @*/
5464: PetscErrorCode DMPlexGetDepthStratum(DM dm, PetscInt depth, PetscInt *start, PetscInt *end)
5465: {
5466:   DM_Plex *mesh = (DM_Plex *)dm->data;
5467:   DMLabel  label;
5468:   PetscInt pStart, pEnd;

5470:   PetscFunctionBegin;
5472:   if (start) {
5473:     PetscAssertPointer(start, 3);
5474:     *start = 0;
5475:   }
5476:   if (end) {
5477:     PetscAssertPointer(end, 4);
5478:     *end = 0;
5479:   }
5480:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5481:   if (pStart == pEnd) PetscFunctionReturn(PETSC_SUCCESS);
5482:   if (depth < 0) {
5483:     if (start) *start = pStart;
5484:     if (end) *end = pEnd;
5485:     PetscFunctionReturn(PETSC_SUCCESS);
5486:   }
5487:   if (mesh->tr) {
5488:     PetscCall(DMPlexTransformGetDepthStratum(mesh->tr, depth, start, end));
5489:   } else {
5490:     PetscCall(DMPlexGetDepthLabel(dm, &label));
5491:     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "No label named depth was found");
5492:     PetscCall(DMLabelGetStratumBounds(label, depth, start, end));
5493:   }
5494:   PetscFunctionReturn(PETSC_SUCCESS);
5495: }

5497: /*@
5498:   DMPlexGetHeightStratum - Get the bounds [`start`, `end`) for all points at a certain height.

5500:   Not Collective

5502:   Input Parameters:
5503: + dm     - The `DMPLEX` object
5504: - height - The requested height

5506:   Output Parameters:
5507: + start - The first point at this `height`
5508: - end   - One beyond the last point at this `height`

5510:   Level: developer

5512:   Notes:
5513:   Height indexing is related to topological codimension.  Height stratum 0 contains the highest topological dimension
5514:   points, often called "cells" or "elements".  If the mesh is "interpolated" (see `DMPlexInterpolate()`), then height
5515:   stratum 1 contains the boundary of these "cells", often called "faces" or "facets".

5517: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetCellTypeStratum()`, `DMPlexGetDepth()`, `DMPlexGetPointHeight()`
5518: @*/
5519: PetscErrorCode DMPlexGetHeightStratum(DM dm, PetscInt height, PetscInt *start, PetscInt *end)
5520: {
5521:   DMLabel  label;
5522:   PetscInt depth, pStart, pEnd;

5524:   PetscFunctionBegin;
5526:   if (start) {
5527:     PetscAssertPointer(start, 3);
5528:     *start = 0;
5529:   }
5530:   if (end) {
5531:     PetscAssertPointer(end, 4);
5532:     *end = 0;
5533:   }
5534:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5535:   if (pStart == pEnd) PetscFunctionReturn(PETSC_SUCCESS);
5536:   if (height < 0) {
5537:     if (start) *start = pStart;
5538:     if (end) *end = pEnd;
5539:     PetscFunctionReturn(PETSC_SUCCESS);
5540:   }
5541:   PetscCall(DMPlexGetDepthLabel(dm, &label));
5542:   if (label) PetscCall(DMLabelGetNumValues(label, &depth));
5543:   else PetscCall(DMGetDimension(dm, &depth));
5544:   PetscCheck(depth >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Depth not yet computed");
5545:   PetscCall(DMPlexGetDepthStratum(dm, depth - 1 - height, start, end));
5546:   PetscFunctionReturn(PETSC_SUCCESS);
5547: }

5549: /*@
5550:   DMPlexGetPointDepth - Get the `depth` of a given point

5552:   Not Collective

5554:   Input Parameters:
5555: + dm    - The `DMPLEX` object
5556: - point - The point

5558:   Output Parameter:
5559: . depth - The depth of the `point`

5561:   Level: intermediate

5563: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellType()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`, `DMPlexGetPointHeight()`
5564: @*/
5565: PetscErrorCode DMPlexGetPointDepth(DM dm, PetscInt point, PetscInt *depth)
5566: {
5567:   PetscFunctionBegin;
5569:   PetscAssertPointer(depth, 3);
5570:   PetscCall(DMLabelGetValue(dm->depthLabel, point, depth));
5571:   PetscFunctionReturn(PETSC_SUCCESS);
5572: }

5574: /*@
5575:   DMPlexGetPointHeight - Get the `height` of a given point

5577:   Not Collective

5579:   Input Parameters:
5580: + dm    - The `DMPLEX` object
5581: - point - The point

5583:   Output Parameter:
5584: . height - The height of the `point`

5586:   Level: intermediate

5588: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellType()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`, `DMPlexGetPointDepth()`
5589: @*/
5590: PetscErrorCode DMPlexGetPointHeight(DM dm, PetscInt point, PetscInt *height)
5591: {
5592:   PetscInt n, pDepth;

5594:   PetscFunctionBegin;
5596:   PetscAssertPointer(height, 3);
5597:   PetscCall(DMLabelGetNumValues(dm->depthLabel, &n));
5598:   PetscCall(DMLabelGetValue(dm->depthLabel, point, &pDepth));
5599:   *height = n - 1 - pDepth; /* DAG depth is n-1 */
5600:   PetscFunctionReturn(PETSC_SUCCESS);
5601: }

5603: /*@
5604:   DMPlexGetCellTypeLabel - Get the `DMLabel` recording the polytope type of each cell

5606:   Not Collective

5608:   Input Parameter:
5609: . dm - The `DMPLEX` object

5611:   Output Parameter:
5612: . celltypeLabel - The `DMLabel` recording cell polytope type

5614:   Level: developer

5616:   Note:
5617:   This function will trigger automatica computation of cell types. This can be disabled by calling
5618:   `DMCreateLabel`(dm, "celltype") beforehand.

5620: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellType()`, `DMPlexGetDepthLabel()`, `DMCreateLabel()`
5621: @*/
5622: PetscErrorCode DMPlexGetCellTypeLabel(DM dm, DMLabel *celltypeLabel)
5623: {
5624:   PetscFunctionBegin;
5626:   PetscAssertPointer(celltypeLabel, 2);
5627:   if (!dm->celltypeLabel) PetscCall(DMPlexComputeCellTypes(dm));
5628:   *celltypeLabel = dm->celltypeLabel;
5629:   PetscFunctionReturn(PETSC_SUCCESS);
5630: }

5632: /*@
5633:   DMPlexGetCellType - Get the polytope type of a given cell

5635:   Not Collective

5637:   Input Parameters:
5638: + dm   - The `DMPLEX` object
5639: - cell - The cell

5641:   Output Parameter:
5642: . celltype - The polytope type of the cell

5644:   Level: intermediate

5646: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPolytopeType`, `DMPlexGetCellTypeLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`
5647: @*/
5648: PetscErrorCode DMPlexGetCellType(DM dm, PetscInt cell, DMPolytopeType *celltype)
5649: {
5650:   DM_Plex *mesh = (DM_Plex *)dm->data;
5651:   DMLabel  label;
5652:   PetscInt ct;

5654:   PetscFunctionBegin;
5656:   PetscAssertPointer(celltype, 3);
5657:   if (mesh->tr) {
5658:     PetscCall(DMPlexTransformGetCellType(mesh->tr, cell, celltype));
5659:   } else {
5660:     PetscInt pStart, pEnd;

5662:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, NULL));
5663:     if (!mesh->cellTypes) { /* XXX remove? optimize? */
5664:       PetscCall(PetscSectionGetChart(mesh->coneSection, NULL, &pEnd));
5665:       PetscCall(PetscMalloc1(pEnd - pStart, &mesh->cellTypes));
5666:       PetscCall(DMPlexGetCellTypeLabel(dm, &label));
5667:       for (PetscInt p = pStart; p < pEnd; p++) {
5668:         PetscCall(DMLabelGetValue(label, p, &ct));
5669:         mesh->cellTypes[p - pStart].value_as_uint8 = (uint8_t)ct;
5670:       }
5671:     }
5672:     *celltype = (DMPolytopeType)mesh->cellTypes[cell - pStart].value_as_uint8;
5673:     if (PetscDefined(USE_DEBUG)) {
5674:       PetscCall(DMPlexGetCellTypeLabel(dm, &label));
5675:       PetscCall(DMLabelGetValue(label, cell, &ct));
5676:       PetscCheck(ct >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " has not been assigned a cell type", cell);
5677:       PetscCheck(ct == (PetscInt)*celltype, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid cellType for %" PetscInt_FMT ": %d != %" PetscInt_FMT, cell, (int)*celltype, ct);
5678:     }
5679:   }
5680:   PetscFunctionReturn(PETSC_SUCCESS);
5681: }

5683: /*@
5684:   DMPlexSetCellType - Set the polytope type of a given cell

5686:   Not Collective

5688:   Input Parameters:
5689: + dm       - The `DMPLEX` object
5690: . cell     - The cell
5691: - celltype - The polytope type of the cell

5693:   Level: advanced

5695:   Note:
5696:   By default, cell types will be automatically computed using `DMPlexComputeCellTypes()` before this function
5697:   is executed. This function will override the computed type. However, if automatic classification will not succeed
5698:   and a user wants to manually specify all types, the classification must be disabled by calling
5699:   DMCreateLabel(dm, "celltype") before getting or setting any cell types.

5701: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellTypeLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`, `DMPlexComputeCellTypes()`, `DMCreateLabel()`
5702: @*/
5703: PetscErrorCode DMPlexSetCellType(DM dm, PetscInt cell, DMPolytopeType celltype)
5704: {
5705:   DM_Plex *mesh = (DM_Plex *)dm->data;
5706:   DMLabel  label;
5707:   PetscInt pStart, pEnd;

5709:   PetscFunctionBegin;
5711:   PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
5712:   PetscCall(DMPlexGetCellTypeLabel(dm, &label));
5713:   PetscCall(DMLabelSetValue(label, cell, celltype));
5714:   if (!mesh->cellTypes) PetscCall(PetscMalloc1(pEnd - pStart, &mesh->cellTypes));
5715:   mesh->cellTypes[cell - pStart].value_as_uint8 = (uint8_t)celltype;
5716:   PetscFunctionReturn(PETSC_SUCCESS);
5717: }

5719: PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm)
5720: {
5721:   PetscSection section;
5722:   PetscInt     maxHeight;
5723:   const char  *prefix;

5725:   PetscFunctionBegin;
5726:   PetscCall(DMClone(dm, cdm));
5727:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
5728:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix));
5729:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_"));
5730:   PetscCall(DMPlexGetMaxProjectionHeight(dm, &maxHeight));
5731:   PetscCall(DMPlexSetMaxProjectionHeight(*cdm, maxHeight));
5732:   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
5733:   PetscCall(DMSetLocalSection(*cdm, section));
5734:   PetscCall(PetscSectionDestroy(&section));

5736:   PetscCall(DMSetNumFields(*cdm, 1));
5737:   PetscCall(DMCreateDS(*cdm));
5738:   (*cdm)->cloneOpts = PETSC_TRUE;
5739:   if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(*cdm));
5740:   PetscFunctionReturn(PETSC_SUCCESS);
5741: }

5743: PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field)
5744: {
5745:   Vec coordsLocal, cellCoordsLocal;
5746:   DM  coordsDM, cellCoordsDM;

5748:   PetscFunctionBegin;
5749:   *field = NULL;
5750:   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
5751:   PetscCall(DMGetCoordinateDM(dm, &coordsDM));
5752:   PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoordsLocal));
5753:   PetscCall(DMGetCellCoordinateDM(dm, &cellCoordsDM));
5754:   if (coordsLocal && coordsDM) {
5755:     if (cellCoordsLocal && cellCoordsDM) PetscCall(DMFieldCreateDSWithDG(coordsDM, cellCoordsDM, 0, coordsLocal, cellCoordsLocal, field));
5756:     else PetscCall(DMFieldCreateDS(coordsDM, 0, coordsLocal, field));
5757:   }
5758:   PetscFunctionReturn(PETSC_SUCCESS);
5759: }

5761: /*@
5762:   DMPlexGetConeSection - Return a section which describes the layout of cone data

5764:   Not Collective

5766:   Input Parameter:
5767: . dm - The `DMPLEX` object

5769:   Output Parameter:
5770: . section - The `PetscSection` object

5772:   Level: developer

5774: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetSupportSection()`, `DMPlexGetCones()`, `DMPlexGetConeOrientations()`, `PetscSection`
5775: @*/
5776: PetscErrorCode DMPlexGetConeSection(DM dm, PetscSection *section)
5777: {
5778:   DM_Plex *mesh = (DM_Plex *)dm->data;

5780:   PetscFunctionBegin;
5782:   if (section) *section = mesh->coneSection;
5783:   PetscFunctionReturn(PETSC_SUCCESS);
5784: }

5786: /*@
5787:   DMPlexGetSupportSection - Return a section which describes the layout of support data

5789:   Not Collective

5791:   Input Parameter:
5792: . dm - The `DMPLEX` object

5794:   Output Parameter:
5795: . section - The `PetscSection` object

5797:   Level: developer

5799: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSection()`, `PetscSection`
5800: @*/
5801: PetscErrorCode DMPlexGetSupportSection(DM dm, PetscSection *section)
5802: {
5803:   DM_Plex *mesh = (DM_Plex *)dm->data;

5805:   PetscFunctionBegin;
5807:   if (section) *section = mesh->supportSection;
5808:   PetscFunctionReturn(PETSC_SUCCESS);
5809: }

5811: /*@C
5812:   DMPlexGetCones - Return cone data

5814:   Not Collective

5816:   Input Parameter:
5817: . dm - The `DMPLEX` object

5819:   Output Parameter:
5820: . cones - The cone for each point

5822:   Level: developer

5824: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSection()`
5825: @*/
5826: PetscErrorCode DMPlexGetCones(DM dm, PetscInt *cones[])
5827: {
5828:   DM_Plex *mesh = (DM_Plex *)dm->data;

5830:   PetscFunctionBegin;
5832:   if (cones) *cones = mesh->cones;
5833:   PetscFunctionReturn(PETSC_SUCCESS);
5834: }

5836: /*@C
5837:   DMPlexGetConeOrientations - Return cone orientation data

5839:   Not Collective

5841:   Input Parameter:
5842: . dm - The `DMPLEX` object

5844:   Output Parameter:
5845: . coneOrientations - The array of cone orientations for all points

5847:   Level: developer

5849:   Notes:
5850:   The `PetscSection` returned by `DMPlexGetConeSection()` partitions coneOrientations into cone orientations of particular points
5851:   as returned by `DMPlexGetConeOrientation()`.

5853:   The meaning of coneOrientations values is detailed in `DMPlexGetConeOrientation()`.

5855: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSection()`, `DMPlexGetConeOrientation()`, `PetscSection`
5856: @*/
5857: PetscErrorCode DMPlexGetConeOrientations(DM dm, PetscInt *coneOrientations[])
5858: {
5859:   DM_Plex *mesh = (DM_Plex *)dm->data;

5861:   PetscFunctionBegin;
5863:   if (coneOrientations) *coneOrientations = mesh->coneOrientations;
5864:   PetscFunctionReturn(PETSC_SUCCESS);
5865: }

5867: /* FEM Support */

5869: PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
5870: {
5871:   PetscInt depth;

5873:   PetscFunctionBegin;
5874:   PetscCall(DMPlexGetDepth(plex, &depth));
5875:   PetscCall(DMGetStratumIS(plex, "dim", depth, cellIS));
5876:   if (!*cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, cellIS));
5877:   PetscFunctionReturn(PETSC_SUCCESS);
5878: }

5880: PetscErrorCode DMPlexGetAllFaces_Internal(DM plex, IS *faceIS)
5881: {
5882:   PetscInt depth;

5884:   PetscFunctionBegin;
5885:   PetscCall(DMPlexGetDepth(plex, &depth));
5886:   PetscCall(DMGetStratumIS(plex, "dim", depth - 1, faceIS));
5887:   if (!*faceIS) PetscCall(DMGetStratumIS(plex, "depth", depth - 1, faceIS));
5888:   PetscFunctionReturn(PETSC_SUCCESS);
5889: }

5891: /*
5892:  Returns number of components and tensor degree for the field.  For interpolated meshes, line should be a point
5893:  representing a line in the section.
5894: */
5895: static PetscErrorCode PetscSectionFieldGetTensorDegree_Private(DM dm, PetscSection section, PetscInt field, PetscInt line, PetscInt *Nc, PetscInt *k, PetscBool *continuous, PetscBool *tensor)
5896: {
5897:   PetscObject  obj;
5898:   PetscClassId id;
5899:   PetscFE      fe = NULL;

5901:   PetscFunctionBeginHot;
5902:   PetscCall(PetscSectionGetFieldComponents(section, field, Nc));
5903:   PetscCall(DMGetField(dm, field, NULL, &obj));
5904:   PetscCall(PetscObjectGetClassId(obj, &id));
5905:   if (id == PETSCFE_CLASSID) fe = (PetscFE)obj;

5907:   if (!fe) {
5908:     /* Assume the full interpolated mesh is in the chart; lines in particular */
5909:     /* An order k SEM disc has k-1 dofs on an edge */
5910:     PetscCall(PetscSectionGetFieldDof(section, line, field, k));
5911:     *k = *k / *Nc + 1;
5912:   } else {
5913:     PetscInt       dual_space_size, dim;
5914:     PetscDualSpace dsp;

5916:     PetscCall(DMGetDimension(dm, &dim));
5917:     PetscCall(PetscFEGetDualSpace(fe, &dsp));
5918:     PetscCall(PetscDualSpaceGetDimension(dsp, &dual_space_size));
5919:     *k = (PetscInt)PetscCeilReal(PetscPowReal(dual_space_size / *Nc, 1.0 / dim)) - 1;
5920:     PetscCall(PetscDualSpaceLagrangeGetContinuity(dsp, continuous));
5921:     PetscCall(PetscDualSpaceLagrangeGetTensor(dsp, tensor));
5922:   }
5923:   PetscFunctionReturn(PETSC_SUCCESS);
5924: }

5926: static PetscErrorCode GetFieldSize_Private(PetscInt dim, PetscInt k, PetscBool tensor, PetscInt *dof)
5927: {
5928:   PetscFunctionBeginHot;
5929:   if (tensor) {
5930:     *dof = PetscPowInt(k + 1, dim);
5931:   } else {
5932:     switch (dim) {
5933:     case 1:
5934:       *dof = k + 1;
5935:       break;
5936:     case 2:
5937:       *dof = ((k + 1) * (k + 2)) / 2;
5938:       break;
5939:     case 3:
5940:       *dof = ((k + 1) * (k + 2) * (k + 3)) / 6;
5941:       break;
5942:     default:
5943:       *dof = 0;
5944:     }
5945:   }
5946:   PetscFunctionReturn(PETSC_SUCCESS);
5947: }

5949: /*@
5950:   DMPlexSetClosurePermutationTensor - Create a permutation from the default (BFS) point ordering in the closure, to a
5951:   lexicographic ordering over the tensor product cell (i.e., line, quad, hex, etc.), and set this permutation in the
5952:   section provided (or the section of the `DM`).

5954:   Input Parameters:
5955: + dm      - The `DM`
5956: . point   - Either a cell (highest dim point) or an edge (dim 1 point), or `PETSC_DETERMINE`
5957: - section - The `PetscSection` to reorder, or `NULL` for the default section

5959:   Example:
5960:   A typical interpolated single-quad mesh might order points as
5961: .vb
5962:   [c0, v1, v2, v3, v4, e5, e6, e7, e8]

5964:   v4 -- e6 -- v3
5965:   |           |
5966:   e7    c0    e8
5967:   |           |
5968:   v1 -- e5 -- v2
5969: .ve

5971:   (There is no significance to the ordering described here.)  The default section for a Q3 quad might typically assign
5972:   dofs in the order of points, e.g.,
5973: .vb
5974:     c0 -> [0,1,2,3]
5975:     v1 -> [4]
5976:     ...
5977:     e5 -> [8, 9]
5978: .ve

5980:   which corresponds to the dofs
5981: .vb
5982:     6   10  11  7
5983:     13  2   3   15
5984:     12  0   1   14
5985:     4   8   9   5
5986: .ve

5988:   The closure in BFS ordering works through height strata (cells, edges, vertices) to produce the ordering
5989: .vb
5990:   0 1 2 3 8 9 14 15 11 10 13 12 4 5 7 6
5991: .ve

5993:   After calling DMPlexSetClosurePermutationTensor(), the closure will be ordered lexicographically,
5994: .vb
5995:    4 8 9 5 12 0 1 14 13 2 3 15 6 10 11 7
5996: .ve

5998:   Level: developer

6000:   Notes:
6001:   The point is used to determine the number of dofs/field on an edge. For SEM, this is related to the polynomial
6002:   degree of the basis.

6004:   This is required to run with libCEED.

6006: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMGetLocalSection()`, `PetscSectionSetClosurePermutation()`, `DMSetGlobalSection()`
6007: @*/
6008: PetscErrorCode DMPlexSetClosurePermutationTensor(DM dm, PetscInt point, PetscSection section)
6009: {
6010:   DMLabel   label;
6011:   PetscInt  dim, depth = -1, eStart = -1, Nf;
6012:   PetscBool continuous = PETSC_TRUE, tensor = PETSC_TRUE;

6014:   PetscFunctionBegin;
6015:   PetscCall(DMGetDimension(dm, &dim));
6016:   if (dim < 1) PetscFunctionReturn(PETSC_SUCCESS);
6017:   if (point < 0) {
6018:     PetscInt sStart, sEnd;

6020:     PetscCall(DMPlexGetDepthStratum(dm, 1, &sStart, &sEnd));
6021:     point = sEnd - sStart ? sStart : point;
6022:   }
6023:   PetscCall(DMPlexGetDepthLabel(dm, &label));
6024:   if (point >= 0) PetscCall(DMLabelGetValue(label, point, &depth));
6025:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
6026:   if (depth == 1) {
6027:     eStart = point;
6028:   } else if (depth == dim) {
6029:     const PetscInt *cone;

6031:     PetscCall(DMPlexGetCone(dm, point, &cone));
6032:     if (dim == 2) eStart = cone[0];
6033:     else if (dim == 3) {
6034:       const PetscInt *cone2;
6035:       PetscCall(DMPlexGetCone(dm, cone[0], &cone2));
6036:       eStart = cone2[0];
6037:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " of depth %" PetscInt_FMT " cannot be used to bootstrap spectral ordering for dim %" PetscInt_FMT, point, depth, dim);
6038:   } else PetscCheck(depth < 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " of depth %" PetscInt_FMT " cannot be used to bootstrap spectral ordering for dim %" PetscInt_FMT, point, depth, dim);

6040:   PetscCall(PetscSectionGetNumFields(section, &Nf));
6041:   for (PetscInt d = 1; d <= dim; d++) {
6042:     PetscInt  k, f, Nc, c, i, j, size = 0, offset = 0, foffset = 0;
6043:     PetscInt *perm;

6045:     for (f = 0; f < Nf; ++f) {
6046:       PetscInt dof;

6048:       PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
6049:       PetscCheck(dim == 1 || tensor || !continuous, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Continuous field %" PetscInt_FMT " must have a tensor product discretization", f);
6050:       if (!continuous && d < dim) continue;
6051:       PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
6052:       size += dof * Nc;
6053:     }
6054:     PetscCall(PetscMalloc1(size, &perm));
6055:     for (f = 0; f < Nf; ++f) {
6056:       switch (d) {
6057:       case 1:
6058:         PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
6059:         if (!continuous && d < dim) continue;
6060:         /*
6061:          Original ordering is [ edge of length k-1; vtx0; vtx1 ]
6062:          We want              [ vtx0; edge of length k-1; vtx1 ]
6063:          */
6064:         if (continuous) {
6065:           for (c = 0; c < Nc; c++, offset++) perm[offset] = (k - 1) * Nc + c + foffset;
6066:           for (i = 0; i < k - 1; i++)
6067:             for (c = 0; c < Nc; c++, offset++) perm[offset] = i * Nc + c + foffset;
6068:           for (c = 0; c < Nc; c++, offset++) perm[offset] = k * Nc + c + foffset;
6069:           foffset = offset;
6070:         } else {
6071:           PetscInt dof;

6073:           PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
6074:           for (i = 0; i < dof * Nc; ++i, ++offset) perm[offset] = i + foffset;
6075:           foffset = offset;
6076:         }
6077:         break;
6078:       case 2:
6079:         /* The original quad closure is oriented clockwise, {f, e_b, e_r, e_t, e_l, v_lb, v_rb, v_tr, v_tl} */
6080:         PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
6081:         if (!continuous && d < dim) continue;
6082:         /* The SEM order is

6084:          v_lb, {e_b}, v_rb,
6085:          e^{(k-1)-i}_l, {f^{i*(k-1)}}, e^i_r,
6086:          v_lt, reverse {e_t}, v_rt
6087:          */
6088:         if (continuous) {
6089:           const PetscInt of   = 0;
6090:           const PetscInt oeb  = of + PetscSqr(k - 1);
6091:           const PetscInt oer  = oeb + (k - 1);
6092:           const PetscInt oet  = oer + (k - 1);
6093:           const PetscInt oel  = oet + (k - 1);
6094:           const PetscInt ovlb = oel + (k - 1);
6095:           const PetscInt ovrb = ovlb + 1;
6096:           const PetscInt ovrt = ovrb + 1;
6097:           const PetscInt ovlt = ovrt + 1;
6098:           PetscInt       o;

6100:           /* bottom */
6101:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlb * Nc + c + foffset;
6102:           for (o = oeb; o < oer; ++o)
6103:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6104:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrb * Nc + c + foffset;
6105:           /* middle */
6106:           for (i = 0; i < k - 1; ++i) {
6107:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oel + (k - 2) - i) * Nc + c + foffset;
6108:             for (o = of + (k - 1) * i; o < of + (k - 1) * (i + 1); ++o)
6109:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6110:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oer + i) * Nc + c + foffset;
6111:           }
6112:           /* top */
6113:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlt * Nc + c + foffset;
6114:           for (o = oel - 1; o >= oet; --o)
6115:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6116:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrt * Nc + c + foffset;
6117:           foffset = offset;
6118:         } else {
6119:           PetscInt dof;

6121:           PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
6122:           for (i = 0; i < dof * Nc; ++i, ++offset) perm[offset] = i + foffset;
6123:           foffset = offset;
6124:         }
6125:         break;
6126:       case 3:
6127:         /* The original hex closure is

6129:          {c,
6130:          f_b, f_t, f_f, f_b, f_r, f_l,
6131:          e_bl, e_bb, e_br, e_bf,  e_tf, e_tr, e_tb, e_tl,  e_rf, e_lf, e_lb, e_rb,
6132:          v_blf, v_blb, v_brb, v_brf, v_tlf, v_trf, v_trb, v_tlb}
6133:          */
6134:         PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
6135:         if (!continuous && d < dim) continue;
6136:         /* The SEM order is
6137:          Bottom Slice
6138:          v_blf, {e^{(k-1)-n}_bf}, v_brf,
6139:          e^{i}_bl, f^{n*(k-1)+(k-1)-i}_b, e^{(k-1)-i}_br,
6140:          v_blb, {e_bb}, v_brb,

6142:          Middle Slice (j)
6143:          {e^{(k-1)-j}_lf}, {f^{j*(k-1)+n}_f}, e^j_rf,
6144:          f^{i*(k-1)+j}_l, {c^{(j*(k-1) + i)*(k-1)+n}_t}, f^{j*(k-1)+i}_r,
6145:          e^j_lb, {f^{j*(k-1)+(k-1)-n}_b}, e^{(k-1)-j}_rb,

6147:          Top Slice
6148:          v_tlf, {e_tf}, v_trf,
6149:          e^{(k-1)-i}_tl, {f^{i*(k-1)}_t}, e^{i}_tr,
6150:          v_tlb, {e^{(k-1)-n}_tb}, v_trb,
6151:          */
6152:         if (continuous) {
6153:           const PetscInt oc    = 0;
6154:           const PetscInt ofb   = oc + PetscSqr(k - 1) * (k - 1);
6155:           const PetscInt oft   = ofb + PetscSqr(k - 1);
6156:           const PetscInt off   = oft + PetscSqr(k - 1);
6157:           const PetscInt ofk   = off + PetscSqr(k - 1);
6158:           const PetscInt ofr   = ofk + PetscSqr(k - 1);
6159:           const PetscInt ofl   = ofr + PetscSqr(k - 1);
6160:           const PetscInt oebl  = ofl + PetscSqr(k - 1);
6161:           const PetscInt oebb  = oebl + (k - 1);
6162:           const PetscInt oebr  = oebb + (k - 1);
6163:           const PetscInt oebf  = oebr + (k - 1);
6164:           const PetscInt oetf  = oebf + (k - 1);
6165:           const PetscInt oetr  = oetf + (k - 1);
6166:           const PetscInt oetb  = oetr + (k - 1);
6167:           const PetscInt oetl  = oetb + (k - 1);
6168:           const PetscInt oerf  = oetl + (k - 1);
6169:           const PetscInt oelf  = oerf + (k - 1);
6170:           const PetscInt oelb  = oelf + (k - 1);
6171:           const PetscInt oerb  = oelb + (k - 1);
6172:           const PetscInt ovblf = oerb + (k - 1);
6173:           const PetscInt ovblb = ovblf + 1;
6174:           const PetscInt ovbrb = ovblb + 1;
6175:           const PetscInt ovbrf = ovbrb + 1;
6176:           const PetscInt ovtlf = ovbrf + 1;
6177:           const PetscInt ovtrf = ovtlf + 1;
6178:           const PetscInt ovtrb = ovtrf + 1;
6179:           const PetscInt ovtlb = ovtrb + 1;
6180:           PetscInt       o, n;

6182:           /* Bottom Slice */
6183:           /*   bottom */
6184:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblf * Nc + c + foffset;
6185:           for (o = oetf - 1; o >= oebf; --o)
6186:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6187:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrf * Nc + c + foffset;
6188:           /*   middle */
6189:           for (i = 0; i < k - 1; ++i) {
6190:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebl + i) * Nc + c + foffset;
6191:             for (n = 0; n < k - 1; ++n) {
6192:               o = ofb + n * (k - 1) + i;
6193:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6194:             }
6195:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebr + (k - 2) - i) * Nc + c + foffset;
6196:           }
6197:           /*   top */
6198:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblb * Nc + c + foffset;
6199:           for (o = oebb; o < oebr; ++o)
6200:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6201:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrb * Nc + c + foffset;

6203:           /* Middle Slice */
6204:           for (j = 0; j < k - 1; ++j) {
6205:             /*   bottom */
6206:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelf + (k - 2) - j) * Nc + c + foffset;
6207:             for (o = off + j * (k - 1); o < off + (j + 1) * (k - 1); ++o)
6208:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6209:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerf + j) * Nc + c + foffset;
6210:             /*   middle */
6211:             for (i = 0; i < k - 1; ++i) {
6212:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofl + i * (k - 1) + j) * Nc + c + foffset;
6213:               for (n = 0; n < k - 1; ++n)
6214:                 for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oc + (j * (k - 1) + i) * (k - 1) + n) * Nc + c + foffset;
6215:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofr + j * (k - 1) + i) * Nc + c + foffset;
6216:             }
6217:             /*   top */
6218:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelb + j) * Nc + c + foffset;
6219:             for (o = ofk + j * (k - 1) + (k - 2); o >= ofk + j * (k - 1); --o)
6220:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6221:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerb + (k - 2) - j) * Nc + c + foffset;
6222:           }

6224:           /* Top Slice */
6225:           /*   bottom */
6226:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlf * Nc + c + foffset;
6227:           for (o = oetf; o < oetr; ++o)
6228:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6229:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrf * Nc + c + foffset;
6230:           /*   middle */
6231:           for (i = 0; i < k - 1; ++i) {
6232:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetl + (k - 2) - i) * Nc + c + foffset;
6233:             for (n = 0; n < k - 1; ++n)
6234:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oft + i * (k - 1) + n) * Nc + c + foffset;
6235:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetr + i) * Nc + c + foffset;
6236:           }
6237:           /*   top */
6238:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlb * Nc + c + foffset;
6239:           for (o = oetl - 1; o >= oetb; --o)
6240:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6241:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrb * Nc + c + foffset;

6243:           foffset = offset;
6244:         } else {
6245:           PetscInt dof;

6247:           PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
6248:           for (i = 0; i < dof * Nc; ++i, ++offset) perm[offset] = i + foffset;
6249:           foffset = offset;
6250:         }
6251:         break;
6252:       default:
6253:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No spectral ordering for dimension %" PetscInt_FMT, d);
6254:       }
6255:     }
6256:     PetscCheck(offset == size, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Number of permutation entries %" PetscInt_FMT " != %" PetscInt_FMT, offset, size);
6257:     /* Check permutation */
6258:     {
6259:       PetscInt *check;

6261:       PetscCall(PetscMalloc1(size, &check));
6262:       for (i = 0; i < size; ++i) {
6263:         check[i] = -1;
6264:         PetscCheck(perm[i] >= 0 && perm[i] < size, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid permutation index p[%" PetscInt_FMT "] = %" PetscInt_FMT, i, perm[i]);
6265:       }
6266:       for (i = 0; i < size; ++i) check[perm[i]] = i;
6267:       for (i = 0; i < size; ++i) PetscCheck(check[i] >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing permutation index %" PetscInt_FMT, i);
6268:       PetscCall(PetscFree(check));
6269:     }
6270:     PetscCall(PetscSectionSetClosurePermutation_Internal(section, (PetscObject)dm, d, size, PETSC_OWN_POINTER, perm));
6271:     if (d == dim) { // Add permutation for localized (in case this is a coordinate DM)
6272:       PetscInt *loc_perm;
6273:       PetscCall(PetscMalloc1(size * 2, &loc_perm));
6274:       for (PetscInt i = 0; i < size; i++) {
6275:         loc_perm[i]        = perm[i];
6276:         loc_perm[size + i] = size + perm[i];
6277:       }
6278:       PetscCall(PetscSectionSetClosurePermutation_Internal(section, (PetscObject)dm, d, size * 2, PETSC_OWN_POINTER, loc_perm));
6279:     }
6280:   }
6281:   PetscFunctionReturn(PETSC_SUCCESS);
6282: }

6284: PetscErrorCode DMPlexGetPointDualSpaceFEM(DM dm, PetscInt point, PetscInt field, PetscDualSpace *dspace)
6285: {
6286:   PetscDS  prob;
6287:   PetscInt depth, Nf, h;
6288:   DMLabel  label;

6290:   PetscFunctionBeginHot;
6291:   PetscCall(DMGetDS(dm, &prob));
6292:   Nf      = prob->Nf;
6293:   label   = dm->depthLabel;
6294:   *dspace = NULL;
6295:   if (field < Nf) {
6296:     PetscObject disc = prob->disc[field];

6298:     if (disc->classid == PETSCFE_CLASSID) {
6299:       PetscDualSpace dsp;

6301:       PetscCall(PetscFEGetDualSpace((PetscFE)disc, &dsp));
6302:       PetscCall(DMLabelGetNumValues(label, &depth));
6303:       PetscCall(DMLabelGetValue(label, point, &h));
6304:       h = depth - 1 - h;
6305:       if (h) {
6306:         PetscCall(PetscDualSpaceGetHeightSubspace(dsp, h, dspace));
6307:       } else {
6308:         *dspace = dsp;
6309:       }
6310:     }
6311:   }
6312:   PetscFunctionReturn(PETSC_SUCCESS);
6313: }

6315: static inline PetscErrorCode DMPlexVecGetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
6316: {
6317:   PetscScalar       *array;
6318:   const PetscScalar *vArray;
6319:   const PetscInt    *cone, *coneO;
6320:   PetscInt           pStart, pEnd, p, numPoints, size = 0, offset = 0;

6322:   PetscFunctionBeginHot;
6323:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6324:   PetscCall(DMPlexGetConeSize(dm, point, &numPoints));
6325:   PetscCall(DMPlexGetCone(dm, point, &cone));
6326:   PetscCall(DMPlexGetConeOrientation(dm, point, &coneO));
6327:   if (!values || !*values) {
6328:     if ((point >= pStart) && (point < pEnd)) {
6329:       PetscInt dof;

6331:       PetscCall(PetscSectionGetDof(section, point, &dof));
6332:       size += dof;
6333:     }
6334:     for (p = 0; p < numPoints; ++p) {
6335:       const PetscInt cp = cone[p];
6336:       PetscInt       dof;

6338:       if ((cp < pStart) || (cp >= pEnd)) continue;
6339:       PetscCall(PetscSectionGetDof(section, cp, &dof));
6340:       size += dof;
6341:     }
6342:     if (!values) {
6343:       if (csize) *csize = size;
6344:       PetscFunctionReturn(PETSC_SUCCESS);
6345:     }
6346:     PetscCall(DMGetWorkArray(dm, size, MPIU_SCALAR, &array));
6347:   } else {
6348:     array = *values;
6349:   }
6350:   size = 0;
6351:   PetscCall(VecGetArrayRead(v, &vArray));
6352:   if ((point >= pStart) && (point < pEnd)) {
6353:     PetscInt           dof, off, d;
6354:     const PetscScalar *varr;

6356:     PetscCall(PetscSectionGetDof(section, point, &dof));
6357:     PetscCall(PetscSectionGetOffset(section, point, &off));
6358:     varr = PetscSafePointerPlusOffset(vArray, off);
6359:     for (d = 0; d < dof; ++d, ++offset) array[offset] = varr[d];
6360:     size += dof;
6361:   }
6362:   for (p = 0; p < numPoints; ++p) {
6363:     const PetscInt     cp = cone[p];
6364:     PetscInt           o  = coneO[p];
6365:     PetscInt           dof, off, d;
6366:     const PetscScalar *varr;

6368:     if ((cp < pStart) || (cp >= pEnd)) continue;
6369:     PetscCall(PetscSectionGetDof(section, cp, &dof));
6370:     PetscCall(PetscSectionGetOffset(section, cp, &off));
6371:     varr = PetscSafePointerPlusOffset(vArray, off);
6372:     if (o >= 0) {
6373:       for (d = 0; d < dof; ++d, ++offset) array[offset] = varr[d];
6374:     } else {
6375:       for (d = dof - 1; d >= 0; --d, ++offset) array[offset] = varr[d];
6376:     }
6377:     size += dof;
6378:   }
6379:   PetscCall(VecRestoreArrayRead(v, &vArray));
6380:   if (!*values) {
6381:     if (csize) *csize = size;
6382:     *values = array;
6383:   } else {
6384:     PetscCheck(size <= *csize, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %" PetscInt_FMT " < actual size %" PetscInt_FMT, *csize, size);
6385:     *csize = size;
6386:   }
6387:   PetscFunctionReturn(PETSC_SUCCESS);
6388: }

6390: /* Compress out points not in the section */
6391: static inline PetscErrorCode CompressPoints_Private(PetscSection section, PetscInt *numPoints, PetscInt points[])
6392: {
6393:   const PetscInt np = *numPoints;
6394:   PetscInt       pStart, pEnd, p, q;

6396:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6397:   for (p = 0, q = 0; p < np; ++p) {
6398:     const PetscInt r = points[p * 2];
6399:     if ((r >= pStart) && (r < pEnd)) {
6400:       points[q * 2]     = r;
6401:       points[q * 2 + 1] = points[p * 2 + 1];
6402:       ++q;
6403:     }
6404:   }
6405:   *numPoints = q;
6406:   return PETSC_SUCCESS;
6407: }

6409: /* Compressed closure does not apply closure permutation */
6410: PetscErrorCode DMPlexGetCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt ornt, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
6411: {
6412:   const PetscInt *cla = NULL;
6413:   PetscInt        np, *pts = NULL;

6415:   PetscFunctionBeginHot;
6416:   PetscCall(PetscSectionGetClosureIndex(section, (PetscObject)dm, clSec, clPoints));
6417:   if (!ornt && *clPoints) {
6418:     PetscInt dof, off;

6420:     PetscCall(PetscSectionGetDof(*clSec, point, &dof));
6421:     PetscCall(PetscSectionGetOffset(*clSec, point, &off));
6422:     PetscCall(ISGetIndices(*clPoints, &cla));
6423:     np  = dof / 2;
6424:     pts = PetscSafePointerPlusOffset((PetscInt *)cla, off);
6425:   } else {
6426:     PetscCall(DMPlexGetTransitiveClosure_Internal(dm, point, ornt, PETSC_TRUE, &np, &pts));
6427:     PetscCall(CompressPoints_Private(section, &np, pts));
6428:   }
6429:   *numPoints = np;
6430:   *points    = pts;
6431:   *clp       = cla;
6432:   PetscFunctionReturn(PETSC_SUCCESS);
6433: }

6435: PetscErrorCode DMPlexRestoreCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
6436: {
6437:   PetscFunctionBeginHot;
6438:   if (!*clPoints) {
6439:     PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, numPoints, points));
6440:   } else {
6441:     PetscCall(ISRestoreIndices(*clPoints, clp));
6442:   }
6443:   *numPoints = 0;
6444:   *points    = NULL;
6445:   *clSec     = NULL;
6446:   *clPoints  = NULL;
6447:   *clp       = NULL;
6448:   PetscFunctionReturn(PETSC_SUCCESS);
6449: }

6451: static inline PetscErrorCode DMPlexVecGetClosure_Static(DM dm, PetscSection section, PetscInt numPoints, const PetscInt points[], const PetscInt clperm[], const PetscScalar vArray[], PetscInt *size, PetscScalar array[])
6452: {
6453:   PetscInt            offset = 0, p;
6454:   const PetscInt    **perms  = NULL;
6455:   const PetscScalar **flips  = NULL;

6457:   PetscFunctionBeginHot;
6458:   *size = 0;
6459:   PetscCall(PetscSectionGetPointSyms(section, numPoints, points, &perms, &flips));
6460:   for (p = 0; p < numPoints; p++) {
6461:     const PetscInt     point = points[2 * p];
6462:     const PetscInt    *perm  = perms ? perms[p] : NULL;
6463:     const PetscScalar *flip  = flips ? flips[p] : NULL;
6464:     PetscInt           dof, off, d;
6465:     const PetscScalar *varr;

6467:     PetscCall(PetscSectionGetDof(section, point, &dof));
6468:     PetscCall(PetscSectionGetOffset(section, point, &off));
6469:     varr = PetscSafePointerPlusOffset(vArray, off);
6470:     if (clperm) {
6471:       if (perm) {
6472:         for (d = 0; d < dof; d++) array[clperm[offset + perm[d]]] = varr[d];
6473:       } else {
6474:         for (d = 0; d < dof; d++) array[clperm[offset + d]] = varr[d];
6475:       }
6476:       if (flip) {
6477:         for (d = 0; d < dof; d++) array[clperm[offset + d]] *= flip[d];
6478:       }
6479:     } else {
6480:       if (perm) {
6481:         for (d = 0; d < dof; d++) array[offset + perm[d]] = varr[d];
6482:       } else {
6483:         for (d = 0; d < dof; d++) array[offset + d] = varr[d];
6484:       }
6485:       if (flip) {
6486:         for (d = 0; d < dof; d++) array[offset + d] *= flip[d];
6487:       }
6488:     }
6489:     offset += dof;
6490:   }
6491:   PetscCall(PetscSectionRestorePointSyms(section, numPoints, points, &perms, &flips));
6492:   *size = offset;
6493:   PetscFunctionReturn(PETSC_SUCCESS);
6494: }

6496: static inline PetscErrorCode DMPlexVecGetClosure_Fields_Static(DM dm, PetscSection section, PetscInt numPoints, const PetscInt points[], PetscInt numFields, const PetscInt clperm[], const PetscScalar vArray[], PetscInt *size, PetscScalar array[])
6497: {
6498:   PetscInt offset = 0, f;

6500:   PetscFunctionBeginHot;
6501:   *size = 0;
6502:   for (f = 0; f < numFields; ++f) {
6503:     PetscInt            p;
6504:     const PetscInt    **perms = NULL;
6505:     const PetscScalar **flips = NULL;

6507:     PetscCall(PetscSectionGetFieldPointSyms(section, f, numPoints, points, &perms, &flips));
6508:     for (p = 0; p < numPoints; p++) {
6509:       const PetscInt     point = points[2 * p];
6510:       PetscInt           fdof, foff, b;
6511:       const PetscScalar *varr;
6512:       const PetscInt    *perm = perms ? perms[p] : NULL;
6513:       const PetscScalar *flip = flips ? flips[p] : NULL;

6515:       PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
6516:       PetscCall(PetscSectionGetFieldOffset(section, point, f, &foff));
6517:       varr = &vArray[foff];
6518:       if (clperm) {
6519:         if (perm) {
6520:           for (b = 0; b < fdof; b++) array[clperm[offset + perm[b]]] = varr[b];
6521:         } else {
6522:           for (b = 0; b < fdof; b++) array[clperm[offset + b]] = varr[b];
6523:         }
6524:         if (flip) {
6525:           for (b = 0; b < fdof; b++) array[clperm[offset + b]] *= flip[b];
6526:         }
6527:       } else {
6528:         if (perm) {
6529:           for (b = 0; b < fdof; b++) array[offset + perm[b]] = varr[b];
6530:         } else {
6531:           for (b = 0; b < fdof; b++) array[offset + b] = varr[b];
6532:         }
6533:         if (flip) {
6534:           for (b = 0; b < fdof; b++) array[offset + b] *= flip[b];
6535:         }
6536:       }
6537:       offset += fdof;
6538:     }
6539:     PetscCall(PetscSectionRestoreFieldPointSyms(section, f, numPoints, points, &perms, &flips));
6540:   }
6541:   *size = offset;
6542:   PetscFunctionReturn(PETSC_SUCCESS);
6543: }

6545: PetscErrorCode DMPlexVecGetOrientedClosure_Internal(DM dm, PetscSection section, PetscBool useClPerm, Vec v, PetscInt point, PetscInt ornt, PetscInt *csize, PetscScalar *values[])
6546: {
6547:   PetscSection    clSection;
6548:   IS              clPoints;
6549:   PetscInt       *points = NULL;
6550:   const PetscInt *clp, *perm = NULL;
6551:   PetscInt        depth, numFields, numPoints, asize;

6553:   PetscFunctionBeginHot;
6555:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
6558:   PetscCall(DMPlexGetDepth(dm, &depth));
6559:   PetscCall(PetscSectionGetNumFields(section, &numFields));
6560:   if (depth == 1 && numFields < 2) {
6561:     PetscCall(DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values));
6562:     PetscFunctionReturn(PETSC_SUCCESS);
6563:   }
6564:   /* Get points */
6565:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, ornt, &numPoints, &points, &clSection, &clPoints, &clp));
6566:   /* Get sizes */
6567:   asize = 0;
6568:   for (PetscInt p = 0; p < numPoints * 2; p += 2) {
6569:     PetscInt dof;
6570:     PetscCall(PetscSectionGetDof(section, points[p], &dof));
6571:     asize += dof;
6572:   }
6573:   if (values) {
6574:     const PetscScalar *vArray;
6575:     PetscInt           size;

6577:     if (*values) {
6578:       PetscCheck(*csize >= asize, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided array size %" PetscInt_FMT " not sufficient to hold closure size %" PetscInt_FMT, *csize, asize);
6579:     } else PetscCall(DMGetWorkArray(dm, asize, MPIU_SCALAR, values));
6580:     if (useClPerm) PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, asize, &perm));
6581:     PetscCall(VecGetArrayRead(v, &vArray));
6582:     /* Get values */
6583:     if (numFields > 0) PetscCall(DMPlexVecGetClosure_Fields_Static(dm, section, numPoints, points, numFields, perm, vArray, &size, *values));
6584:     else PetscCall(DMPlexVecGetClosure_Static(dm, section, numPoints, points, perm, vArray, &size, *values));
6585:     PetscCheck(asize == size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Section size %" PetscInt_FMT " does not match Vec closure size %" PetscInt_FMT, asize, size);
6586:     /* Cleanup array */
6587:     PetscCall(VecRestoreArrayRead(v, &vArray));
6588:   }
6589:   if (csize) *csize = asize;
6590:   /* Cleanup points */
6591:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
6592:   PetscFunctionReturn(PETSC_SUCCESS);
6593: }

6595: /*@C
6596:   DMPlexVecGetClosure - Get an array of the values on the closure of 'point'

6598:   Not collective

6600:   Input Parameters:
6601: + dm      - The `DM`
6602: . section - The section describing the layout in `v`, or `NULL` to use the default section
6603: . v       - The local vector
6604: - point   - The point in the `DM`

6606:   Input/Output Parameters:
6607: + csize  - The size of the input values array, or `NULL`; on output the number of values in the closure
6608: - values - An array to use for the values, or *values = `NULL` to have it allocated automatically;
6609:            if the user provided `NULL`, it is a borrowed array and should not be freed, use  `DMPlexVecRestoreClosure()` to return it

6611:   Level: intermediate

6613:   Notes:
6614:   `DMPlexVecGetClosure()`/`DMPlexVecRestoreClosure()` only allocates the values array if it set to `NULL` in the
6615:   calling function. This is because `DMPlexVecGetClosure()` is typically called in the inner loop of a `Vec` or `Mat`
6616:   assembly function, and a user may already have allocated storage for this operation.

6618:   A typical use could be
6619: .vb
6620:    values = NULL;
6621:    PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values));
6622:    for (cl = 0; cl < clSize; ++cl) {
6623:      <Compute on closure>
6624:    }
6625:    PetscCall(DMPlexVecRestoreClosure(dm, NULL, v, p, &clSize, &values));
6626: .ve
6627:   or
6628: .vb
6629:    PetscMalloc1(clMaxSize, &values);
6630:    for (p = pStart; p < pEnd; ++p) {
6631:      clSize = clMaxSize;
6632:      PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values));
6633:      for (cl = 0; cl < clSize; ++cl) {
6634:        <Compute on closure>
6635:      }
6636:    }
6637:    PetscFree(values);
6638: .ve

6640:   Fortran Notes:
6641:   The `csize` argument is not present in the Fortran binding.

6643:   `values` must be declared with
6644: .vb
6645:   PetscScalar,dimension(:),pointer   :: values
6646: .ve
6647:   and it will be allocated internally by PETSc to hold the values returned

6649: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecRestoreClosure()`, `DMPlexVecSetClosure()`, `DMPlexMatSetClosure()`
6650: @*/
6651: PetscErrorCode DMPlexVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
6652: {
6653:   PetscFunctionBeginHot;
6654:   PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, section, PETSC_TRUE, v, point, 0, csize, values));
6655:   PetscFunctionReturn(PETSC_SUCCESS);
6656: }

6658: PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt depth, PetscInt *csize, PetscScalar *values[])
6659: {
6660:   DMLabel            depthLabel;
6661:   PetscSection       clSection;
6662:   IS                 clPoints;
6663:   PetscScalar       *array;
6664:   const PetscScalar *vArray;
6665:   PetscInt          *points = NULL;
6666:   const PetscInt    *clp, *perm = NULL;
6667:   PetscInt           mdepth, numFields, numPoints, Np = 0, p, clsize, size;

6669:   PetscFunctionBeginHot;
6671:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
6674:   PetscCall(DMPlexGetDepth(dm, &mdepth));
6675:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
6676:   PetscCall(PetscSectionGetNumFields(section, &numFields));
6677:   if (mdepth == 1 && numFields < 2) {
6678:     PetscCall(DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values));
6679:     PetscFunctionReturn(PETSC_SUCCESS);
6680:   }
6681:   /* Get points */
6682:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &numPoints, &points, &clSection, &clPoints, &clp));
6683:   for (clsize = 0, p = 0; p < Np; p++) {
6684:     PetscInt dof;
6685:     PetscCall(PetscSectionGetDof(section, points[2 * p], &dof));
6686:     clsize += dof;
6687:   }
6688:   PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, clsize, &perm));
6689:   /* Filter points */
6690:   for (p = 0; p < numPoints * 2; p += 2) {
6691:     PetscInt dep;

6693:     PetscCall(DMLabelGetValue(depthLabel, points[p], &dep));
6694:     if (dep != depth) continue;
6695:     points[Np * 2 + 0] = points[p];
6696:     points[Np * 2 + 1] = points[p + 1];
6697:     ++Np;
6698:   }
6699:   /* Get array */
6700:   if (!values || !*values) {
6701:     PetscInt asize = 0, dof;

6703:     for (p = 0; p < Np * 2; p += 2) {
6704:       PetscCall(PetscSectionGetDof(section, points[p], &dof));
6705:       asize += dof;
6706:     }
6707:     if (!values) {
6708:       PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
6709:       if (csize) *csize = asize;
6710:       PetscFunctionReturn(PETSC_SUCCESS);
6711:     }
6712:     PetscCall(DMGetWorkArray(dm, asize, MPIU_SCALAR, &array));
6713:   } else {
6714:     array = *values;
6715:   }
6716:   PetscCall(VecGetArrayRead(v, &vArray));
6717:   /* Get values */
6718:   if (numFields > 0) PetscCall(DMPlexVecGetClosure_Fields_Static(dm, section, Np, points, numFields, perm, vArray, &size, array));
6719:   else PetscCall(DMPlexVecGetClosure_Static(dm, section, Np, points, perm, vArray, &size, array));
6720:   /* Cleanup points */
6721:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
6722:   /* Cleanup array */
6723:   PetscCall(VecRestoreArrayRead(v, &vArray));
6724:   if (!*values) {
6725:     if (csize) *csize = size;
6726:     *values = array;
6727:   } else {
6728:     PetscCheck(size <= *csize, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %" PetscInt_FMT " < actual size %" PetscInt_FMT, *csize, size);
6729:     *csize = size;
6730:   }
6731:   PetscFunctionReturn(PETSC_SUCCESS);
6732: }

6734: /*@C
6735:   DMPlexVecRestoreClosure - Restore the array of the values on the closure of 'point' obtained with `DMPlexVecGetClosure()`

6737:   Not collective

6739:   Input Parameters:
6740: + dm      - The `DM`
6741: . section - The section describing the layout in `v`, or `NULL` to use the default section
6742: . v       - The local vector
6743: . point   - The point in the `DM`
6744: . csize   - The number of values in the closure, or `NULL`
6745: - values  - The array of values

6747:   Level: intermediate

6749:   Note:
6750:   The array values are discarded and not copied back into `v`. In order to copy values back to `v`, use `DMPlexVecSetClosure()`

6752:   Fortran Note:
6753:   The `csize` argument is not present in the Fortran binding since it is internal to the array.

6755: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecGetClosure()`, `DMPlexVecSetClosure()`, `DMPlexMatSetClosure()`
6756: @*/
6757: PetscErrorCode DMPlexVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
6758: {
6759:   PetscInt size = 0;

6761:   PetscFunctionBegin;
6762:   /* Should work without recalculating size */
6763:   PetscCall(DMRestoreWorkArray(dm, size, MPIU_SCALAR, (void *)values));
6764:   *values = NULL;
6765:   PetscFunctionReturn(PETSC_SUCCESS);
6766: }

6768: static inline void add(PetscScalar *x, PetscScalar y)
6769: {
6770:   *x += y;
6771: }
6772: static inline void insert(PetscScalar *x, PetscScalar y)
6773: {
6774:   *x = y;
6775: }

6777: static inline PetscErrorCode updatePoint_private(PetscSection section, PetscInt point, PetscInt dof, void (*fuse)(PetscScalar *, PetscScalar), PetscBool setBC, const PetscInt perm[], const PetscScalar flip[], const PetscInt clperm[], const PetscScalar values[], PetscInt offset, PetscScalar array[])
6778: {
6779:   PetscInt        cdof;  /* The number of constraints on this point */
6780:   const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6781:   PetscScalar    *a;
6782:   PetscInt        off, cind = 0, k;

6784:   PetscFunctionBegin;
6785:   PetscCall(PetscSectionGetConstraintDof(section, point, &cdof));
6786:   PetscCall(PetscSectionGetOffset(section, point, &off));
6787:   a = &array[off];
6788:   if (!cdof || setBC) {
6789:     if (clperm) {
6790:       if (perm) {
6791:         for (k = 0; k < dof; ++k) fuse(&a[k], values[clperm[offset + perm[k]]] * (flip ? flip[perm[k]] : 1.));
6792:       } else {
6793:         for (k = 0; k < dof; ++k) fuse(&a[k], values[clperm[offset + k]] * (flip ? flip[k] : 1.));
6794:       }
6795:     } else {
6796:       if (perm) {
6797:         for (k = 0; k < dof; ++k) fuse(&a[k], values[offset + perm[k]] * (flip ? flip[perm[k]] : 1.));
6798:       } else {
6799:         for (k = 0; k < dof; ++k) fuse(&a[k], values[offset + k] * (flip ? flip[k] : 1.));
6800:       }
6801:     }
6802:   } else {
6803:     PetscCall(PetscSectionGetConstraintIndices(section, point, &cdofs));
6804:     if (clperm) {
6805:       if (perm) {
6806:         for (k = 0; k < dof; ++k) {
6807:           if ((cind < cdof) && (k == cdofs[cind])) {
6808:             ++cind;
6809:             continue;
6810:           }
6811:           fuse(&a[k], values[clperm[offset + perm[k]]] * (flip ? flip[perm[k]] : 1.));
6812:         }
6813:       } else {
6814:         for (k = 0; k < dof; ++k) {
6815:           if ((cind < cdof) && (k == cdofs[cind])) {
6816:             ++cind;
6817:             continue;
6818:           }
6819:           fuse(&a[k], values[clperm[offset + k]] * (flip ? flip[k] : 1.));
6820:         }
6821:       }
6822:     } else {
6823:       if (perm) {
6824:         for (k = 0; k < dof; ++k) {
6825:           if ((cind < cdof) && (k == cdofs[cind])) {
6826:             ++cind;
6827:             continue;
6828:           }
6829:           fuse(&a[k], values[offset + perm[k]] * (flip ? flip[perm[k]] : 1.));
6830:         }
6831:       } else {
6832:         for (k = 0; k < dof; ++k) {
6833:           if ((cind < cdof) && (k == cdofs[cind])) {
6834:             ++cind;
6835:             continue;
6836:           }
6837:           fuse(&a[k], values[offset + k] * (flip ? flip[k] : 1.));
6838:         }
6839:       }
6840:     }
6841:   }
6842:   PetscFunctionReturn(PETSC_SUCCESS);
6843: }

6845: static inline PetscErrorCode updatePointBC_private(PetscSection section, PetscInt point, PetscInt dof, void (*fuse)(PetscScalar *, PetscScalar), const PetscInt perm[], const PetscScalar flip[], const PetscInt clperm[], const PetscScalar values[], PetscInt offset, PetscScalar array[])
6846: {
6847:   PetscInt        cdof;  /* The number of constraints on this point */
6848:   const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6849:   PetscScalar    *a;
6850:   PetscInt        off, cind = 0, k;

6852:   PetscFunctionBegin;
6853:   PetscCall(PetscSectionGetConstraintDof(section, point, &cdof));
6854:   PetscCall(PetscSectionGetOffset(section, point, &off));
6855:   a = &array[off];
6856:   if (cdof) {
6857:     PetscCall(PetscSectionGetConstraintIndices(section, point, &cdofs));
6858:     if (clperm) {
6859:       if (perm) {
6860:         for (k = 0; k < dof; ++k) {
6861:           if ((cind < cdof) && (k == cdofs[cind])) {
6862:             fuse(&a[k], values[clperm[offset + perm[k]]] * (flip ? flip[perm[k]] : 1.));
6863:             cind++;
6864:           }
6865:         }
6866:       } else {
6867:         for (k = 0; k < dof; ++k) {
6868:           if ((cind < cdof) && (k == cdofs[cind])) {
6869:             fuse(&a[k], values[clperm[offset + k]] * (flip ? flip[k] : 1.));
6870:             cind++;
6871:           }
6872:         }
6873:       }
6874:     } else {
6875:       if (perm) {
6876:         for (k = 0; k < dof; ++k) {
6877:           if ((cind < cdof) && (k == cdofs[cind])) {
6878:             fuse(&a[k], values[offset + perm[k]] * (flip ? flip[perm[k]] : 1.));
6879:             cind++;
6880:           }
6881:         }
6882:       } else {
6883:         for (k = 0; k < dof; ++k) {
6884:           if ((cind < cdof) && (k == cdofs[cind])) {
6885:             fuse(&a[k], values[offset + k] * (flip ? flip[k] : 1.));
6886:             cind++;
6887:           }
6888:         }
6889:       }
6890:     }
6891:   }
6892:   PetscFunctionReturn(PETSC_SUCCESS);
6893: }

6895: static inline PetscErrorCode updatePointFields_private(PetscSection section, PetscInt point, const PetscInt *perm, const PetscScalar *flip, PetscInt f, void (*fuse)(PetscScalar *, PetscScalar), PetscBool setBC, const PetscInt clperm[], const PetscScalar values[], PetscInt *offset, PetscScalar array[])
6896: {
6897:   PetscScalar    *a;
6898:   PetscInt        fdof, foff, fcdof, foffset = *offset;
6899:   const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
6900:   PetscInt        cind = 0, b;

6902:   PetscFunctionBegin;
6903:   PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
6904:   PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &fcdof));
6905:   PetscCall(PetscSectionGetFieldOffset(section, point, f, &foff));
6906:   a = &array[foff];
6907:   if (!fcdof || setBC) {
6908:     if (clperm) {
6909:       if (perm) {
6910:         for (b = 0; b < fdof; b++) fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6911:       } else {
6912:         for (b = 0; b < fdof; b++) fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
6913:       }
6914:     } else {
6915:       if (perm) {
6916:         for (b = 0; b < fdof; b++) fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
6917:       } else {
6918:         for (b = 0; b < fdof; b++) fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
6919:       }
6920:     }
6921:   } else {
6922:     PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
6923:     if (clperm) {
6924:       if (perm) {
6925:         for (b = 0; b < fdof; b++) {
6926:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6927:             ++cind;
6928:             continue;
6929:           }
6930:           fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6931:         }
6932:       } else {
6933:         for (b = 0; b < fdof; b++) {
6934:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6935:             ++cind;
6936:             continue;
6937:           }
6938:           fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
6939:         }
6940:       }
6941:     } else {
6942:       if (perm) {
6943:         for (b = 0; b < fdof; b++) {
6944:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6945:             ++cind;
6946:             continue;
6947:           }
6948:           fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
6949:         }
6950:       } else {
6951:         for (b = 0; b < fdof; b++) {
6952:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6953:             ++cind;
6954:             continue;
6955:           }
6956:           fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
6957:         }
6958:       }
6959:     }
6960:   }
6961:   *offset += fdof;
6962:   PetscFunctionReturn(PETSC_SUCCESS);
6963: }

6965: static inline PetscErrorCode updatePointFieldsBC_private(PetscSection section, PetscInt point, const PetscInt perm[], const PetscScalar flip[], PetscInt f, PetscInt Ncc, const PetscInt comps[], void (*fuse)(PetscScalar *, PetscScalar), const PetscInt clperm[], const PetscScalar values[], PetscInt *offset, PetscScalar array[])
6966: {
6967:   PetscScalar    *a;
6968:   PetscInt        fdof, foff, fcdof, foffset = *offset;
6969:   const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
6970:   PetscInt        Nc, cind = 0, ncind = 0, b;
6971:   PetscBool       ncSet, fcSet;

6973:   PetscFunctionBegin;
6974:   PetscCall(PetscSectionGetFieldComponents(section, f, &Nc));
6975:   PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
6976:   PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &fcdof));
6977:   PetscCall(PetscSectionGetFieldOffset(section, point, f, &foff));
6978:   a = &array[foff];
6979:   if (fcdof) {
6980:     /* We just override fcdof and fcdofs with Ncc and comps */
6981:     PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
6982:     if (clperm) {
6983:       if (perm) {
6984:         if (comps) {
6985:           for (b = 0; b < fdof; b++) {
6986:             ncSet = fcSet = PETSC_FALSE;
6987:             if (b % Nc == comps[ncind]) {
6988:               ncind = (ncind + 1) % Ncc;
6989:               ncSet = PETSC_TRUE;
6990:             }
6991:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6992:               ++cind;
6993:               fcSet = PETSC_TRUE;
6994:             }
6995:             if (ncSet && fcSet) fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6996:           }
6997:         } else {
6998:           for (b = 0; b < fdof; b++) {
6999:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7000:               fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
7001:               ++cind;
7002:             }
7003:           }
7004:         }
7005:       } else {
7006:         if (comps) {
7007:           for (b = 0; b < fdof; b++) {
7008:             ncSet = fcSet = PETSC_FALSE;
7009:             if (b % Nc == comps[ncind]) {
7010:               ncind = (ncind + 1) % Ncc;
7011:               ncSet = PETSC_TRUE;
7012:             }
7013:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7014:               ++cind;
7015:               fcSet = PETSC_TRUE;
7016:             }
7017:             if (ncSet && fcSet) fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
7018:           }
7019:         } else {
7020:           for (b = 0; b < fdof; b++) {
7021:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7022:               fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
7023:               ++cind;
7024:             }
7025:           }
7026:         }
7027:       }
7028:     } else {
7029:       if (perm) {
7030:         if (comps) {
7031:           for (b = 0; b < fdof; b++) {
7032:             ncSet = fcSet = PETSC_FALSE;
7033:             if (b % Nc == comps[ncind]) {
7034:               ncind = (ncind + 1) % Ncc;
7035:               ncSet = PETSC_TRUE;
7036:             }
7037:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7038:               ++cind;
7039:               fcSet = PETSC_TRUE;
7040:             }
7041:             if (ncSet && fcSet) fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
7042:           }
7043:         } else {
7044:           for (b = 0; b < fdof; b++) {
7045:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7046:               fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
7047:               ++cind;
7048:             }
7049:           }
7050:         }
7051:       } else {
7052:         if (comps) {
7053:           for (b = 0; b < fdof; b++) {
7054:             ncSet = fcSet = PETSC_FALSE;
7055:             if (b % Nc == comps[ncind]) {
7056:               ncind = (ncind + 1) % Ncc;
7057:               ncSet = PETSC_TRUE;
7058:             }
7059:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7060:               ++cind;
7061:               fcSet = PETSC_TRUE;
7062:             }
7063:             if (ncSet && fcSet) fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
7064:           }
7065:         } else {
7066:           for (b = 0; b < fdof; b++) {
7067:             if ((cind < fcdof) && (b == fcdofs[cind])) {
7068:               fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
7069:               ++cind;
7070:             }
7071:           }
7072:         }
7073:       }
7074:     }
7075:   }
7076:   *offset += fdof;
7077:   PetscFunctionReturn(PETSC_SUCCESS);
7078: }

7080: static inline PetscErrorCode DMPlexVecSetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
7081: {
7082:   PetscScalar    *array;
7083:   const PetscInt *cone, *coneO;
7084:   PetscInt        pStart, pEnd, p, numPoints, off, dof;

7086:   PetscFunctionBeginHot;
7087:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
7088:   PetscCall(DMPlexGetConeSize(dm, point, &numPoints));
7089:   PetscCall(DMPlexGetCone(dm, point, &cone));
7090:   PetscCall(DMPlexGetConeOrientation(dm, point, &coneO));
7091:   PetscCall(VecGetArray(v, &array));
7092:   for (p = 0, off = 0; p <= numPoints; ++p, off += dof) {
7093:     const PetscInt cp = !p ? point : cone[p - 1];
7094:     const PetscInt o  = !p ? 0 : coneO[p - 1];

7096:     if ((cp < pStart) || (cp >= pEnd)) {
7097:       dof = 0;
7098:       continue;
7099:     }
7100:     PetscCall(PetscSectionGetDof(section, cp, &dof));
7101:     /* ADD_VALUES */
7102:     {
7103:       const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
7104:       PetscScalar    *a;
7105:       PetscInt        cdof, coff, cind = 0, k;

7107:       PetscCall(PetscSectionGetConstraintDof(section, cp, &cdof));
7108:       PetscCall(PetscSectionGetOffset(section, cp, &coff));
7109:       a = &array[coff];
7110:       if (!cdof) {
7111:         if (o >= 0) {
7112:           for (k = 0; k < dof; ++k) a[k] += values[off + k];
7113:         } else {
7114:           for (k = 0; k < dof; ++k) a[k] += values[off + dof - k - 1];
7115:         }
7116:       } else {
7117:         PetscCall(PetscSectionGetConstraintIndices(section, cp, &cdofs));
7118:         if (o >= 0) {
7119:           for (k = 0; k < dof; ++k) {
7120:             if ((cind < cdof) && (k == cdofs[cind])) {
7121:               ++cind;
7122:               continue;
7123:             }
7124:             a[k] += values[off + k];
7125:           }
7126:         } else {
7127:           for (k = 0; k < dof; ++k) {
7128:             if ((cind < cdof) && (k == cdofs[cind])) {
7129:               ++cind;
7130:               continue;
7131:             }
7132:             a[k] += values[off + dof - k - 1];
7133:           }
7134:         }
7135:       }
7136:     }
7137:   }
7138:   PetscCall(VecRestoreArray(v, &array));
7139:   PetscFunctionReturn(PETSC_SUCCESS);
7140: }

7142: /*@C
7143:   DMPlexVecSetClosure - Set an array of the values on the closure of `point`

7145:   Not collective

7147:   Input Parameters:
7148: + dm      - The `DM`
7149: . section - The section describing the layout in `v`, or `NULL` to use the default section
7150: . v       - The local vector
7151: . point   - The point in the `DM`
7152: . values  - The array of values
7153: - mode    - The insert mode. One of `INSERT_ALL_VALUES`, `ADD_ALL_VALUES`, `INSERT_VALUES`, `ADD_VALUES`, `INSERT_BC_VALUES`, and `ADD_BC_VALUES`,
7154:             where `INSERT_ALL_VALUES` and `ADD_ALL_VALUES` also overwrite boundary conditions.

7156:   Level: intermediate

7158:   Note:
7159:   Usually the input arrays were obtained with `DMPlexVecGetClosure()`

7161:   Fortran Note:
7162:   `values` must be declared with
7163: .vb
7164:   PetscScalar,dimension(:),pointer   :: values
7165: .ve

7167: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecGetClosure()`, `DMPlexMatSetClosure()`
7168: @*/
7169: PetscErrorCode DMPlexVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
7170: {
7171:   PetscSection    clSection;
7172:   IS              clPoints;
7173:   PetscScalar    *array;
7174:   PetscInt       *points = NULL;
7175:   const PetscInt *clp, *clperm = NULL;
7176:   PetscInt        depth, numFields, numPoints, p, clsize;

7178:   PetscFunctionBeginHot;
7180:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
7183:   PetscCall(DMPlexGetDepth(dm, &depth));
7184:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7185:   if (depth == 1 && numFields < 2 && mode == ADD_VALUES) {
7186:     PetscCall(DMPlexVecSetClosure_Depth1_Static(dm, section, v, point, values, mode));
7187:     PetscFunctionReturn(PETSC_SUCCESS);
7188:   }
7189:   /* Get points */
7190:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &numPoints, &points, &clSection, &clPoints, &clp));
7191:   for (clsize = 0, p = 0; p < numPoints; p++) {
7192:     PetscInt dof;
7193:     PetscCall(PetscSectionGetDof(section, points[2 * p], &dof));
7194:     clsize += dof;
7195:   }
7196:   PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, clsize, &clperm));
7197:   /* Get array */
7198:   PetscCall(VecGetArray(v, &array));
7199:   /* Get values */
7200:   if (numFields > 0) {
7201:     PetscInt offset = 0, f;
7202:     for (f = 0; f < numFields; ++f) {
7203:       const PetscInt    **perms = NULL;
7204:       const PetscScalar **flips = NULL;

7206:       PetscCall(PetscSectionGetFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7207:       switch (mode) {
7208:       case INSERT_VALUES:
7209:         for (p = 0; p < numPoints; p++) {
7210:           const PetscInt     point = points[2 * p];
7211:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7212:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7213:           PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, clperm, values, &offset, array));
7214:         }
7215:         break;
7216:       case INSERT_ALL_VALUES:
7217:         for (p = 0; p < numPoints; p++) {
7218:           const PetscInt     point = points[2 * p];
7219:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7220:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7221:           PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, clperm, values, &offset, array));
7222:         }
7223:         break;
7224:       case INSERT_BC_VALUES:
7225:         for (p = 0; p < numPoints; p++) {
7226:           const PetscInt     point = points[2 * p];
7227:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7228:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7229:           PetscCall(updatePointFieldsBC_private(section, point, perm, flip, f, -1, NULL, insert, clperm, values, &offset, array));
7230:         }
7231:         break;
7232:       case ADD_VALUES:
7233:         for (p = 0; p < numPoints; p++) {
7234:           const PetscInt     point = points[2 * p];
7235:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7236:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7237:           PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, clperm, values, &offset, array));
7238:         }
7239:         break;
7240:       case ADD_ALL_VALUES:
7241:         for (p = 0; p < numPoints; p++) {
7242:           const PetscInt     point = points[2 * p];
7243:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7244:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7245:           PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, clperm, values, &offset, array));
7246:         }
7247:         break;
7248:       case ADD_BC_VALUES:
7249:         for (p = 0; p < numPoints; p++) {
7250:           const PetscInt     point = points[2 * p];
7251:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7252:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7253:           PetscCall(updatePointFieldsBC_private(section, point, perm, flip, f, -1, NULL, add, clperm, values, &offset, array));
7254:         }
7255:         break;
7256:       default:
7257:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
7258:       }
7259:       PetscCall(PetscSectionRestoreFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7260:     }
7261:   } else {
7262:     PetscInt            dof, off;
7263:     const PetscInt    **perms = NULL;
7264:     const PetscScalar **flips = NULL;

7266:     PetscCall(PetscSectionGetPointSyms(section, numPoints, points, &perms, &flips));
7267:     switch (mode) {
7268:     case INSERT_VALUES:
7269:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7270:         const PetscInt     point = points[2 * p];
7271:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7272:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7273:         PetscCall(PetscSectionGetDof(section, point, &dof));
7274:         PetscCall(updatePoint_private(section, point, dof, insert, PETSC_FALSE, perm, flip, clperm, values, off, array));
7275:       }
7276:       break;
7277:     case INSERT_ALL_VALUES:
7278:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7279:         const PetscInt     point = points[2 * p];
7280:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7281:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7282:         PetscCall(PetscSectionGetDof(section, point, &dof));
7283:         PetscCall(updatePoint_private(section, point, dof, insert, PETSC_TRUE, perm, flip, clperm, values, off, array));
7284:       }
7285:       break;
7286:     case INSERT_BC_VALUES:
7287:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7288:         const PetscInt     point = points[2 * p];
7289:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7290:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7291:         PetscCall(PetscSectionGetDof(section, point, &dof));
7292:         PetscCall(updatePointBC_private(section, point, dof, insert, perm, flip, clperm, values, off, array));
7293:       }
7294:       break;
7295:     case ADD_VALUES:
7296:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7297:         const PetscInt     point = points[2 * p];
7298:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7299:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7300:         PetscCall(PetscSectionGetDof(section, point, &dof));
7301:         PetscCall(updatePoint_private(section, point, dof, add, PETSC_FALSE, perm, flip, clperm, values, off, array));
7302:       }
7303:       break;
7304:     case ADD_ALL_VALUES:
7305:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7306:         const PetscInt     point = points[2 * p];
7307:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7308:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7309:         PetscCall(PetscSectionGetDof(section, point, &dof));
7310:         PetscCall(updatePoint_private(section, point, dof, add, PETSC_TRUE, perm, flip, clperm, values, off, array));
7311:       }
7312:       break;
7313:     case ADD_BC_VALUES:
7314:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7315:         const PetscInt     point = points[2 * p];
7316:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7317:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7318:         PetscCall(PetscSectionGetDof(section, point, &dof));
7319:         PetscCall(updatePointBC_private(section, point, dof, add, perm, flip, clperm, values, off, array));
7320:       }
7321:       break;
7322:     default:
7323:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
7324:     }
7325:     PetscCall(PetscSectionRestorePointSyms(section, numPoints, points, &perms, &flips));
7326:   }
7327:   /* Cleanup points */
7328:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
7329:   /* Cleanup array */
7330:   PetscCall(VecRestoreArray(v, &array));
7331:   PetscFunctionReturn(PETSC_SUCCESS);
7332: }

7334: /* Check whether the given point is in the label. If not, update the offset to skip this point */
7335: static inline PetscErrorCode CheckPoint_Private(DMLabel label, PetscInt labelId, PetscSection section, PetscInt point, PetscInt f, PetscInt *offset, PetscBool *contains)
7336: {
7337:   PetscFunctionBegin;
7338:   *contains = PETSC_TRUE;
7339:   if (label) {
7340:     PetscInt fdof;

7342:     PetscCall(DMLabelStratumHasPoint(label, labelId, point, contains));
7343:     if (!*contains) {
7344:       PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
7345:       *offset += fdof;
7346:       PetscFunctionReturn(PETSC_SUCCESS);
7347:     }
7348:   }
7349:   PetscFunctionReturn(PETSC_SUCCESS);
7350: }

7352: /* Unlike DMPlexVecSetClosure(), this uses plex-native closure permutation, not a user-specified permutation such as DMPlexSetClosurePermutationTensor(). */
7353: PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM dm, PetscSection section, Vec v, PetscBool fieldActive[], PetscInt point, PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt labelId, const PetscScalar values[], InsertMode mode)
7354: {
7355:   PetscSection    clSection;
7356:   IS              clPoints;
7357:   PetscScalar    *array;
7358:   PetscInt       *points = NULL;
7359:   const PetscInt *clp;
7360:   PetscInt        numFields, numPoints, p;
7361:   PetscInt        offset = 0, f;

7363:   PetscFunctionBeginHot;
7365:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
7368:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7369:   /* Get points */
7370:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &numPoints, &points, &clSection, &clPoints, &clp));
7371:   /* Get array */
7372:   PetscCall(VecGetArray(v, &array));
7373:   /* Get values */
7374:   for (f = 0; f < numFields; ++f) {
7375:     const PetscInt    **perms = NULL;
7376:     const PetscScalar **flips = NULL;
7377:     PetscBool           contains;

7379:     if (!fieldActive[f]) {
7380:       for (p = 0; p < numPoints * 2; p += 2) {
7381:         PetscInt fdof;
7382:         PetscCall(PetscSectionGetFieldDof(section, points[p], f, &fdof));
7383:         offset += fdof;
7384:       }
7385:       continue;
7386:     }
7387:     PetscCall(PetscSectionGetFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7388:     switch (mode) {
7389:     case INSERT_VALUES:
7390:       for (p = 0; p < numPoints; p++) {
7391:         const PetscInt     point = points[2 * p];
7392:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7393:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7394:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7395:         if (!contains) continue;
7396:         PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, NULL, values, &offset, array));
7397:       }
7398:       break;
7399:     case INSERT_ALL_VALUES:
7400:       for (p = 0; p < numPoints; p++) {
7401:         const PetscInt     point = points[2 * p];
7402:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7403:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7404:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7405:         if (!contains) continue;
7406:         PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, NULL, values, &offset, array));
7407:       }
7408:       break;
7409:     case INSERT_BC_VALUES:
7410:       for (p = 0; p < numPoints; p++) {
7411:         const PetscInt     point = points[2 * p];
7412:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7413:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7414:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7415:         if (!contains) continue;
7416:         PetscCall(updatePointFieldsBC_private(section, point, perm, flip, f, Ncc, comps, insert, NULL, values, &offset, array));
7417:       }
7418:       break;
7419:     case ADD_VALUES:
7420:       for (p = 0; p < numPoints; p++) {
7421:         const PetscInt     point = points[2 * p];
7422:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7423:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7424:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7425:         if (!contains) continue;
7426:         PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, NULL, values, &offset, array));
7427:       }
7428:       break;
7429:     case ADD_ALL_VALUES:
7430:       for (p = 0; p < numPoints; p++) {
7431:         const PetscInt     point = points[2 * p];
7432:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7433:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7434:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7435:         if (!contains) continue;
7436:         PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, NULL, values, &offset, array));
7437:       }
7438:       break;
7439:     default:
7440:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
7441:     }
7442:     PetscCall(PetscSectionRestoreFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7443:   }
7444:   /* Cleanup points */
7445:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
7446:   /* Cleanup array */
7447:   PetscCall(VecRestoreArray(v, &array));
7448:   PetscFunctionReturn(PETSC_SUCCESS);
7449: }

7451: static PetscErrorCode DMPlexPrintMatSetValues(PetscViewer viewer, Mat A, PetscInt point, PetscInt numRIndices, const PetscInt rindices[], PetscInt numCIndices, const PetscInt cindices[], const PetscScalar values[])
7452: {
7453:   PetscMPIInt rank;
7454:   PetscInt    i, j;

7456:   PetscFunctionBegin;
7457:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
7458:   PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]mat for point %" PetscInt_FMT "\n", rank, point));
7459:   for (i = 0; i < numRIndices; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]mat row indices[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, i, rindices[i]));
7460:   for (i = 0; i < numCIndices; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]mat col indices[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, i, cindices[i]));
7461:   numCIndices = numCIndices ? numCIndices : numRIndices;
7462:   if (!values) PetscFunctionReturn(PETSC_SUCCESS);
7463:   for (i = 0; i < numRIndices; i++) {
7464:     PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]", rank));
7465:     for (j = 0; j < numCIndices; j++) {
7466: #if defined(PETSC_USE_COMPLEX)
7467:       PetscCall(PetscViewerASCIIPrintf(viewer, " (%g,%g)", (double)PetscRealPart(values[i * numCIndices + j]), (double)PetscImaginaryPart(values[i * numCIndices + j])));
7468: #else
7469:       PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)values[i * numCIndices + j]));
7470: #endif
7471:     }
7472:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
7473:   }
7474:   PetscFunctionReturn(PETSC_SUCCESS);
7475: }

7477: /*
7478:   DMPlexGetIndicesPoint_Internal - Add the indices for dofs on a point to an index array

7480:   Input Parameters:
7481: + section - The section for this data layout
7482: . islocal - Is the section (and thus indices being requested) local or global?
7483: . point   - The point contributing dofs with these indices
7484: . off     - The global offset of this point
7485: . loff    - The local offset of each field
7486: . setBC   - The flag determining whether to include indices of boundary values
7487: . perm    - A permutation of the dofs on this point, or NULL
7488: - indperm - A permutation of the entire indices array, or NULL

7490:   Output Parameter:
7491: . indices - Indices for dofs on this point

7493:   Level: developer

7495:   Note: The indices could be local or global, depending on the value of 'off'.
7496: */
7497: PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection section, PetscBool islocal, PetscInt point, PetscInt off, PetscInt *loff, PetscBool setBC, const PetscInt perm[], const PetscInt indperm[], PetscInt indices[])
7498: {
7499:   PetscInt        dof;   /* The number of unknowns on this point */
7500:   PetscInt        cdof;  /* The number of constraints on this point */
7501:   const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
7502:   PetscInt        cind = 0, k;

7504:   PetscFunctionBegin;
7505:   PetscCheck(islocal || !setBC, PetscObjectComm((PetscObject)section), PETSC_ERR_ARG_INCOMP, "setBC incompatible with global indices; use a local section or disable setBC");
7506:   PetscCall(PetscSectionGetDof(section, point, &dof));
7507:   PetscCall(PetscSectionGetConstraintDof(section, point, &cdof));
7508:   if (!cdof || setBC) {
7509:     for (k = 0; k < dof; ++k) {
7510:       const PetscInt preind = perm ? *loff + perm[k] : *loff + k;
7511:       const PetscInt ind    = indperm ? indperm[preind] : preind;

7513:       indices[ind] = off + k;
7514:     }
7515:   } else {
7516:     PetscCall(PetscSectionGetConstraintIndices(section, point, &cdofs));
7517:     for (k = 0; k < dof; ++k) {
7518:       const PetscInt preind = perm ? *loff + perm[k] : *loff + k;
7519:       const PetscInt ind    = indperm ? indperm[preind] : preind;

7521:       if ((cind < cdof) && (k == cdofs[cind])) {
7522:         /* Insert check for returning constrained indices */
7523:         indices[ind] = -(off + k + 1);
7524:         ++cind;
7525:       } else {
7526:         indices[ind] = off + k - (islocal ? 0 : cind);
7527:       }
7528:     }
7529:   }
7530:   *loff += dof;
7531:   PetscFunctionReturn(PETSC_SUCCESS);
7532: }

7534: /*
7535:  DMPlexGetIndicesPointFields_Internal - gets section indices for a point in its canonical ordering.

7537:  Input Parameters:
7538: + section - a section (global or local)
7539: - islocal - `PETSC_TRUE` if requesting local indices (i.e., section is local); `PETSC_FALSE` for global
7540: . point - point within section
7541: . off - The offset of this point in the (local or global) indexed space - should match islocal and (usually) the section
7542: . foffs - array of length numFields containing the offset in canonical point ordering (the location in indices) of each field
7543: . setBC - identify constrained (boundary condition) points via involution.
7544: . perms - perms[f][permsoff][:] is a permutation of dofs within each field
7545: . permsoff - offset
7546: - indperm - index permutation

7548:  Output Parameter:
7549: . foffs - each entry is incremented by the number of (unconstrained if setBC=FALSE) dofs in that field
7550: . indices - array to hold indices (as defined by section) of each dof associated with point

7552:  Notes:
7553:  If section is local and setBC=true, there is no distinction between constrained and unconstrained dofs.
7554:  If section is local and setBC=false, the indices for constrained points are the involution -(i+1) of their position
7555:  in the local vector.

7557:  If section is global and setBC=false, the indices for constrained points are negative (and their value is not
7558:  significant).  It is invalid to call with a global section and setBC=true.

7560:  Developer Note:
7561:  The section is only used for field layout, so islocal is technically a statement about the offset (off).  At some point
7562:  in the future, global sections may have fields set, in which case we could pass the global section and obtain the
7563:  offset could be obtained from the section instead of passing it explicitly as we do now.

7565:  Example:
7566:  Suppose a point contains one field with three components, and for which the unconstrained indices are {10, 11, 12}.
7567:  When the middle component is constrained, we get the array {10, -12, 12} for (islocal=TRUE, setBC=FALSE).
7568:  Note that -12 is the involution of 11, so the user can involute negative indices to recover local indices.
7569:  The global vector does not store constrained dofs, so when this function returns global indices, say {110, -112, 111}, the value of -112 is an arbitrary flag that should not be interpreted beyond its sign.

7571:  Level: developer
7572: */
7573: PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection section, PetscBool islocal, PetscInt point, PetscInt off, PetscInt foffs[], PetscBool setBC, const PetscInt ***perms, PetscInt permsoff, const PetscInt indperm[], PetscInt indices[])
7574: {
7575:   PetscInt numFields, foff, f;

7577:   PetscFunctionBegin;
7578:   PetscCheck(islocal || !setBC, PetscObjectComm((PetscObject)section), PETSC_ERR_ARG_INCOMP, "setBC incompatible with global indices; use a local section or disable setBC");
7579:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7580:   for (f = 0, foff = 0; f < numFields; ++f) {
7581:     PetscInt        fdof, cfdof;
7582:     const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
7583:     PetscInt        cind = 0, b;
7584:     const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;

7586:     PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
7587:     PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &cfdof));
7588:     if (!cfdof || setBC) {
7589:       for (b = 0; b < fdof; ++b) {
7590:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7591:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7593:         indices[ind] = off + foff + b;
7594:       }
7595:     } else {
7596:       PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
7597:       for (b = 0; b < fdof; ++b) {
7598:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7599:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7601:         if ((cind < cfdof) && (b == fcdofs[cind])) {
7602:           indices[ind] = -(off + foff + b + 1);
7603:           ++cind;
7604:         } else {
7605:           indices[ind] = off + foff + b - (islocal ? 0 : cind);
7606:         }
7607:       }
7608:     }
7609:     foff += (setBC || islocal ? fdof : (fdof - cfdof));
7610:     foffs[f] += fdof;
7611:   }
7612:   PetscFunctionReturn(PETSC_SUCCESS);
7613: }

7615: /*
7616:   This version believes the globalSection offsets for each field, rather than just the point offset

7618:  . foffs - The offset into 'indices' for each field, since it is segregated by field

7620:  Notes:
7621:  The semantics of this function relate to that of setBC=FALSE in DMPlexGetIndicesPointFields_Internal.
7622:  Since this function uses global indices, setBC=TRUE would be invalid, so no such argument exists.
7623: */
7624: static PetscErrorCode DMPlexGetIndicesPointFieldsSplit_Internal(PetscSection section, PetscSection globalSection, PetscInt point, PetscInt foffs[], const PetscInt ***perms, PetscInt permsoff, const PetscInt indperm[], PetscInt indices[])
7625: {
7626:   PetscInt numFields, foff, f;

7628:   PetscFunctionBegin;
7629:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7630:   for (f = 0; f < numFields; ++f) {
7631:     PetscInt        fdof, cfdof;
7632:     const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
7633:     PetscInt        cind = 0, b;
7634:     const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;

7636:     PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
7637:     PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &cfdof));
7638:     PetscCall(PetscSectionGetFieldOffset(globalSection, point, f, &foff));
7639:     if (!cfdof) {
7640:       for (b = 0; b < fdof; ++b) {
7641:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7642:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7644:         indices[ind] = foff + b;
7645:       }
7646:     } else {
7647:       PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
7648:       for (b = 0; b < fdof; ++b) {
7649:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7650:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7652:         if ((cind < cfdof) && (b == fcdofs[cind])) {
7653:           indices[ind] = -(foff + b + 1);
7654:           ++cind;
7655:         } else {
7656:           indices[ind] = foff + b - cind;
7657:         }
7658:       }
7659:     }
7660:     foffs[f] += fdof;
7661:   }
7662:   PetscFunctionReturn(PETSC_SUCCESS);
7663: }

7665: static PetscErrorCode DMPlexAnchorsGetSubMatIndices(PetscInt nPoints, const PetscInt pnts[], PetscSection section, PetscSection cSec, PetscInt tmpIndices[], PetscInt fieldOffsets[], PetscInt indices[], const PetscInt ***perms)
7666: {
7667:   PetscInt numFields, sStart, sEnd, cStart, cEnd;

7669:   PetscFunctionBegin;
7670:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7671:   PetscCall(PetscSectionGetChart(section, &sStart, &sEnd));
7672:   PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
7673:   for (PetscInt p = 0; p < nPoints; p++) {
7674:     PetscInt     b       = pnts[2 * p];
7675:     PetscInt     bSecDof = 0, bOff;
7676:     PetscInt     cSecDof = 0;
7677:     PetscSection indices_section;

7679:     if (b >= sStart && b < sEnd) PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7680:     if (!bSecDof) continue;
7681:     if (b >= cStart && b < cEnd) PetscCall(PetscSectionGetDof(cSec, b, &cSecDof));
7682:     indices_section = cSecDof > 0 ? cSec : section;
7683:     if (numFields) {
7684:       PetscInt fStart[32], fEnd[32];

7686:       fStart[0] = 0;
7687:       fEnd[0]   = 0;
7688:       for (PetscInt f = 0; f < numFields; f++) {
7689:         PetscInt fDof = 0;

7691:         PetscCall(PetscSectionGetFieldDof(indices_section, b, f, &fDof));
7692:         fStart[f + 1] = fStart[f] + fDof;
7693:         fEnd[f + 1]   = fStart[f + 1];
7694:       }
7695:       PetscCall(PetscSectionGetOffset(indices_section, b, &bOff));
7696:       // only apply permutations on one side
7697:       PetscCall(DMPlexGetIndicesPointFields_Internal(indices_section, PETSC_TRUE, b, bOff, fEnd, PETSC_TRUE, perms, perms ? p : -1, NULL, tmpIndices));
7698:       for (PetscInt f = 0; f < numFields; f++) {
7699:         for (PetscInt i = fStart[f]; i < fEnd[f]; i++) { indices[fieldOffsets[f]++] = (cSecDof > 0) ? tmpIndices[i] : -(tmpIndices[i] + 1); }
7700:       }
7701:     } else {
7702:       PetscInt bEnd = 0;

7704:       PetscCall(PetscSectionGetOffset(indices_section, b, &bOff));
7705:       PetscCall(DMPlexGetIndicesPoint_Internal(indices_section, PETSC_TRUE, b, bOff, &bEnd, PETSC_TRUE, (perms && perms[0]) ? perms[0][p] : NULL, NULL, tmpIndices));

7707:       for (PetscInt i = 0; i < bEnd; i++) indices[fieldOffsets[0]++] = (cSecDof > 0) ? tmpIndices[i] : -(tmpIndices[i] + 1);
7708:     }
7709:   }
7710:   PetscFunctionReturn(PETSC_SUCCESS);
7711: }

7713: PETSC_INTERN PetscErrorCode DMPlexAnchorsGetSubMatModification(DM dm, PetscSection section, PetscInt numPoints, PetscInt numIndices, const PetscInt points[], const PetscInt ***perms, PetscInt *outNumPoints, PetscInt *outNumIndices, PetscInt *outPoints[], PetscInt offsets[], PetscScalar *outMat[])
7714: {
7715:   Mat             cMat;
7716:   PetscSection    aSec, cSec;
7717:   IS              aIS;
7718:   PetscInt        aStart = -1, aEnd = -1;
7719:   PetscInt        sStart = -1, sEnd = -1;
7720:   PetscInt        cStart = -1, cEnd = -1;
7721:   const PetscInt *anchors;
7722:   PetscInt        numFields, p;
7723:   PetscInt        newNumPoints = 0, newNumIndices = 0;
7724:   PetscInt       *newPoints, *indices, *newIndices, *tmpIndices, *tmpNewIndices;
7725:   PetscInt        oldOffsets[32];
7726:   PetscInt        newOffsets[32];
7727:   PetscInt        oldOffsetsCopy[32];
7728:   PetscInt        newOffsetsCopy[32];
7729:   PetscScalar    *modMat         = NULL;
7730:   PetscBool       anyConstrained = PETSC_FALSE;

7732:   PetscFunctionBegin;
7735:   PetscCall(PetscSectionGetNumFields(section, &numFields));

7737:   PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
7738:   /* if there are point-to-point constraints */
7739:   if (aSec) {
7740:     PetscCall(PetscArrayzero(newOffsets, 32));
7741:     PetscCall(PetscArrayzero(oldOffsets, 32));
7742:     PetscCall(ISGetIndices(aIS, &anchors));
7743:     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
7744:     PetscCall(PetscSectionGetChart(section, &sStart, &sEnd));
7745:     /* figure out how many points are going to be in the new element matrix
7746:      * (we allow double counting, because it's all just going to be summed
7747:      * into the global matrix anyway) */
7748:     for (p = 0; p < 2 * numPoints; p += 2) {
7749:       PetscInt b    = points[p];
7750:       PetscInt bDof = 0, bSecDof = 0;

7752:       if (b >= sStart && b < sEnd) PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7753:       if (!bSecDof) continue;

7755:       for (PetscInt f = 0; f < numFields; f++) {
7756:         PetscInt fDof = 0;

7758:         PetscCall(PetscSectionGetFieldDof(section, b, f, &fDof));
7759:         oldOffsets[f + 1] += fDof;
7760:       }
7761:       if (b >= aStart && b < aEnd) PetscCall(PetscSectionGetDof(aSec, b, &bDof));
7762:       if (bDof) {
7763:         /* this point is constrained */
7764:         /* it is going to be replaced by its anchors */
7765:         PetscInt bOff, q;

7767:         PetscCall(PetscSectionGetOffset(aSec, b, &bOff));
7768:         for (q = 0; q < bDof; q++) {
7769:           PetscInt a    = anchors[bOff + q];
7770:           PetscInt aDof = 0;

7772:           if (a >= sStart && a < sEnd) PetscCall(PetscSectionGetDof(section, a, &aDof));
7773:           if (aDof) {
7774:             anyConstrained = PETSC_TRUE;
7775:             newNumPoints += 1;
7776:           }
7777:           newNumIndices += aDof;
7778:           for (PetscInt f = 0; f < numFields; ++f) {
7779:             PetscInt fDof = 0;

7781:             if (a >= sStart && a < sEnd) PetscCall(PetscSectionGetFieldDof(section, a, f, &fDof));
7782:             newOffsets[f + 1] += fDof;
7783:           }
7784:         }
7785:       } else {
7786:         /* this point is not constrained */
7787:         newNumPoints++;
7788:         newNumIndices += bSecDof;
7789:         for (PetscInt f = 0; f < numFields; ++f) {
7790:           PetscInt fDof;

7792:           PetscCall(PetscSectionGetFieldDof(section, b, f, &fDof));
7793:           newOffsets[f + 1] += fDof;
7794:         }
7795:       }
7796:     }
7797:   }
7798:   if (!anyConstrained) {
7799:     if (outNumPoints) *outNumPoints = 0;
7800:     if (outNumIndices) *outNumIndices = 0;
7801:     if (outPoints) *outPoints = NULL;
7802:     if (outMat) *outMat = NULL;
7803:     if (aSec) PetscCall(ISRestoreIndices(aIS, &anchors));
7804:     PetscFunctionReturn(PETSC_SUCCESS);
7805:   }

7807:   if (outNumPoints) *outNumPoints = newNumPoints;
7808:   if (outNumIndices) *outNumIndices = newNumIndices;

7810:   for (PetscInt f = 0; f < numFields; ++f) newOffsets[f + 1] += newOffsets[f];
7811:   for (PetscInt f = 0; f < numFields; ++f) oldOffsets[f + 1] += oldOffsets[f];

7813:   if (!outPoints && !outMat) {
7814:     if (offsets) {
7815:       for (PetscInt f = 0; f <= numFields; f++) offsets[f] = newOffsets[f];
7816:     }
7817:     if (aSec) PetscCall(ISRestoreIndices(aIS, &anchors));
7818:     PetscFunctionReturn(PETSC_SUCCESS);
7819:   }

7821:   PetscCheck(!numFields || newOffsets[numFields] == newNumIndices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid size for closure %" PetscInt_FMT " should be %" PetscInt_FMT, newOffsets[numFields], newNumIndices);
7822:   PetscCheck(!numFields || oldOffsets[numFields] == numIndices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid size for closure %" PetscInt_FMT " should be %" PetscInt_FMT, oldOffsets[numFields], numIndices);

7824:   PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
7825:   PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));

7827:   /* output arrays */
7828:   PetscCall(DMGetWorkArray(dm, 2 * newNumPoints, MPIU_INT, &newPoints));
7829:   PetscCall(PetscArrayzero(newPoints, 2 * newNumPoints));

7831:   // get the new Points
7832:   for (PetscInt p = 0, newP = 0; p < numPoints; p++) {
7833:     PetscInt b    = points[2 * p];
7834:     PetscInt bDof = 0, bSecDof = 0, bOff;

7836:     if (b >= sStart && b < sEnd) PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7837:     if (!bSecDof) continue;
7838:     if (b >= aStart && b < aEnd) PetscCall(PetscSectionGetDof(aSec, b, &bDof));
7839:     if (bDof) {
7840:       PetscCall(PetscSectionGetOffset(aSec, b, &bOff));
7841:       for (PetscInt q = 0; q < bDof; q++) {
7842:         PetscInt a = anchors[bOff + q], aDof = 0;

7844:         if (a >= sStart && a < sEnd) PetscCall(PetscSectionGetDof(section, a, &aDof));
7845:         if (aDof) {
7846:           newPoints[2 * newP]     = a;
7847:           newPoints[2 * newP + 1] = 0; // orientations are accounted for in constructing the matrix, newly added points are in default orientation
7848:           newP++;
7849:         }
7850:       }
7851:     } else {
7852:       newPoints[2 * newP]     = b;
7853:       newPoints[2 * newP + 1] = points[2 * p + 1];
7854:       newP++;
7855:     }
7856:   }

7858:   if (outMat) {
7859:     PetscScalar *tmpMat;
7860:     PetscCall(PetscArraycpy(oldOffsetsCopy, oldOffsets, 32));
7861:     PetscCall(PetscArraycpy(newOffsetsCopy, newOffsets, 32));

7863:     PetscCall(DMGetWorkArray(dm, numIndices, MPIU_INT, &indices));
7864:     PetscCall(DMGetWorkArray(dm, numIndices, MPIU_INT, &tmpIndices));
7865:     PetscCall(DMGetWorkArray(dm, newNumIndices, MPIU_INT, &newIndices));
7866:     PetscCall(DMGetWorkArray(dm, newNumIndices, MPIU_INT, &tmpNewIndices));

7868:     for (PetscInt i = 0; i < numIndices; i++) indices[i] = -1;
7869:     for (PetscInt i = 0; i < newNumIndices; i++) newIndices[i] = -1;

7871:     PetscCall(DMPlexAnchorsGetSubMatIndices(numPoints, points, section, cSec, tmpIndices, oldOffsetsCopy, indices, perms));
7872:     PetscCall(DMPlexAnchorsGetSubMatIndices(newNumPoints, newPoints, section, section, tmpNewIndices, newOffsetsCopy, newIndices, NULL));

7874:     PetscCall(DMGetWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &modMat));
7875:     PetscCall(DMGetWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &tmpMat));
7876:     PetscCall(PetscArrayzero(modMat, newNumIndices * numIndices));
7877:     // for each field, insert the anchor modification into modMat
7878:     for (PetscInt f = 0; f < PetscMax(1, numFields); f++) {
7879:       PetscInt fStart    = oldOffsets[f];
7880:       PetscInt fNewStart = newOffsets[f];
7881:       for (PetscInt p = 0, newP = 0, o = fStart, oNew = fNewStart; p < numPoints; p++) {
7882:         PetscInt b    = points[2 * p];
7883:         PetscInt bDof = 0, bSecDof = 0, bOff;

7885:         if (b >= sStart && b < sEnd) {
7886:           if (numFields) {
7887:             PetscCall(PetscSectionGetFieldDof(section, b, f, &bSecDof));
7888:           } else {
7889:             PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7890:           }
7891:         }
7892:         if (!bSecDof) continue;
7893:         if (b >= aStart && b < aEnd) PetscCall(PetscSectionGetDof(aSec, b, &bDof));
7894:         if (bDof) {
7895:           PetscCall(PetscSectionGetOffset(aSec, b, &bOff));
7896:           for (PetscInt q = 0; q < bDof; q++, newP++) {
7897:             PetscInt a = anchors[bOff + q], aDof = 0;

7899:             if (a >= sStart && a < sEnd) {
7900:               if (numFields) {
7901:                 PetscCall(PetscSectionGetFieldDof(section, a, f, &aDof));
7902:               } else {
7903:                 PetscCall(PetscSectionGetDof(section, a, &aDof));
7904:               }
7905:             }
7906:             if (aDof) {
7907:               PetscCall(MatGetValues(cMat, bSecDof, &indices[o], aDof, &newIndices[oNew], tmpMat));
7908:               for (Pe