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_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 DMPlexGetFieldType_Internal(DM dm, PetscSection section, PetscInt field, PetscInt *sStart, PetscInt *sEnd, PetscViewerVTKFieldType *ft)
125: {
126:   PetscInt cdim, pStart, pEnd, vStart, vEnd, cStart, cEnd;
127:   PetscInt vcdof[2] = {0, 0}, globalvcdof[2];

129:   PetscFunctionBegin;
130:   *ft = PETSC_VTK_INVALID;
131:   PetscCall(DMGetCoordinateDim(dm, &cdim));
132:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
133:   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
134:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
135:   if (field >= 0) {
136:     if ((vStart >= pStart) && (vStart < pEnd)) PetscCall(PetscSectionGetFieldDof(section, vStart, field, &vcdof[0]));
137:     if ((cStart >= pStart) && (cStart < pEnd)) PetscCall(PetscSectionGetFieldDof(section, cStart, field, &vcdof[1]));
138:   } else {
139:     if ((vStart >= pStart) && (vStart < pEnd)) PetscCall(PetscSectionGetDof(section, vStart, &vcdof[0]));
140:     if ((cStart >= pStart) && (cStart < pEnd)) PetscCall(PetscSectionGetDof(section, cStart, &vcdof[1]));
141:   }
142:   PetscCallMPI(MPIU_Allreduce(vcdof, globalvcdof, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
143:   if (globalvcdof[0]) {
144:     *sStart = vStart;
145:     *sEnd   = vEnd;
146:     if (globalvcdof[0] == cdim) *ft = PETSC_VTK_POINT_VECTOR_FIELD;
147:     else *ft = PETSC_VTK_POINT_FIELD;
148:   } else if (globalvcdof[1]) {
149:     *sStart = cStart;
150:     *sEnd   = cEnd;
151:     if (globalvcdof[1] == cdim) *ft = PETSC_VTK_CELL_VECTOR_FIELD;
152:     else *ft = PETSC_VTK_CELL_FIELD;
153:   } else {
154:     if (field >= 0) {
155:       const char *fieldname;

157:       PetscCall(PetscSectionGetFieldName(section, field, &fieldname));
158:       PetscCall(PetscInfo((PetscObject)dm, "Could not classify VTK output type of section field %" PetscInt_FMT " \"%s\"\n", field, fieldname));
159:     } else {
160:       PetscCall(PetscInfo((PetscObject)dm, "Could not classify VTK output type of section\n"));
161:     }
162:   }
163:   PetscFunctionReturn(PETSC_SUCCESS);
164: }

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

169:   Collective

171:   Input Parameters:
172: + dm     - The `DMPLEX` object
173: . n      - The number of vectors
174: . u      - The array of local vectors
175: - viewer - The `PetscViewer`

177:   Level: advanced

179: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `VecViewFromOptions()`, `VecView()`
180: @*/
181: PetscErrorCode DMPlexVecView1D(DM dm, PetscInt n, Vec u[], PetscViewer viewer)
182: {
183:   PetscDS            ds;
184:   PetscDraw          draw = NULL;
185:   PetscDrawLG        lg;
186:   Vec                coordinates;
187:   const PetscScalar *coords, **sol;
188:   PetscReal         *vals;
189:   PetscInt          *Nc;
190:   PetscInt           Nf, f, c, Nl, l, i, vStart, vEnd, v;
191:   char             **names;

193:   PetscFunctionBegin;
194:   PetscCall(DMGetDS(dm, &ds));
195:   PetscCall(PetscDSGetNumFields(ds, &Nf));
196:   PetscCall(PetscDSGetTotalComponents(ds, &Nl));
197:   PetscCall(PetscDSGetComponents(ds, &Nc));

199:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
200:   if (!draw) PetscFunctionReturn(PETSC_SUCCESS);
201:   PetscCall(PetscDrawLGCreate(draw, n * Nl, &lg));

203:   PetscCall(PetscMalloc3(n, &sol, n * Nl, &names, n * Nl, &vals));
204:   for (i = 0, l = 0; i < n; ++i) {
205:     const char *vname;

207:     PetscCall(PetscObjectGetName((PetscObject)u[i], &vname));
208:     for (f = 0; f < Nf; ++f) {
209:       PetscObject disc;
210:       const char *fname;
211:       char        tmpname[PETSC_MAX_PATH_LEN];

213:       PetscCall(PetscDSGetDiscretization(ds, f, &disc));
214:       /* TODO Create names for components */
215:       for (c = 0; c < Nc[f]; ++c, ++l) {
216:         PetscCall(PetscObjectGetName(disc, &fname));
217:         PetscCall(PetscStrncpy(tmpname, vname, sizeof(tmpname)));
218:         PetscCall(PetscStrlcat(tmpname, ":", sizeof(tmpname)));
219:         PetscCall(PetscStrlcat(tmpname, fname, sizeof(tmpname)));
220:         PetscCall(PetscStrallocpy(tmpname, &names[l]));
221:       }
222:     }
223:   }
224:   PetscCall(PetscDrawLGSetLegend(lg, (const char *const *)names));
225:   /* Just add P_1 support for now */
226:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
227:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
228:   PetscCall(VecGetArrayRead(coordinates, &coords));
229:   for (i = 0; i < n; ++i) PetscCall(VecGetArrayRead(u[i], &sol[i]));
230:   for (v = vStart; v < vEnd; ++v) {
231:     PetscScalar *x, *svals;

233:     PetscCall(DMPlexPointLocalRead(dm, v, coords, &x));
234:     for (i = 0; i < n; ++i) {
235:       PetscCall(DMPlexPointLocalRead(dm, v, sol[i], &svals));
236:       for (l = 0; l < Nl; ++l) vals[i * Nl + l] = PetscRealPart(svals[l]);
237:     }
238:     PetscCall(PetscDrawLGAddCommonPoint(lg, PetscRealPart(x[0]), vals));
239:   }
240:   PetscCall(VecRestoreArrayRead(coordinates, &coords));
241:   for (i = 0; i < n; ++i) PetscCall(VecRestoreArrayRead(u[i], &sol[i]));
242:   for (l = 0; l < n * Nl; ++l) PetscCall(PetscFree(names[l]));
243:   PetscCall(PetscFree3(sol, names, vals));

245:   PetscCall(PetscDrawLGDraw(lg));
246:   PetscCall(PetscDrawLGDestroy(&lg));
247:   PetscFunctionReturn(PETSC_SUCCESS);
248: }

250: static PetscErrorCode VecView_Plex_Local_Draw_1D(Vec u, PetscViewer viewer)
251: {
252:   DM dm;

254:   PetscFunctionBegin;
255:   PetscCall(VecGetDM(u, &dm));
256:   PetscCall(DMPlexVecView1D(dm, 1, &u, viewer));
257:   PetscFunctionReturn(PETSC_SUCCESS);
258: }

260: static PetscErrorCode VecView_Plex_Local_Draw_2D(Vec v, PetscViewer viewer)
261: {
262:   DM                 dm;
263:   PetscSection       s;
264:   PetscDraw          draw, popup;
265:   DM                 cdm;
266:   PetscSection       coordSection;
267:   Vec                coordinates;
268:   const PetscScalar *array;
269:   PetscReal          lbound[3], ubound[3];
270:   PetscReal          vbound[2], time;
271:   PetscBool          flg;
272:   PetscInt           dim, Nf, f, Nc, comp, vStart, vEnd, cStart, cEnd, c, N, level, step, w = 0;
273:   const char        *name;
274:   char               title[PETSC_MAX_PATH_LEN];

276:   PetscFunctionBegin;
277:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
278:   PetscCall(VecGetDM(v, &dm));
279:   PetscCall(DMGetCoordinateDim(dm, &dim));
280:   PetscCall(DMGetLocalSection(dm, &s));
281:   PetscCall(PetscSectionGetNumFields(s, &Nf));
282:   PetscCall(DMGetCoarsenLevel(dm, &level));
283:   PetscCall(DMGetCoordinateDM(dm, &cdm));
284:   PetscCall(DMGetLocalSection(cdm, &coordSection));
285:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
286:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
287:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));

289:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
290:   PetscCall(DMGetOutputSequenceNumber(dm, &step, &time));

292:   PetscCall(VecGetLocalSize(coordinates, &N));
293:   PetscCall(DMGetBoundingBox(dm, lbound, ubound));
294:   PetscCall(PetscDrawClear(draw));

296:   /* Could implement something like DMDASelectFields() */
297:   for (f = 0; f < Nf; ++f) {
298:     DM          fdm = dm;
299:     Vec         fv  = v;
300:     IS          fis;
301:     char        prefix[PETSC_MAX_PATH_LEN];
302:     const char *fname;

304:     PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
305:     PetscCall(PetscSectionGetFieldName(s, f, &fname));

307:     if (v->hdr.prefix) PetscCall(PetscStrncpy(prefix, v->hdr.prefix, sizeof(prefix)));
308:     else prefix[0] = '\0';
309:     if (Nf > 1) {
310:       PetscCall(DMCreateSubDM(dm, 1, &f, &fis, &fdm));
311:       PetscCall(VecGetSubVector(v, fis, &fv));
312:       PetscCall(PetscStrlcat(prefix, fname, sizeof(prefix)));
313:       PetscCall(PetscStrlcat(prefix, "_", sizeof(prefix)));
314:     }
315:     for (comp = 0; comp < Nc; ++comp, ++w) {
316:       PetscInt nmax = 2;

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

323:       /* TODO Get max and min only for this component */
324:       PetscCall(PetscOptionsGetRealArray(NULL, prefix, "-vec_view_bounds", vbound, &nmax, &flg));
325:       if (!flg) {
326:         PetscCall(VecMin(fv, NULL, &vbound[0]));
327:         PetscCall(VecMax(fv, NULL, &vbound[1]));
328:         if (vbound[1] <= vbound[0]) vbound[1] = vbound[0] + 1.0;
329:       }

331:       PetscCall(PetscDrawGetPopup(draw, &popup));
332:       PetscCall(PetscDrawScalePopup(popup, vbound[0], vbound[1]));
333:       PetscCall(PetscDrawSetCoordinates(draw, lbound[0], lbound[1], ubound[0], ubound[1]));
334:       PetscCall(VecGetArrayRead(fv, &array));
335:       for (c = cStart; c < cEnd; ++c) {
336:         DMPolytopeType     ct;
337:         PetscScalar       *coords = NULL, *a = NULL;
338:         const PetscScalar *coords_arr;
339:         PetscBool          isDG;
340:         PetscInt           numCoords;
341:         int                color[4] = {-1, -1, -1, -1};

343:         PetscCall(DMPlexGetCellType(dm, c, &ct));
344:         PetscCall(DMPlexPointLocalRead(fdm, c, array, &a));
345:         if (a) {
346:           color[0] = PetscDrawRealToColor(PetscRealPart(a[comp]), vbound[0], vbound[1]);
347:           color[1] = color[2] = color[3] = color[0];
348:         } else {
349:           PetscScalar *vals = NULL;
350:           PetscInt     numVals, va;

352:           PetscCall(DMPlexVecGetClosure(fdm, NULL, fv, c, &numVals, &vals));
353:           if (!numVals) {
354:             PetscCall(DMPlexVecRestoreClosure(fdm, NULL, fv, c, &numVals, &vals));
355:             continue;
356:           }
357:           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);
358:           switch (numVals / Nc) {
359:           case 1: /* P1 Clamped Segment Prism */
360:           case 2: /* P1 Segment Prism, P2 Clamped Segment Prism */
361:             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]);
362:             for (va = 0; va < numVals / Nc; ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va * Nc + comp]), vbound[0], vbound[1]);
363:             break;
364:           case 3: /* P1 Triangle */
365:           case 4: /* P1 Quadrangle */
366:             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]);
367:             for (va = 0; va < numVals / Nc; ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va * Nc + comp]), vbound[0], vbound[1]);
368:             break;
369:           case 6: /* P2 Triangle */
370:           case 8: /* P2 Quadrangle */
371:             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]);
372:             for (va = 0; va < numVals / (Nc * 2); ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va * Nc + comp + numVals / (Nc * 2)]), vbound[0], vbound[1]);
373:             break;
374:           default:
375:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of values for cell closure %" PetscInt_FMT " cannot be handled", numVals / Nc);
376:           }
377:           PetscCall(DMPlexVecRestoreClosure(fdm, NULL, fv, c, &numVals, &vals));
378:         }
379:         PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
380:         switch (numCoords) {
381:         case 6:
382:         case 12: /* Localized triangle */
383:           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]));
384:           break;
385:         case 8:
386:         case 16: /* Localized quadrilateral */
387:           if (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) {
388:             PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscMax(color[0], color[1])));
389:           } else {
390:             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]));
391:             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]));
392:           }
393:           break;
394:         default:
395:           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells with %" PetscInt_FMT " coordinates", numCoords);
396:         }
397:         PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
398:       }
399:       PetscCall(VecRestoreArrayRead(fv, &array));
400:       PetscCall(PetscDrawFlush(draw));
401:       PetscCall(PetscDrawPause(draw));
402:       PetscCall(PetscDrawSave(draw));
403:     }
404:     if (Nf > 1) {
405:       PetscCall(VecRestoreSubVector(v, fis, &fv));
406:       PetscCall(ISDestroy(&fis));
407:       PetscCall(DMDestroy(&fdm));
408:     }
409:   }
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }

413: static PetscErrorCode VecView_Plex_Local_Draw(Vec v, PetscViewer viewer)
414: {
415:   DM        dm;
416:   PetscDraw draw;
417:   PetscInt  dim;
418:   PetscBool isnull;

420:   PetscFunctionBegin;
421:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
422:   PetscCall(PetscDrawIsNull(draw, &isnull));
423:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

425:   PetscCall(VecGetDM(v, &dm));
426:   PetscCall(DMGetCoordinateDim(dm, &dim));
427:   switch (dim) {
428:   case 1:
429:     PetscCall(VecView_Plex_Local_Draw_1D(v, viewer));
430:     break;
431:   case 2:
432:     PetscCall(VecView_Plex_Local_Draw_2D(v, viewer));
433:     break;
434:   default:
435:     SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot draw meshes of dimension %" PetscInt_FMT ". Try PETSCVIEWERGLVIS", dim);
436:   }
437:   PetscFunctionReturn(PETSC_SUCCESS);
438: }

440: static PetscErrorCode VecView_Plex_Local_VTK(Vec v, PetscViewer viewer)
441: {
442:   DM                      dm;
443:   Vec                     locv;
444:   const char             *name;
445:   PetscSection            section;
446:   PetscInt                pStart, pEnd;
447:   PetscInt                numFields;
448:   PetscViewerVTKFieldType ft;

450:   PetscFunctionBegin;
451:   PetscCall(VecGetDM(v, &dm));
452:   PetscCall(DMCreateLocalVector(dm, &locv)); /* VTK viewer requires exclusive ownership of the vector */
453:   PetscCall(PetscObjectGetName((PetscObject)v, &name));
454:   PetscCall(PetscObjectSetName((PetscObject)locv, name));
455:   PetscCall(VecCopy(v, locv));
456:   PetscCall(DMGetLocalSection(dm, &section));
457:   PetscCall(PetscSectionGetNumFields(section, &numFields));
458:   if (!numFields) {
459:     PetscCall(DMPlexGetFieldType_Internal(dm, section, PETSC_DETERMINE, &pStart, &pEnd, &ft));
460:     PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)dm, DMPlexVTKWriteAll, PETSC_DEFAULT, ft, PETSC_TRUE, (PetscObject)locv));
461:   } else {
462:     PetscInt f;

464:     for (f = 0; f < numFields; f++) {
465:       PetscCall(DMPlexGetFieldType_Internal(dm, section, f, &pStart, &pEnd, &ft));
466:       if (ft == PETSC_VTK_INVALID) continue;
467:       PetscCall(PetscObjectReference((PetscObject)locv));
468:       PetscCall(PetscViewerVTKAddField(viewer, (PetscObject)dm, DMPlexVTKWriteAll, f, ft, PETSC_TRUE, (PetscObject)locv));
469:     }
470:     PetscCall(VecDestroy(&locv));
471:   }
472:   PetscFunctionReturn(PETSC_SUCCESS);
473: }

475: PetscErrorCode VecView_Plex_Local(Vec v, PetscViewer viewer)
476: {
477:   DM        dm;
478:   PetscBool isvtk, ishdf5, isdraw, isglvis, iscgns;

480:   PetscFunctionBegin;
481:   PetscCall(VecGetDM(v, &dm));
482:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
483:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
484:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
485:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
486:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
487:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
488:   if (isvtk || ishdf5 || isdraw || isglvis || iscgns) {
489:     PetscInt    i, numFields;
490:     PetscObject fe;
491:     PetscBool   fem  = PETSC_FALSE;
492:     Vec         locv = v;
493:     const char *name;
494:     PetscInt    step;
495:     PetscReal   time;

497:     PetscCall(DMGetNumFields(dm, &numFields));
498:     for (i = 0; i < numFields; i++) {
499:       PetscCall(DMGetField(dm, i, NULL, &fe));
500:       if (fe->classid == PETSCFE_CLASSID) {
501:         fem = PETSC_TRUE;
502:         break;
503:       }
504:     }
505:     if (fem) {
506:       PetscObject isZero;

508:       PetscCall(DMGetLocalVector(dm, &locv));
509:       PetscCall(PetscObjectGetName((PetscObject)v, &name));
510:       PetscCall(PetscObjectSetName((PetscObject)locv, name));
511:       PetscCall(PetscObjectQuery((PetscObject)v, "__Vec_bc_zero__", &isZero));
512:       PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", isZero));
513:       PetscCall(VecCopy(v, locv));
514:       PetscCall(DMGetOutputSequenceNumber(dm, NULL, &time));
515:       PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL));
516:     }
517:     if (isvtk) {
518:       PetscCall(VecView_Plex_Local_VTK(locv, viewer));
519:     } else if (ishdf5) {
520: #if defined(PETSC_HAVE_HDF5)
521:       PetscCall(VecView_Plex_Local_HDF5_Internal(locv, viewer));
522: #else
523:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
524: #endif
525:     } else if (isdraw) {
526:       PetscCall(VecView_Plex_Local_Draw(locv, viewer));
527:     } else if (isglvis) {
528:       PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
529:       PetscCall(PetscViewerGLVisSetSnapId(viewer, step));
530:       PetscCall(VecView_GLVis(locv, viewer));
531:     } else if (iscgns) {
532: #if defined(PETSC_HAVE_CGNS)
533:       PetscCall(VecView_Plex_Local_CGNS(locv, viewer));
534: #else
535:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CGNS not supported in this build.\nPlease reconfigure using --download-cgns");
536: #endif
537:     }
538:     if (fem) {
539:       PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", NULL));
540:       PetscCall(DMRestoreLocalVector(dm, &locv));
541:     }
542:   } else {
543:     PetscBool isseq;

545:     PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
546:     if (isseq) PetscCall(VecView_Seq(v, viewer));
547:     else PetscCall(VecView_MPI(v, viewer));
548:   }
549:   PetscFunctionReturn(PETSC_SUCCESS);
550: }

552: PetscErrorCode VecView_Plex(Vec v, PetscViewer viewer)
553: {
554:   DM        dm;
555:   PetscBool isvtk, ishdf5, isdraw, isglvis, isexodusii, iscgns;

557:   PetscFunctionBegin;
558:   PetscCall(VecGetDM(v, &dm));
559:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
560:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
561:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
562:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
563:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
564:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
565:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWEREXODUSII, &isexodusii));
566:   if (isvtk || isdraw || isglvis || iscgns) {
567:     Vec         locv;
568:     PetscObject isZero;
569:     const char *name;

571:     PetscCall(DMGetLocalVector(dm, &locv));
572:     PetscCall(PetscObjectGetName((PetscObject)v, &name));
573:     PetscCall(PetscObjectSetName((PetscObject)locv, name));
574:     PetscCall(DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv));
575:     PetscCall(DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv));
576:     PetscCall(PetscObjectQuery((PetscObject)v, "__Vec_bc_zero__", &isZero));
577:     PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", isZero));
578:     PetscCall(VecView_Plex_Local(locv, viewer));
579:     PetscCall(PetscObjectCompose((PetscObject)locv, "__Vec_bc_zero__", NULL));
580:     PetscCall(DMRestoreLocalVector(dm, &locv));
581:   } else if (ishdf5) {
582: #if defined(PETSC_HAVE_HDF5)
583:     PetscCall(VecView_Plex_HDF5_Internal(v, viewer));
584: #else
585:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
586: #endif
587:   } else if (isexodusii) {
588: #if defined(PETSC_HAVE_EXODUSII)
589:     PetscCall(VecView_PlexExodusII_Internal(v, viewer));
590: #else
591:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "ExodusII not supported in this build.\nPlease reconfigure using --download-exodusii");
592: #endif
593:   } else {
594:     PetscBool isseq;

596:     PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
597:     if (isseq) PetscCall(VecView_Seq(v, viewer));
598:     else PetscCall(VecView_MPI(v, viewer));
599:   }
600:   PetscFunctionReturn(PETSC_SUCCESS);
601: }

603: PetscErrorCode VecView_Plex_Native(Vec originalv, PetscViewer viewer)
604: {
605:   DM                dm;
606:   MPI_Comm          comm;
607:   PetscViewerFormat format;
608:   Vec               v;
609:   PetscBool         isvtk, ishdf5;

611:   PetscFunctionBegin;
612:   PetscCall(VecGetDM(originalv, &dm));
613:   PetscCall(PetscObjectGetComm((PetscObject)originalv, &comm));
614:   PetscCheck(dm, comm, PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
615:   PetscCall(PetscViewerGetFormat(viewer, &format));
616:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
617:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
618:   if (format == PETSC_VIEWER_NATIVE) {
619:     /* Natural ordering is the common case for DMDA, NATIVE means plain vector, for PLEX is the opposite */
620:     /* this need a better fix */
621:     if (dm->useNatural) {
622:       if (dm->sfNatural) {
623:         const char *vecname;
624:         PetscInt    n, nroots;

626:         PetscCall(VecGetLocalSize(originalv, &n));
627:         PetscCall(PetscSFGetGraph(dm->sfNatural, &nroots, NULL, NULL, NULL));
628:         if (n == nroots) {
629:           PetscCall(DMPlexCreateNaturalVector(dm, &v));
630:           PetscCall(DMPlexGlobalToNaturalBegin(dm, originalv, v));
631:           PetscCall(DMPlexGlobalToNaturalEnd(dm, originalv, v));
632:           PetscCall(PetscObjectGetName((PetscObject)originalv, &vecname));
633:           PetscCall(PetscObjectSetName((PetscObject)v, vecname));
634:         } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "DM global to natural SF only handles global vectors");
635:       } else SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created");
636:     } else v = originalv;
637:   } else v = originalv;

639:   if (ishdf5) {
640: #if defined(PETSC_HAVE_HDF5)
641:     PetscCall(VecView_Plex_HDF5_Native_Internal(v, viewer));
642: #else
643:     SETERRQ(comm, PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
644: #endif
645:   } else if (isvtk) {
646:     SETERRQ(comm, PETSC_ERR_SUP, "VTK format does not support viewing in natural order. Please switch to HDF5.");
647:   } else {
648:     PetscBool isseq;

650:     PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
651:     if (isseq) PetscCall(VecView_Seq(v, viewer));
652:     else PetscCall(VecView_MPI(v, viewer));
653:   }
654:   if (v != originalv) PetscCall(VecDestroy(&v));
655:   PetscFunctionReturn(PETSC_SUCCESS);
656: }

658: PetscErrorCode VecLoad_Plex_Local(Vec v, PetscViewer viewer)
659: {
660:   DM        dm;
661:   PetscBool ishdf5;

663:   PetscFunctionBegin;
664:   PetscCall(VecGetDM(v, &dm));
665:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
666:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
667:   if (ishdf5) {
668:     DM          dmBC;
669:     Vec         gv;
670:     const char *name;

672:     PetscCall(DMGetOutputDM(dm, &dmBC));
673:     PetscCall(DMGetGlobalVector(dmBC, &gv));
674:     PetscCall(PetscObjectGetName((PetscObject)v, &name));
675:     PetscCall(PetscObjectSetName((PetscObject)gv, name));
676:     PetscCall(VecLoad_Default(gv, viewer));
677:     PetscCall(DMGlobalToLocalBegin(dmBC, gv, INSERT_VALUES, v));
678:     PetscCall(DMGlobalToLocalEnd(dmBC, gv, INSERT_VALUES, v));
679:     PetscCall(DMRestoreGlobalVector(dmBC, &gv));
680:   } else PetscCall(VecLoad_Default(v, viewer));
681:   PetscFunctionReturn(PETSC_SUCCESS);
682: }

684: PetscErrorCode VecLoad_Plex(Vec v, PetscViewer viewer)
685: {
686:   DM        dm;
687:   PetscBool ishdf5, isexodusii, iscgns;

689:   PetscFunctionBegin;
690:   PetscCall(VecGetDM(v, &dm));
691:   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
692:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
693:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWEREXODUSII, &isexodusii));
694:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
695:   if (ishdf5) {
696: #if defined(PETSC_HAVE_HDF5)
697:     PetscCall(VecLoad_Plex_HDF5_Internal(v, viewer));
698: #else
699:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
700: #endif
701:   } else if (isexodusii) {
702: #if defined(PETSC_HAVE_EXODUSII)
703:     PetscCall(VecLoad_PlexExodusII_Internal(v, viewer));
704: #else
705:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "ExodusII not supported in this build.\nPlease reconfigure using --download-exodusii");
706: #endif
707:   } else if (iscgns) {
708: #if defined(PETSC_HAVE_CGNS)
709:     PetscCall(VecLoad_Plex_CGNS_Internal(v, viewer));
710: #else
711:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CGNS not supported in this build.\nPlease reconfigure using --download-cgns");
712: #endif
713:   } else PetscCall(VecLoad_Default(v, viewer));
714:   PetscFunctionReturn(PETSC_SUCCESS);
715: }

717: PetscErrorCode VecLoad_Plex_Native(Vec originalv, PetscViewer viewer)
718: {
719:   DM                dm;
720:   PetscViewerFormat format;
721:   PetscBool         ishdf5;

723:   PetscFunctionBegin;
724:   PetscCall(VecGetDM(originalv, &dm));
725:   PetscCheck(dm, PetscObjectComm((PetscObject)originalv), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
726:   PetscCall(PetscViewerGetFormat(viewer, &format));
727:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
728:   if (format == PETSC_VIEWER_NATIVE) {
729:     if (dm->useNatural) {
730:       if (dm->sfNatural) {
731:         if (ishdf5) {
732: #if defined(PETSC_HAVE_HDF5)
733:           Vec         v;
734:           const char *vecname;

736:           PetscCall(DMPlexCreateNaturalVector(dm, &v));
737:           PetscCall(PetscObjectGetName((PetscObject)originalv, &vecname));
738:           PetscCall(PetscObjectSetName((PetscObject)v, vecname));
739:           PetscCall(VecLoad_Plex_HDF5_Native_Internal(v, viewer));
740:           PetscCall(DMPlexNaturalToGlobalBegin(dm, v, originalv));
741:           PetscCall(DMPlexNaturalToGlobalEnd(dm, v, originalv));
742:           PetscCall(VecDestroy(&v));
743: #else
744:           SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
745: #endif
746:         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Reading in natural order is not supported for anything but HDF5.");
747:       }
748:     } else PetscCall(VecLoad_Default(originalv, viewer));
749:   }
750:   PetscFunctionReturn(PETSC_SUCCESS);
751: }

753: PETSC_UNUSED static PetscErrorCode DMPlexView_Ascii_Geometry(DM dm, PetscViewer viewer)
754: {
755:   PetscSection       coordSection;
756:   Vec                coordinates;
757:   DMLabel            depthLabel, celltypeLabel;
758:   const char        *name[4];
759:   const PetscScalar *a;
760:   PetscInt           dim, pStart, pEnd, cStart, cEnd, c;

762:   PetscFunctionBegin;
763:   PetscCall(DMGetDimension(dm, &dim));
764:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
765:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
766:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
767:   PetscCall(DMPlexGetCellTypeLabel(dm, &celltypeLabel));
768:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
769:   PetscCall(PetscSectionGetChart(coordSection, &pStart, &pEnd));
770:   PetscCall(VecGetArrayRead(coordinates, &a));
771:   name[0]       = "vertex";
772:   name[1]       = "edge";
773:   name[dim - 1] = "face";
774:   name[dim]     = "cell";
775:   for (c = cStart; c < cEnd; ++c) {
776:     PetscInt *closure = NULL;
777:     PetscInt  closureSize, cl, ct;

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

786:       if ((point < pStart) || (point >= pEnd)) continue;
787:       PetscCall(PetscSectionGetDof(coordSection, point, &dof));
788:       if (!dof) continue;
789:       PetscCall(DMLabelGetValue(depthLabel, point, &depth));
790:       PetscCall(PetscSectionGetOffset(coordSection, point, &off));
791:       PetscCall(PetscViewerASCIIPrintf(viewer, "%s %" PetscInt_FMT " coords:", name[depth], point));
792:       for (p = 0; p < dof / dim; ++p) {
793:         PetscCall(PetscViewerASCIIPrintf(viewer, " ("));
794:         for (d = 0; d < dim; ++d) {
795:           if (d > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
796:           PetscCall(PetscViewerASCIIPrintf(viewer, "%g", (double)PetscRealPart(a[off + p * dim + d])));
797:         }
798:         PetscCall(PetscViewerASCIIPrintf(viewer, ")"));
799:       }
800:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
801:     }
802:     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
803:     PetscCall(PetscViewerASCIIPopTab(viewer));
804:   }
805:   PetscCall(VecRestoreArrayRead(coordinates, &a));
806:   PetscFunctionReturn(PETSC_SUCCESS);
807: }

809: typedef enum {
810:   CS_CARTESIAN,
811:   CS_POLAR,
812:   CS_CYLINDRICAL,
813:   CS_SPHERICAL
814: } CoordSystem;
815: const char *CoordSystems[] = {"cartesian", "polar", "cylindrical", "spherical", "CoordSystem", "CS_", NULL};

817: static PetscErrorCode DMPlexView_Ascii_Coordinates(PetscViewer viewer, CoordSystem cs, PetscInt dim, const PetscScalar x[])
818: {
819:   PetscInt i;

821:   PetscFunctionBegin;
822:   if (dim > 3) {
823:     for (i = 0; i < dim; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %g", (double)PetscRealPart(x[i])));
824:   } else {
825:     PetscReal coords[3], trcoords[3] = {0., 0., 0.};

827:     for (i = 0; i < dim; ++i) coords[i] = PetscRealPart(x[i]);
828:     switch (cs) {
829:     case CS_CARTESIAN:
830:       for (i = 0; i < dim; ++i) trcoords[i] = coords[i];
831:       break;
832:     case CS_POLAR:
833:       PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Polar coordinates are for 2 dimension, not %" PetscInt_FMT, dim);
834:       trcoords[0] = PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1]));
835:       trcoords[1] = PetscAtan2Real(coords[1], coords[0]);
836:       break;
837:     case CS_CYLINDRICAL:
838:       PetscCheck(dim == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cylindrical coordinates are for 3 dimension, not %" PetscInt_FMT, dim);
839:       trcoords[0] = PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1]));
840:       trcoords[1] = PetscAtan2Real(coords[1], coords[0]);
841:       trcoords[2] = coords[2];
842:       break;
843:     case CS_SPHERICAL:
844:       PetscCheck(dim == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Spherical coordinates are for 3 dimension, not %" PetscInt_FMT, dim);
845:       trcoords[0] = PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1]) + PetscSqr(coords[2]));
846:       trcoords[1] = PetscAtan2Real(PetscSqrtReal(PetscSqr(coords[0]) + PetscSqr(coords[1])), coords[2]);
847:       trcoords[2] = PetscAtan2Real(coords[1], coords[0]);
848:       break;
849:     }
850:     for (i = 0; i < dim; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %g", (double)trcoords[i]));
851:   }
852:   PetscFunctionReturn(PETSC_SUCCESS);
853: }

855: static PetscErrorCode DMPlexView_Ascii(DM dm, PetscViewer viewer)
856: {
857:   DM_Plex          *mesh = (DM_Plex *)dm->data;
858:   DM                cdm, cdmCell;
859:   PetscSection      coordSection, coordSectionCell;
860:   Vec               coordinates, coordinatesCell;
861:   PetscViewerFormat format;

863:   PetscFunctionBegin;
864:   PetscCall(PetscViewerGetFormat(viewer, &format));
865:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
866:     const char *name;
867:     PetscInt    dim, cellHeight, maxConeSize, maxSupportSize;
868:     PetscInt    pStart, pEnd, p, numLabels, l;
869:     PetscMPIInt rank, size;

871:     PetscCall(DMGetCoordinateDM(dm, &cdm));
872:     PetscCall(DMGetCoordinateSection(dm, &coordSection));
873:     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
874:     PetscCall(DMGetCellCoordinateDM(dm, &cdmCell));
875:     PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell));
876:     PetscCall(DMGetCellCoordinatesLocal(dm, &coordinatesCell));
877:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
878:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
879:     PetscCall(PetscObjectGetName((PetscObject)dm, &name));
880:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
881:     PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
882:     PetscCall(DMGetDimension(dm, &dim));
883:     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
884:     if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
885:     else PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
886:     if (cellHeight) PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells are at height %" PetscInt_FMT "\n", cellHeight));
887:     PetscCall(PetscViewerASCIIPrintf(viewer, "Supports:\n"));
888:     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
889:     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max support size: %" PetscInt_FMT "\n", rank, maxSupportSize));
890:     for (p = pStart; p < pEnd; ++p) {
891:       PetscInt dof, off, s;

893:       PetscCall(PetscSectionGetDof(mesh->supportSection, p, &dof));
894:       PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
895:       for (s = off; s < off + dof; ++s) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", rank, p, mesh->supports[s]));
896:     }
897:     PetscCall(PetscViewerFlush(viewer));
898:     PetscCall(PetscViewerASCIIPrintf(viewer, "Cones:\n"));
899:     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max cone size: %" PetscInt_FMT "\n", rank, maxConeSize));
900:     for (p = pStart; p < pEnd; ++p) {
901:       PetscInt dof, off, c;

903:       PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
904:       PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
905:       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]));
906:     }
907:     PetscCall(PetscViewerFlush(viewer));
908:     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
909:     if (coordSection && coordinates) {
910:       CoordSystem        cs = CS_CARTESIAN;
911:       const PetscScalar *array, *arrayCell = NULL;
912:       PetscInt           Nf, Nc, pvStart, pvEnd, pcStart = PETSC_INT_MAX, pcEnd = PETSC_INT_MIN, pStart, pEnd, p;
913:       PetscMPIInt        rank;
914:       const char        *name;

916:       PetscCall(PetscOptionsGetEnum(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_coord_system", CoordSystems, (PetscEnum *)&cs, NULL));
917:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
918:       PetscCall(PetscSectionGetNumFields(coordSection, &Nf));
919:       PetscCheck(Nf == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coordinate section should have 1 field, not %" PetscInt_FMT, Nf);
920:       PetscCall(PetscSectionGetFieldComponents(coordSection, 0, &Nc));
921:       PetscCall(PetscSectionGetChart(coordSection, &pvStart, &pvEnd));
922:       if (coordSectionCell) PetscCall(PetscSectionGetChart(coordSectionCell, &pcStart, &pcEnd));
923:       pStart = PetscMin(pvStart, pcStart);
924:       pEnd   = PetscMax(pvEnd, pcEnd);
925:       PetscCall(PetscObjectGetName((PetscObject)coordinates, &name));
926:       PetscCall(PetscViewerASCIIPrintf(viewer, "%s with %" PetscInt_FMT " fields\n", name, Nf));
927:       PetscCall(PetscViewerASCIIPrintf(viewer, "  field 0 with %" PetscInt_FMT " components\n", Nc));
928:       if (cs != CS_CARTESIAN) PetscCall(PetscViewerASCIIPrintf(viewer, "  output coordinate system: %s\n", CoordSystems[cs]));

930:       PetscCall(VecGetArrayRead(coordinates, &array));
931:       if (coordinatesCell) PetscCall(VecGetArrayRead(coordinatesCell, &arrayCell));
932:       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
933:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank));
934:       for (p = pStart; p < pEnd; ++p) {
935:         PetscInt dof, off;

937:         if (p >= pvStart && p < pvEnd) {
938:           PetscCall(PetscSectionGetDof(coordSection, p, &dof));
939:           PetscCall(PetscSectionGetOffset(coordSection, p, &off));
940:           if (dof) {
941:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  (%4" PetscInt_FMT ") dof %2" PetscInt_FMT " offset %3" PetscInt_FMT, p, dof, off));
942:             PetscCall(DMPlexView_Ascii_Coordinates(viewer, cs, dof, &array[off]));
943:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
944:           }
945:         }
946:         if (cdmCell && p >= pcStart && p < pcEnd) {
947:           PetscCall(PetscSectionGetDof(coordSectionCell, p, &dof));
948:           PetscCall(PetscSectionGetOffset(coordSectionCell, p, &off));
949:           if (dof) {
950:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  (%4" PetscInt_FMT ") dof %2" PetscInt_FMT " offset %3" PetscInt_FMT, p, dof, off));
951:             PetscCall(DMPlexView_Ascii_Coordinates(viewer, cs, dof, &arrayCell[off]));
952:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
953:           }
954:         }
955:       }
956:       PetscCall(PetscViewerFlush(viewer));
957:       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
958:       PetscCall(VecRestoreArrayRead(coordinates, &array));
959:       if (coordinatesCell) PetscCall(VecRestoreArrayRead(coordinatesCell, &arrayCell));
960:     }
961:     PetscCall(DMGetNumLabels(dm, &numLabels));
962:     if (numLabels) PetscCall(PetscViewerASCIIPrintf(viewer, "Labels:\n"));
963:     for (l = 0; l < numLabels; ++l) {
964:       DMLabel     label;
965:       PetscBool   isdepth;
966:       const char *name;

968:       PetscCall(DMGetLabelName(dm, l, &name));
969:       PetscCall(PetscStrcmp(name, "depth", &isdepth));
970:       if (isdepth) continue;
971:       PetscCall(DMGetLabel(dm, name, &label));
972:       PetscCall(DMLabelView(label, viewer));
973:     }
974:     if (size > 1) {
975:       PetscSF sf;

977:       PetscCall(DMGetPointSF(dm, &sf));
978:       PetscCall(PetscSFView(sf, viewer));
979:     }
980:     if (mesh->periodic.face_sfs)
981:       for (PetscInt i = 0; i < mesh->periodic.num_face_sfs; i++) PetscCall(PetscSFView(mesh->periodic.face_sfs[i], viewer));
982:     PetscCall(PetscViewerFlush(viewer));
983:   } else if (format == PETSC_VIEWER_ASCII_LATEX) {
984:     const char  *name, *color;
985:     const char  *defcolors[3]  = {"gray", "orange", "green"};
986:     const char  *deflcolors[4] = {"blue", "cyan", "red", "magenta"};
987:     char         lname[PETSC_MAX_PATH_LEN];
988:     PetscReal    scale      = 2.0;
989:     PetscReal    tikzscale  = 1.0;
990:     PetscBool    useNumbers = PETSC_TRUE, drawNumbers[4], drawColors[4], useLabels, useColors, plotEdges, drawHasse = PETSC_FALSE;
991:     double       tcoords[3];
992:     PetscScalar *coords;
993:     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;
994:     PetscMPIInt  rank, size;
995:     char       **names, **colors, **lcolors;
996:     PetscBool    flg, lflg;
997:     PetscBT      wp = NULL;
998:     PetscInt     pEnd, pStart;

1000:     PetscCall(DMGetCoordinateDM(dm, &cdm));
1001:     PetscCall(DMGetCoordinateSection(dm, &coordSection));
1002:     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1003:     PetscCall(DMGetCellCoordinateDM(dm, &cdmCell));
1004:     PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell));
1005:     PetscCall(DMGetCellCoordinatesLocal(dm, &coordinatesCell));
1006:     PetscCall(DMGetDimension(dm, &dim));
1007:     PetscCall(DMPlexGetDepth(dm, &depth));
1008:     PetscCall(DMGetNumLabels(dm, &numLabels));
1009:     numLabels  = PetscMax(numLabels, 10);
1010:     numColors  = 10;
1011:     numLColors = 10;
1012:     PetscCall(PetscCalloc3(numLabels, &names, numColors, &colors, numLColors, &lcolors));
1013:     PetscCall(PetscOptionsGetReal(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_scale", &scale, NULL));
1014:     PetscCall(PetscOptionsGetReal(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_tikzscale", &tikzscale, NULL));
1015:     PetscCall(PetscOptionsGetBool(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_numbers", &useNumbers, NULL));
1016:     for (d = 0; d < 4; ++d) drawNumbers[d] = useNumbers;
1017:     for (d = 0; d < 4; ++d) drawColors[d] = PETSC_TRUE;
1018:     n = 4;
1019:     PetscCall(PetscOptionsGetBoolArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_numbers_depth", drawNumbers, &n, &flg));
1020:     PetscCheck(!flg || n == dim + 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of flags %" PetscInt_FMT " != %" PetscInt_FMT " dim+1", n, dim + 1);
1021:     n = 4;
1022:     PetscCall(PetscOptionsGetBoolArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_colors_depth", drawColors, &n, &flg));
1023:     PetscCheck(!flg || n == dim + 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of flags %" PetscInt_FMT " != %" PetscInt_FMT " dim+1", n, dim + 1);
1024:     PetscCall(PetscOptionsGetStringArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_labels", names, &numLabels, &useLabels));
1025:     if (!useLabels) numLabels = 0;
1026:     PetscCall(PetscOptionsGetStringArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_colors", colors, &numColors, &useColors));
1027:     if (!useColors) {
1028:       numColors = 3;
1029:       for (c = 0; c < numColors; ++c) PetscCall(PetscStrallocpy(defcolors[c], &colors[c]));
1030:     }
1031:     PetscCall(PetscOptionsGetStringArray(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_lcolors", lcolors, &numLColors, &useColors));
1032:     if (!useColors) {
1033:       numLColors = 4;
1034:       for (c = 0; c < numLColors; ++c) PetscCall(PetscStrallocpy(deflcolors[c], &lcolors[c]));
1035:     }
1036:     PetscCall(PetscOptionsGetString(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_label_filter", lname, sizeof(lname), &lflg));
1037:     plotEdges = (PetscBool)(depth > 1 && drawNumbers[1] && dim < 3);
1038:     PetscCall(PetscOptionsGetBool(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_edges", &plotEdges, &flg));
1039:     PetscCheck(!flg || !plotEdges || depth >= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh must be interpolated");
1040:     if (depth < dim) plotEdges = PETSC_FALSE;
1041:     PetscCall(PetscOptionsGetBool(((PetscObject)viewer)->options, ((PetscObject)viewer)->prefix, "-dm_plex_view_hasse", &drawHasse, NULL));

1043:     /* filter points with labelvalue != labeldefaultvalue */
1044:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1045:     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1046:     PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
1047:     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1048:     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
1049:     if (lflg) {
1050:       DMLabel lbl;

1052:       PetscCall(DMGetLabel(dm, lname, &lbl));
1053:       if (lbl) {
1054:         PetscInt val, defval;

1056:         PetscCall(DMLabelGetDefaultValue(lbl, &defval));
1057:         PetscCall(PetscBTCreate(pEnd - pStart, &wp));
1058:         for (c = pStart; c < pEnd; c++) {
1059:           PetscInt *closure = NULL;
1060:           PetscInt  closureSize;

1062:           PetscCall(DMLabelGetValue(lbl, c, &val));
1063:           if (val == defval) continue;

1065:           PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1066:           for (p = 0; p < closureSize * 2; p += 2) PetscCall(PetscBTSet(wp, closure[p] - pStart));
1067:           PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1068:         }
1069:       }
1070:     }

1072:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1073:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
1074:     PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1075:     PetscCall(PetscViewerASCIIPrintf(viewer, "\
1076: \\documentclass[tikz]{standalone}\n\n\
1077: \\usepackage{pgflibraryshapes}\n\
1078: \\usetikzlibrary{backgrounds}\n\
1079: \\usetikzlibrary{arrows}\n\
1080: \\begin{document}\n"));
1081:     if (size > 1) {
1082:       PetscCall(PetscViewerASCIIPrintf(viewer, "%s for process ", name));
1083:       for (p = 0; p < size; ++p) {
1084:         if (p) PetscCall(PetscViewerASCIIPrintf(viewer, (p == size - 1) ? ", and " : ", "));
1085:         PetscCall(PetscViewerASCIIPrintf(viewer, "{\\textcolor{%s}%" PetscInt_FMT "}", colors[p % numColors], p));
1086:       }
1087:       PetscCall(PetscViewerASCIIPrintf(viewer, ".\n\n\n"));
1088:     }
1089:     if (drawHasse) {
1090:       PetscInt maxStratum = PetscMax(vEnd - vStart, PetscMax(eEnd - eStart, PetscMax(fEnd - fStart, cEnd - cStart)));

1092:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vStart}{%" PetscInt_FMT "}\n", vStart));
1093:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vEnd}{%" PetscInt_FMT "}\n", vEnd - 1));
1094:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numVertices}{%" PetscInt_FMT "}\n", vEnd - vStart));
1095:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vShift}{%.2f}\n", 3 + (maxStratum - (vEnd - vStart)) / 2.));
1096:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eStart}{%" PetscInt_FMT "}\n", eStart));
1097:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eEnd}{%" PetscInt_FMT "}\n", eEnd - 1));
1098:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eShift}{%.2f}\n", 3 + (maxStratum - (eEnd - eStart)) / 2.));
1099:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numEdges}{%" PetscInt_FMT "}\n", eEnd - eStart));
1100:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\fStart}{%" PetscInt_FMT "}\n", fStart));
1101:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\fEnd}{%" PetscInt_FMT "}\n", fEnd - 1));
1102:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\fShift}{%.2f}\n", 3 + (maxStratum - (fEnd - fStart)) / 2.));
1103:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numFaces}{%" PetscInt_FMT "}\n", fEnd - fStart));
1104:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cStart}{%" PetscInt_FMT "}\n", cStart));
1105:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cEnd}{%" PetscInt_FMT "}\n", cEnd - 1));
1106:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numCells}{%" PetscInt_FMT "}\n", cEnd - cStart));
1107:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cShift}{%.2f}\n", 3 + (maxStratum - (cEnd - cStart)) / 2.));
1108:     }
1109:     PetscCall(PetscViewerASCIIPrintf(viewer, "\\begin{tikzpicture}[scale = %g,font=\\fontsize{8}{8}\\selectfont]\n", (double)tikzscale));

1111:     /* Plot vertices */
1112:     PetscCall(VecGetArray(coordinates, &coords));
1113:     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1114:     for (v = vStart; v < vEnd; ++v) {
1115:       PetscInt  off, dof, d;
1116:       PetscBool isLabeled = PETSC_FALSE;

1118:       if (wp && !PetscBTLookup(wp, v - pStart)) continue;
1119:       PetscCall(PetscSectionGetDof(coordSection, v, &dof));
1120:       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
1121:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\path ("));
1122:       PetscCheck(dof <= 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "coordSection vertex %" PetscInt_FMT " has dof %" PetscInt_FMT " > 3", v, dof);
1123:       for (d = 0; d < dof; ++d) {
1124:         tcoords[d] = (double)(scale * PetscRealPart(coords[off + d]));
1125:         tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1126:       }
1127:       /* Rotate coordinates since PGF makes z point out of the page instead of up */
1128:       if (dim == 3) {
1129:         PetscReal tmp = tcoords[1];
1130:         tcoords[1]    = tcoords[2];
1131:         tcoords[2]    = -tmp;
1132:       }
1133:       for (d = 0; d < dof; ++d) {
1134:         if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ","));
1135:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double)tcoords[d]));
1136:       }
1137:       if (drawHasse) color = colors[0 % numColors];
1138:       else color = colors[rank % numColors];
1139:       for (l = 0; l < numLabels; ++l) {
1140:         PetscInt val;
1141:         PetscCall(DMGetLabelValue(dm, names[l], v, &val));
1142:         if (val >= 0) {
1143:           color     = lcolors[l % numLColors];
1144:           isLabeled = PETSC_TRUE;
1145:           break;
1146:         }
1147:       }
1148:       if (drawNumbers[0]) {
1149:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [draw,shape=circle,color=%s] {%" PetscInt_FMT "};\n", v, rank, color, v));
1150:       } else if (drawColors[0]) {
1151:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", v, rank, !isLabeled ? 1 : 2, color));
1152:       } else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [] {};\n", v, rank));
1153:     }
1154:     PetscCall(VecRestoreArray(coordinates, &coords));
1155:     PetscCall(PetscViewerFlush(viewer));
1156:     /* Plot edges */
1157:     if (plotEdges) {
1158:       PetscCall(VecGetArray(coordinates, &coords));
1159:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\path\n"));
1160:       for (e = eStart; e < eEnd; ++e) {
1161:         const PetscInt *cone;
1162:         PetscInt        coneSize, offA, offB, dof, d;

1164:         if (wp && !PetscBTLookup(wp, e - pStart)) continue;
1165:         PetscCall(DMPlexGetConeSize(dm, e, &coneSize));
1166:         PetscCheck(coneSize == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " cone should have two vertices, not %" PetscInt_FMT, e, coneSize);
1167:         PetscCall(DMPlexGetCone(dm, e, &cone));
1168:         PetscCall(PetscSectionGetDof(coordSection, cone[0], &dof));
1169:         PetscCall(PetscSectionGetOffset(coordSection, cone[0], &offA));
1170:         PetscCall(PetscSectionGetOffset(coordSection, cone[1], &offB));
1171:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "("));
1172:         for (d = 0; d < dof; ++d) {
1173:           tcoords[d] = (double)(0.5 * scale * PetscRealPart(coords[offA + d] + coords[offB + d]));
1174:           tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1175:         }
1176:         /* Rotate coordinates since PGF makes z point out of the page instead of up */
1177:         if (dim == 3) {
1178:           PetscReal tmp = tcoords[1];
1179:           tcoords[1]    = tcoords[2];
1180:           tcoords[2]    = -tmp;
1181:         }
1182:         for (d = 0; d < dof; ++d) {
1183:           if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ","));
1184:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double)tcoords[d]));
1185:         }
1186:         if (drawHasse) color = colors[1 % numColors];
1187:         else color = colors[rank % numColors];
1188:         for (l = 0; l < numLabels; ++l) {
1189:           PetscInt val;
1190:           PetscCall(DMGetLabelValue(dm, names[l], e, &val));
1191:           if (val >= 0) {
1192:             color = lcolors[l % numLColors];
1193:             break;
1194:           }
1195:         }
1196:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [draw,shape=circle,color=%s] {%" PetscInt_FMT "} --\n", e, rank, color, e));
1197:       }
1198:       PetscCall(VecRestoreArray(coordinates, &coords));
1199:       PetscCall(PetscViewerFlush(viewer));
1200:       PetscCall(PetscViewerASCIIPrintf(viewer, "(0,0);\n"));
1201:     }
1202:     /* Plot cells */
1203:     if (dim == 3 || !drawNumbers[1]) {
1204:       for (e = eStart; e < eEnd; ++e) {
1205:         const PetscInt *cone;

1207:         if (wp && !PetscBTLookup(wp, e - pStart)) continue;
1208:         color = colors[rank % numColors];
1209:         for (l = 0; l < numLabels; ++l) {
1210:           PetscInt val;
1211:           PetscCall(DMGetLabelValue(dm, names[l], e, &val));
1212:           if (val >= 0) {
1213:             color = lcolors[l % numLColors];
1214:             break;
1215:           }
1216:         }
1217:         PetscCall(DMPlexGetCone(dm, e, &cone));
1218:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] (%" PetscInt_FMT "_%d) -- (%" PetscInt_FMT "_%d);\n", color, cone[0], rank, cone[1], rank));
1219:       }
1220:     } else {
1221:       DMPolytopeType ct;

1223:       /* Drawing a 2D polygon */
1224:       for (c = cStart; c < cEnd; ++c) {
1225:         if (wp && !PetscBTLookup(wp, c - pStart)) continue;
1226:         PetscCall(DMPlexGetCellType(dm, c, &ct));
1227:         if (DMPolytopeTypeIsHybrid(ct)) {
1228:           const PetscInt *cone;
1229:           PetscInt        coneSize, e;

1231:           PetscCall(DMPlexGetCone(dm, c, &cone));
1232:           PetscCall(DMPlexGetConeSize(dm, c, &coneSize));
1233:           for (e = 0; e < coneSize; ++e) {
1234:             const PetscInt *econe;

1236:             PetscCall(DMPlexGetCone(dm, cone[e], &econe));
1237:             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));
1238:           }
1239:         } else {
1240:           PetscInt *closure = NULL;
1241:           PetscInt  closureSize, Nv = 0, v;

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

1247:             if ((point >= vStart) && (point < vEnd)) closure[Nv++] = point;
1248:           }
1249:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] ", colors[rank % numColors]));
1250:           for (v = 0; v <= Nv; ++v) {
1251:             const PetscInt vertex = closure[v % Nv];

1253:             if (v > 0) {
1254:               if (plotEdges) {
1255:                 const PetscInt *edge;
1256:                 PetscInt        endpoints[2], ne;

1258:                 endpoints[0] = closure[v - 1];
1259:                 endpoints[1] = vertex;
1260:                 PetscCall(DMPlexGetJoin(dm, 2, endpoints, &ne, &edge));
1261:                 PetscCheck(ne == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find edge for vertices %" PetscInt_FMT ", %" PetscInt_FMT, endpoints[0], endpoints[1]);
1262:                 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " -- (%" PetscInt_FMT "_%d) -- ", edge[0], rank));
1263:                 PetscCall(DMPlexRestoreJoin(dm, 2, endpoints, &ne, &edge));
1264:               } else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " -- "));
1265:             }
1266:             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "(%" PetscInt_FMT "_%d)", vertex, rank));
1267:           }
1268:           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ";\n"));
1269:           PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
1270:         }
1271:       }
1272:     }
1273:     for (c = cStart; c < cEnd; ++c) {
1274:       double             ccoords[3] = {0.0, 0.0, 0.0};
1275:       PetscBool          isLabeled  = PETSC_FALSE;
1276:       PetscScalar       *cellCoords = NULL;
1277:       const PetscScalar *array;
1278:       PetscInt           numCoords, cdim, d;
1279:       PetscBool          isDG;

1281:       if (wp && !PetscBTLookup(wp, c - pStart)) continue;
1282:       PetscCall(DMGetCoordinateDim(dm, &cdim));
1283:       PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &array, &cellCoords));
1284:       PetscCheck(!(numCoords % cdim), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "coordinate dim %" PetscInt_FMT " does not divide numCoords %" PetscInt_FMT, cdim, numCoords);
1285:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\\path ("));
1286:       for (p = 0; p < numCoords / cdim; ++p) {
1287:         for (d = 0; d < cdim; ++d) {
1288:           tcoords[d] = (double)(scale * PetscRealPart(cellCoords[p * cdim + d]));
1289:           tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1290:         }
1291:         /* Rotate coordinates since PGF makes z point out of the page instead of up */
1292:         if (cdim == 3) {
1293:           PetscReal tmp = tcoords[1];
1294:           tcoords[1]    = tcoords[2];
1295:           tcoords[2]    = -tmp;
1296:         }
1297:         for (d = 0; d < dim; ++d) ccoords[d] += tcoords[d];
1298:       }
1299:       for (d = 0; d < cdim; ++d) ccoords[d] /= (numCoords / cdim);
1300:       PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &array, &cellCoords));
1301:       for (d = 0; d < cdim; ++d) {
1302:         if (d > 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ","));
1303:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double)ccoords[d]));
1304:       }
1305:       if (drawHasse) color = colors[depth % numColors];
1306:       else color = colors[rank % numColors];
1307:       for (l = 0; l < numLabels; ++l) {
1308:         PetscInt val;
1309:         PetscCall(DMGetLabelValue(dm, names[l], c, &val));
1310:         if (val >= 0) {
1311:           color     = lcolors[l % numLColors];
1312:           isLabeled = PETSC_TRUE;
1313:           break;
1314:         }
1315:       }
1316:       if (drawNumbers[dim]) {
1317:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [draw,shape=circle,color=%s] {%" PetscInt_FMT "};\n", c, rank, color, c));
1318:       } else if (drawColors[dim]) {
1319:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", c, rank, !isLabeled ? 1 : 2, color));
1320:       } else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ") node(%" PetscInt_FMT "_%d) [] {};\n", c, rank));
1321:     }
1322:     if (drawHasse) {
1323:       int height = 0;

1325:       color = colors[depth % numColors];
1326:       PetscCall(PetscViewerASCIIPrintf(viewer, "%% Cells\n"));
1327:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\c in {\\cStart,...,\\cEnd}\n"));
1328:       PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1329:       PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\c_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\cShift+\\c-\\cStart,%d) {\\c};\n", rank, color, height++));
1330:       PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));

1332:       if (depth > 2) {
1333:         color = colors[1 % numColors];
1334:         PetscCall(PetscViewerASCIIPrintf(viewer, "%% Faces\n"));
1335:         PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\f in {\\fStart,...,\\fEnd}\n"));
1336:         PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1337:         PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\f_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\fShift+\\f-\\fStart,%d) {\\f};\n", rank, color, height++));
1338:         PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));
1339:       }

1341:       color = colors[1 % numColors];
1342:       PetscCall(PetscViewerASCIIPrintf(viewer, "%% Edges\n"));
1343:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\e in {\\eStart,...,\\eEnd}\n"));
1344:       PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1345:       PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\e_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\eShift+\\e-\\eStart,%d) {\\e};\n", rank, color, height++));
1346:       PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));

1348:       color = colors[0 % numColors];
1349:       PetscCall(PetscViewerASCIIPrintf(viewer, "%% Vertices\n"));
1350:       PetscCall(PetscViewerASCIIPrintf(viewer, "\\foreach \\v in {\\vStart,...,\\vEnd}\n"));
1351:       PetscCall(PetscViewerASCIIPrintf(viewer, "{\n"));
1352:       PetscCall(PetscViewerASCIIPrintf(viewer, "  \\node(\\v_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\vShift+\\v-\\vStart,%d) {\\v};\n", rank, color, height++));
1353:       PetscCall(PetscViewerASCIIPrintf(viewer, "}\n"));

1355:       for (p = pStart; p < pEnd; ++p) {
1356:         const PetscInt *cone;
1357:         PetscInt        coneSize, cp;

1359:         PetscCall(DMPlexGetCone(dm, p, &cone));
1360:         PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
1361:         for (cp = 0; cp < coneSize; ++cp) PetscCall(PetscViewerASCIIPrintf(viewer, "\\draw[->, shorten >=1pt] (%" PetscInt_FMT "_%d) -- (%" PetscInt_FMT "_%d);\n", cone[cp], rank, p, rank));
1362:       }
1363:     }
1364:     PetscCall(PetscViewerFlush(viewer));
1365:     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1366:     PetscCall(PetscViewerASCIIPrintf(viewer, "\\end{tikzpicture}\n"));
1367:     PetscCall(PetscViewerASCIIPrintf(viewer, "\\end{document}\n"));
1368:     for (l = 0; l < numLabels; ++l) PetscCall(PetscFree(names[l]));
1369:     for (c = 0; c < numColors; ++c) PetscCall(PetscFree(colors[c]));
1370:     for (c = 0; c < numLColors; ++c) PetscCall(PetscFree(lcolors[c]));
1371:     PetscCall(PetscFree3(names, colors, lcolors));
1372:     PetscCall(PetscBTDestroy(&wp));
1373:   } else if (format == PETSC_VIEWER_LOAD_BALANCE) {
1374:     Vec                    cown, acown;
1375:     VecScatter             sct;
1376:     ISLocalToGlobalMapping g2l;
1377:     IS                     gid, acis;
1378:     MPI_Comm               comm, ncomm = MPI_COMM_NULL;
1379:     MPI_Group              ggroup, ngroup;
1380:     PetscScalar           *array, nid;
1381:     const PetscInt        *idxs;
1382:     PetscInt              *idxs2, *start, *adjacency, *work;
1383:     PetscInt64             lm[3], gm[3];
1384:     PetscInt               i, c, cStart, cEnd, cum, numVertices, ect, ectn, cellHeight;
1385:     PetscMPIInt            d1, d2, rank;

1387:     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1388:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1389: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
1390:     PetscCallMPI(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, rank, MPI_INFO_NULL, &ncomm));
1391: #endif
1392:     if (ncomm != MPI_COMM_NULL) {
1393:       PetscCallMPI(MPI_Comm_group(comm, &ggroup));
1394:       PetscCallMPI(MPI_Comm_group(ncomm, &ngroup));
1395:       d1 = 0;
1396:       PetscCallMPI(MPI_Group_translate_ranks(ngroup, 1, &d1, ggroup, &d2));
1397:       nid = d2;
1398:       PetscCallMPI(MPI_Group_free(&ggroup));
1399:       PetscCallMPI(MPI_Group_free(&ngroup));
1400:       PetscCallMPI(MPI_Comm_free(&ncomm));
1401:     } else nid = 0.0;

1403:     /* Get connectivity */
1404:     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1405:     PetscCall(DMPlexCreatePartitionerGraph(dm, cellHeight, &numVertices, &start, &adjacency, &gid));

1407:     /* filter overlapped local cells */
1408:     PetscCall(DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd));
1409:     PetscCall(ISGetIndices(gid, &idxs));
1410:     PetscCall(ISGetLocalSize(gid, &cum));
1411:     PetscCall(PetscMalloc1(cum, &idxs2));
1412:     for (c = cStart, cum = 0; c < cEnd; c++) {
1413:       if (idxs[c - cStart] < 0) continue;
1414:       idxs2[cum++] = idxs[c - cStart];
1415:     }
1416:     PetscCall(ISRestoreIndices(gid, &idxs));
1417:     PetscCheck(numVertices == cum, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected %" PetscInt_FMT " != %" PetscInt_FMT, numVertices, cum);
1418:     PetscCall(ISDestroy(&gid));
1419:     PetscCall(ISCreateGeneral(comm, numVertices, idxs2, PETSC_OWN_POINTER, &gid));

1421:     /* support for node-aware cell locality */
1422:     PetscCall(ISCreateGeneral(comm, start[numVertices], adjacency, PETSC_USE_POINTER, &acis));
1423:     PetscCall(VecCreateSeq(PETSC_COMM_SELF, start[numVertices], &acown));
1424:     PetscCall(VecCreateMPI(comm, numVertices, PETSC_DECIDE, &cown));
1425:     PetscCall(VecGetArray(cown, &array));
1426:     for (c = 0; c < numVertices; c++) array[c] = nid;
1427:     PetscCall(VecRestoreArray(cown, &array));
1428:     PetscCall(VecScatterCreate(cown, acis, acown, NULL, &sct));
1429:     PetscCall(VecScatterBegin(sct, cown, acown, INSERT_VALUES, SCATTER_FORWARD));
1430:     PetscCall(VecScatterEnd(sct, cown, acown, INSERT_VALUES, SCATTER_FORWARD));
1431:     PetscCall(ISDestroy(&acis));
1432:     PetscCall(VecScatterDestroy(&sct));
1433:     PetscCall(VecDestroy(&cown));

1435:     /* compute edgeCut */
1436:     for (c = 0, cum = 0; c < numVertices; c++) cum = PetscMax(cum, start[c + 1] - start[c]);
1437:     PetscCall(PetscMalloc1(cum, &work));
1438:     PetscCall(ISLocalToGlobalMappingCreateIS(gid, &g2l));
1439:     PetscCall(ISLocalToGlobalMappingSetType(g2l, ISLOCALTOGLOBALMAPPINGHASH));
1440:     PetscCall(ISDestroy(&gid));
1441:     PetscCall(VecGetArray(acown, &array));
1442:     for (c = 0, ect = 0, ectn = 0; c < numVertices; c++) {
1443:       PetscInt totl;

1445:       totl = start[c + 1] - start[c];
1446:       PetscCall(ISGlobalToLocalMappingApply(g2l, IS_GTOLM_MASK, totl, adjacency + start[c], NULL, work));
1447:       for (i = 0; i < totl; i++) {
1448:         if (work[i] < 0) {
1449:           ect += 1;
1450:           ectn += (array[i + start[c]] != nid) ? 0 : 1;
1451:         }
1452:       }
1453:     }
1454:     PetscCall(PetscFree(work));
1455:     PetscCall(VecRestoreArray(acown, &array));
1456:     lm[0] = numVertices > 0 ? numVertices : PETSC_INT_MAX;
1457:     lm[1] = -numVertices;
1458:     PetscCallMPI(MPIU_Allreduce(lm, gm, 2, MPIU_INT64, MPI_MIN, comm));
1459:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Cell balance: %.2f (max %" PetscInt_FMT ", min %" PetscInt_FMT, -((double)gm[1]) / ((double)gm[0]), -(PetscInt)gm[1], (PetscInt)gm[0]));
1460:     lm[0] = ect;                     /* edgeCut */
1461:     lm[1] = ectn;                    /* node-aware edgeCut */
1462:     lm[2] = numVertices > 0 ? 0 : 1; /* empty processes */
1463:     PetscCallMPI(MPIU_Allreduce(lm, gm, 3, MPIU_INT64, MPI_SUM, comm));
1464:     PetscCall(PetscViewerASCIIPrintf(viewer, ", empty %" PetscInt_FMT ")\n", (PetscInt)gm[2]));
1465: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
1466:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge Cut: %" PetscInt_FMT " (on node %.3f)\n", (PetscInt)(gm[0] / 2), gm[0] ? ((double)gm[1]) / ((double)gm[0]) : 1.));
1467: #else
1468:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge Cut: %" PetscInt_FMT " (on node %.3f)\n", (PetscInt)(gm[0] / 2), 0.0));
1469: #endif
1470:     PetscCall(ISLocalToGlobalMappingDestroy(&g2l));
1471:     PetscCall(PetscFree(start));
1472:     PetscCall(PetscFree(adjacency));
1473:     PetscCall(VecDestroy(&acown));
1474:   } else {
1475:     const char    *name;
1476:     PetscInt      *sizes, *hybsizes, *ghostsizes;
1477:     PetscInt       locDepth, depth, cellHeight, dim, d;
1478:     PetscInt       pStart, pEnd, p, gcStart, gcEnd, gcNum;
1479:     PetscInt       numLabels, l, maxSize = 17;
1480:     DMPolytopeType ct0 = DM_POLYTOPE_UNKNOWN;
1481:     MPI_Comm       comm;
1482:     PetscMPIInt    size, rank;

1484:     PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1485:     PetscCallMPI(MPI_Comm_size(comm, &size));
1486:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
1487:     PetscCall(DMGetDimension(dm, &dim));
1488:     PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1489:     PetscCall(PetscObjectGetName((PetscObject)dm, &name));
1490:     if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s"));
1491:     else PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s"));
1492:     if (cellHeight) PetscCall(PetscViewerASCIIPrintf(viewer, "  Cells are at height %" PetscInt_FMT "\n", cellHeight));
1493:     PetscCall(DMPlexGetDepth(dm, &locDepth));
1494:     PetscCallMPI(MPIU_Allreduce(&locDepth, &depth, 1, MPIU_INT, MPI_MAX, comm));
1495:     PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &gcStart, &gcEnd));
1496:     gcNum = gcEnd - gcStart;
1497:     if (size < maxSize) PetscCall(PetscCalloc3(size, &sizes, size, &hybsizes, size, &ghostsizes));
1498:     else PetscCall(PetscCalloc3(3, &sizes, 3, &hybsizes, 3, &ghostsizes));
1499:     for (d = 0; d <= depth; d++) {
1500:       PetscInt Nc[2] = {0, 0}, ict;

1502:       PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
1503:       if (pStart < pEnd) PetscCall(DMPlexGetCellType(dm, pStart, &ct0));
1504:       ict = ct0;
1505:       PetscCallMPI(MPI_Bcast(&ict, 1, MPIU_INT, 0, comm));
1506:       ct0 = (DMPolytopeType)ict;
1507:       for (p = pStart; p < pEnd; ++p) {
1508:         DMPolytopeType ct;

1510:         PetscCall(DMPlexGetCellType(dm, p, &ct));
1511:         if (ct == ct0) ++Nc[0];
1512:         else ++Nc[1];
1513:       }
1514:       if (size < maxSize) {
1515:         PetscCallMPI(MPI_Gather(&Nc[0], 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm));
1516:         PetscCallMPI(MPI_Gather(&Nc[1], 1, MPIU_INT, hybsizes, 1, MPIU_INT, 0, comm));
1517:         if (d == depth) PetscCallMPI(MPI_Gather(&gcNum, 1, MPIU_INT, ghostsizes, 1, MPIU_INT, 0, comm));
1518:         PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of %" PetscInt_FMT "-cells per rank:", (depth == 1) && d ? dim : d));
1519:         for (p = 0; p < size; ++p) {
1520:           if (rank == 0) {
1521:             PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p] + hybsizes[p]));
1522:             if (hybsizes[p] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ")", hybsizes[p]));
1523:             if (ghostsizes[p] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " [%" PetscInt_FMT "]", ghostsizes[p]));
1524:           }
1525:         }
1526:       } else {
1527:         PetscInt locMinMax[2];

1529:         locMinMax[0] = Nc[0] + Nc[1];
1530:         locMinMax[1] = Nc[0] + Nc[1];
1531:         PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes));
1532:         locMinMax[0] = Nc[1];
1533:         locMinMax[1] = Nc[1];
1534:         PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, hybsizes));
1535:         if (d == depth) {
1536:           locMinMax[0] = gcNum;
1537:           locMinMax[1] = gcNum;
1538:           PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, ghostsizes));
1539:         }
1540:         PetscCall(PetscViewerASCIIPrintf(viewer, "  Min/Max of %" PetscInt_FMT "-cells per rank:", (depth == 1) && d ? dim : d));
1541:         PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1]));
1542:         if (hybsizes[0] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT "/%" PetscInt_FMT ")", hybsizes[0], hybsizes[1]));
1543:         if (ghostsizes[0] > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " [%" PetscInt_FMT "/%" PetscInt_FMT "]", ghostsizes[0], ghostsizes[1]));
1544:       }
1545:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1546:     }
1547:     PetscCall(PetscFree3(sizes, hybsizes, ghostsizes));
1548:     {
1549:       const PetscReal *maxCell;
1550:       const PetscReal *L;
1551:       PetscBool        localized;

1553:       PetscCall(DMGetPeriodicity(dm, &maxCell, NULL, &L));
1554:       PetscCall(DMGetCoordinatesLocalized(dm, &localized));
1555:       if (L || localized) {
1556:         PetscCall(PetscViewerASCIIPrintf(viewer, "Periodic mesh"));
1557:         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1558:         if (L) {
1559:           PetscCall(PetscViewerASCIIPrintf(viewer, " ("));
1560:           for (d = 0; d < dim; ++d) {
1561:             if (d > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1562:             PetscCall(PetscViewerASCIIPrintf(viewer, "%s", L[d] > 0.0 ? "PERIODIC" : "NONE"));
1563:           }
1564:           PetscCall(PetscViewerASCIIPrintf(viewer, ")"));
1565:         }
1566:         PetscCall(PetscViewerASCIIPrintf(viewer, " coordinates %s\n", localized ? "localized" : "not localized"));
1567:         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1568:       }
1569:     }
1570:     PetscCall(DMGetNumLabels(dm, &numLabels));
1571:     if (numLabels) PetscCall(PetscViewerASCIIPrintf(viewer, "Labels:\n"));
1572:     for (l = 0; l < numLabels; ++l) {
1573:       DMLabel     label;
1574:       const char *name;
1575:       PetscInt   *values;
1576:       PetscInt    numValues, v;

1578:       PetscCall(DMGetLabelName(dm, l, &name));
1579:       PetscCall(DMGetLabel(dm, name, &label));
1580:       PetscCall(DMLabelGetNumValues(label, &numValues));
1581:       PetscCall(PetscViewerASCIIPrintf(viewer, "  %s: %" PetscInt_FMT " strata with value/size (", name, numValues));

1583:       { // Extract array of DMLabel values so it can be sorted
1584:         IS              is_values;
1585:         const PetscInt *is_values_local = NULL;

1587:         PetscCall(DMLabelGetValueIS(label, &is_values));
1588:         PetscCall(ISGetIndices(is_values, &is_values_local));
1589:         PetscCall(PetscMalloc1(numValues, &values));
1590:         PetscCall(PetscArraycpy(values, is_values_local, numValues));
1591:         PetscCall(PetscSortInt(numValues, values));
1592:         PetscCall(ISRestoreIndices(is_values, &is_values_local));
1593:         PetscCall(ISDestroy(&is_values));
1594:       }
1595:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1596:       for (v = 0; v < numValues; ++v) {
1597:         PetscInt size;

1599:         PetscCall(DMLabelGetStratumSize(label, values[v], &size));
1600:         if (v > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1601:         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " (%" PetscInt_FMT ")", values[v], size));
1602:       }
1603:       PetscCall(PetscViewerASCIIPrintf(viewer, ")\n"));
1604:       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1605:       PetscCall(PetscFree(values));
1606:     }
1607:     {
1608:       char    **labelNames;
1609:       PetscInt  Nl = numLabels;
1610:       PetscBool flg;

1612:       PetscCall(PetscMalloc1(Nl, &labelNames));
1613:       PetscCall(PetscOptionsGetStringArray(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_view_labels", labelNames, &Nl, &flg));
1614:       for (l = 0; l < Nl; ++l) {
1615:         DMLabel label;

1617:         PetscCall(DMHasLabel(dm, labelNames[l], &flg));
1618:         if (flg) {
1619:           PetscCall(DMGetLabel(dm, labelNames[l], &label));
1620:           PetscCall(DMLabelView(label, viewer));
1621:         }
1622:         PetscCall(PetscFree(labelNames[l]));
1623:       }
1624:       PetscCall(PetscFree(labelNames));
1625:     }
1626:     /* If no fields are specified, people do not want to see adjacency */
1627:     if (dm->Nf) {
1628:       PetscInt f;

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

1633:         PetscCall(PetscObjectGetName(dm->fields[f].disc, &name));
1634:         if (numLabels) PetscCall(PetscViewerASCIIPrintf(viewer, "Field %s:\n", name));
1635:         PetscCall(PetscViewerASCIIPushTab(viewer));
1636:         if (dm->fields[f].label) PetscCall(DMLabelView(dm->fields[f].label, viewer));
1637:         if (dm->fields[f].adjacency[0]) {
1638:           if (dm->fields[f].adjacency[1]) PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FVM++\n"));
1639:           else PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FVM\n"));
1640:         } else {
1641:           if (dm->fields[f].adjacency[1]) PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FEM\n"));
1642:           else PetscCall(PetscViewerASCIIPrintf(viewer, "adjacency FUNKY\n"));
1643:         }
1644:         PetscCall(PetscViewerASCIIPopTab(viewer));
1645:       }
1646:     }
1647:     PetscCall(DMGetCoarseDM(dm, &cdm));
1648:     if (cdm) {
1649:       PetscCall(PetscViewerASCIIPushTab(viewer));
1650:       PetscCall(PetscViewerASCIIPrintf(viewer, "Defined by transform from:\n"));
1651:       PetscCall(DMPlexView_Ascii(cdm, viewer));
1652:       PetscCall(PetscViewerASCIIPopTab(viewer));
1653:     }
1654:   }
1655:   PetscFunctionReturn(PETSC_SUCCESS);
1656: }

1658: static PetscErrorCode DMPlexDrawCell(DM dm, PetscDraw draw, PetscInt cell, const PetscScalar coords[])
1659: {
1660:   DMPolytopeType ct;
1661:   PetscMPIInt    rank;
1662:   PetscInt       cdim;

1664:   PetscFunctionBegin;
1665:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1666:   PetscCall(DMPlexGetCellType(dm, cell, &ct));
1667:   PetscCall(DMGetCoordinateDim(dm, &cdim));
1668:   switch (ct) {
1669:   case DM_POLYTOPE_SEGMENT:
1670:   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1671:     switch (cdim) {
1672:     case 1: {
1673:       const PetscReal y  = 0.5;  /* TODO Put it in the middle of the viewport */
1674:       const PetscReal dy = 0.05; /* TODO Make it a fraction of the total length */

1676:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), y, PetscRealPart(coords[1]), y, PETSC_DRAW_BLACK));
1677:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), y + dy, PetscRealPart(coords[0]), y - dy, PETSC_DRAW_BLACK));
1678:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[1]), y + dy, PetscRealPart(coords[1]), y - dy, PETSC_DRAW_BLACK));
1679:     } break;
1680:     case 2: {
1681:       const PetscReal dx = (PetscRealPart(coords[3]) - PetscRealPart(coords[1]));
1682:       const PetscReal dy = (PetscRealPart(coords[2]) - PetscRealPart(coords[0]));
1683:       const PetscReal l  = 0.1 / PetscSqrtReal(dx * dx + dy * dy);

1685:       PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1686:       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));
1687:       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));
1688:     } break;
1689:     default:
1690:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of dimension %" PetscInt_FMT, cdim);
1691:     }
1692:     break;
1693:   case DM_POLYTOPE_TRIANGLE:
1694:     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));
1695:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1696:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK));
1697:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK));
1698:     break;
1699:   case DM_POLYTOPE_QUADRILATERAL:
1700:     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));
1701:     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));
1702:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1703:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK));
1704:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_BLACK));
1705:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK));
1706:     break;
1707:   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1708:     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));
1709:     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));
1710:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK));
1711:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_BLACK));
1712:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK));
1713:     PetscCall(PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK));
1714:     break;
1715:   case DM_POLYTOPE_FV_GHOST:
1716:     break;
1717:   default:
1718:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of type %s", DMPolytopeTypes[ct]);
1719:   }
1720:   PetscFunctionReturn(PETSC_SUCCESS);
1721: }

1723: static PetscErrorCode DrawPolygon_Private(DM dm, PetscDraw draw, PetscInt cell, PetscInt Nv, const PetscReal refVertices[], const PetscScalar coords[], PetscInt edgeDiv, PetscReal refCoords[], PetscReal edgeCoords[])
1724: {
1725:   PetscReal   centroid[2] = {0., 0.};
1726:   PetscMPIInt rank;
1727:   PetscMPIInt fillColor;

1729:   PetscFunctionBegin;
1730:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1731:   fillColor = PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS - 2) + 2;
1732:   for (PetscInt v = 0; v < Nv; ++v) {
1733:     centroid[0] += PetscRealPart(coords[v * 2 + 0]) / Nv;
1734:     centroid[1] += PetscRealPart(coords[v * 2 + 1]) / Nv;
1735:   }
1736:   for (PetscInt e = 0; e < Nv; ++e) {
1737:     refCoords[0] = refVertices[e * 2 + 0];
1738:     refCoords[1] = refVertices[e * 2 + 1];
1739:     for (PetscInt d = 1; d <= edgeDiv; ++d) {
1740:       refCoords[d * 2 + 0] = refCoords[0] + (refVertices[(e + 1) % Nv * 2 + 0] - refCoords[0]) * d / edgeDiv;
1741:       refCoords[d * 2 + 1] = refCoords[1] + (refVertices[(e + 1) % Nv * 2 + 1] - refCoords[1]) * d / edgeDiv;
1742:     }
1743:     PetscCall(DMPlexReferenceToCoordinates(dm, cell, edgeDiv + 1, refCoords, edgeCoords));
1744:     for (PetscInt d = 0; d < edgeDiv; ++d) {
1745:       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));
1746:       PetscCall(PetscDrawLine(draw, edgeCoords[d * 2 + 0], edgeCoords[d * 2 + 1], edgeCoords[(d + 1) * 2 + 0], edgeCoords[(d + 1) * 2 + 1], PETSC_DRAW_BLACK));
1747:     }
1748:   }
1749:   PetscFunctionReturn(PETSC_SUCCESS);
1750: }

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

1756:   PetscFunctionBegin;
1757:   PetscCall(DMPlexGetCellType(dm, cell, &ct));
1758:   switch (ct) {
1759:   case DM_POLYTOPE_TRIANGLE: {
1760:     PetscReal refVertices[6] = {-1., -1., 1., -1., -1., 1.};

1762:     PetscCall(DrawPolygon_Private(dm, draw, cell, 3, refVertices, coords, edgeDiv, refCoords, edgeCoords));
1763:   } break;
1764:   case DM_POLYTOPE_QUADRILATERAL: {
1765:     PetscReal refVertices[8] = {-1., -1., 1., -1., 1., 1., -1., 1.};

1767:     PetscCall(DrawPolygon_Private(dm, draw, cell, 4, refVertices, coords, edgeDiv, refCoords, edgeCoords));
1768:   } break;
1769:   default:
1770:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of type %s", DMPolytopeTypes[ct]);
1771:   }
1772:   PetscFunctionReturn(PETSC_SUCCESS);
1773: }

1775: static PetscErrorCode DMPlexView_Draw(DM dm, PetscViewer viewer)
1776: {
1777:   PetscDraw    draw;
1778:   DM           cdm;
1779:   PetscSection coordSection;
1780:   Vec          coordinates;
1781:   PetscReal    xyl[3], xyr[3];
1782:   PetscReal   *refCoords, *edgeCoords;
1783:   PetscBool    isnull, drawAffine;
1784:   PetscInt     dim, vStart, vEnd, cStart, cEnd, c, cDegree, edgeDiv;

1786:   PetscFunctionBegin;
1787:   PetscCall(DMGetCoordinateDim(dm, &dim));
1788:   PetscCheck(dim <= 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot draw meshes of dimension %" PetscInt_FMT, dim);
1789:   PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree));
1790:   drawAffine = cDegree > 1 ? PETSC_FALSE : PETSC_TRUE;
1791:   edgeDiv    = cDegree + 1;
1792:   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_view_draw_affine", &drawAffine, NULL));
1793:   if (!drawAffine) PetscCall(PetscMalloc2((edgeDiv + 1) * dim, &refCoords, (edgeDiv + 1) * dim, &edgeCoords));
1794:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1795:   PetscCall(DMGetLocalSection(cdm, &coordSection));
1796:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1797:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1798:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));

1800:   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
1801:   PetscCall(PetscDrawIsNull(draw, &isnull));
1802:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1803:   PetscCall(PetscDrawSetTitle(draw, "Mesh"));

1805:   PetscCall(DMGetBoundingBox(dm, xyl, xyr));
1806:   PetscCall(PetscDrawSetCoordinates(draw, xyl[0], xyl[1], xyr[0], xyr[1]));
1807:   PetscCall(PetscDrawClear(draw));

1809:   for (c = cStart; c < cEnd; ++c) {
1810:     PetscScalar       *coords = NULL;
1811:     const PetscScalar *coords_arr;
1812:     PetscInt           numCoords;
1813:     PetscBool          isDG;

1815:     PetscCall(DMPlexGetCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
1816:     if (drawAffine) PetscCall(DMPlexDrawCell(dm, draw, c, coords));
1817:     else PetscCall(DMPlexDrawCellHighOrder(dm, draw, c, coords, edgeDiv, refCoords, edgeCoords));
1818:     PetscCall(DMPlexRestoreCellCoordinates(dm, c, &isDG, &numCoords, &coords_arr, &coords));
1819:   }
1820:   if (!drawAffine) PetscCall(PetscFree2(refCoords, edgeCoords));
1821:   PetscCall(PetscDrawFlush(draw));
1822:   PetscCall(PetscDrawPause(draw));
1823:   PetscCall(PetscDrawSave(draw));
1824:   PetscFunctionReturn(PETSC_SUCCESS);
1825: }

1827: static PetscErrorCode DMPlexCreateHighOrderSurrogate_Internal(DM dm, DM *hdm)
1828: {
1829:   DM           odm = dm, rdm = dm, cdm;
1830:   PetscFE      fe;
1831:   PetscSpace   sp;
1832:   PetscClassId id;
1833:   PetscInt     degree;
1834:   PetscBool    hoView = PETSC_TRUE;

1836:   PetscFunctionBegin;
1837:   PetscObjectOptionsBegin((PetscObject)dm);
1838:   PetscCall(PetscOptionsBool("-dm_plex_high_order_view", "Subsample to view meshes with high order coordinates", "DMPlexCreateHighOrderSurrogate_Internal", hoView, &hoView, NULL));
1839:   PetscOptionsEnd();
1840:   PetscCall(PetscObjectReference((PetscObject)dm));
1841:   *hdm = dm;
1842:   if (!hoView) PetscFunctionReturn(PETSC_SUCCESS);
1843:   PetscCall(DMGetCoordinateDM(dm, &cdm));
1844:   PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
1845:   PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
1846:   if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
1847:   PetscCall(PetscFEGetBasisSpace(fe, &sp));
1848:   PetscCall(PetscSpaceGetDegree(sp, &degree, NULL));
1849:   for (PetscInt r = 0, rd = PetscCeilReal(((PetscReal)degree) / 2.); r < (PetscInt)PetscCeilReal(PetscLog2Real(degree)); ++r, rd = PetscCeilReal(((PetscReal)rd) / 2.)) {
1850:     DM  cdm, rcdm;
1851:     Mat In;
1852:     Vec cl, rcl;

1854:     PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)odm), &rdm));
1855:     PetscCall(DMPlexCreateCoordinateSpace(rdm, rd, PETSC_FALSE, NULL));
1856:     PetscCall(PetscObjectSetName((PetscObject)rdm, "Refined Mesh with Linear Coordinates"));
1857:     PetscCall(DMGetCoordinateDM(odm, &cdm));
1858:     PetscCall(DMGetCoordinateDM(rdm, &rcdm));
1859:     PetscCall(DMGetCoordinatesLocal(odm, &cl));
1860:     PetscCall(DMGetCoordinatesLocal(rdm, &rcl));
1861:     PetscCall(DMSetCoarseDM(rcdm, cdm));
1862:     PetscCall(DMCreateInterpolation(cdm, rcdm, &In, NULL));
1863:     PetscCall(MatMult(In, cl, rcl));
1864:     PetscCall(MatDestroy(&In));
1865:     PetscCall(DMSetCoordinatesLocal(rdm, rcl));
1866:     PetscCall(DMDestroy(&odm));
1867:     odm = rdm;
1868:   }
1869:   *hdm = rdm;
1870:   PetscFunctionReturn(PETSC_SUCCESS);
1871: }

1873: #if defined(PETSC_HAVE_EXODUSII)
1874:   #include <exodusII.h>
1875: #include <petscviewerexodusii.h>
1876: #endif

1878: PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer)
1879: {
1880:   PetscBool iascii, ishdf5, isvtk, isdraw, flg, isglvis, isexodus, iscgns;
1881:   char      name[PETSC_MAX_PATH_LEN];

1883:   PetscFunctionBegin;
1886:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1887:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
1888:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
1889:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1890:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
1891:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWEREXODUSII, &isexodus));
1892:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERCGNS, &iscgns));
1893:   if (iascii) {
1894:     PetscViewerFormat format;
1895:     PetscCall(PetscViewerGetFormat(viewer, &format));
1896:     if (format == PETSC_VIEWER_ASCII_GLVIS) PetscCall(DMPlexView_GLVis(dm, viewer));
1897:     else PetscCall(DMPlexView_Ascii(dm, viewer));
1898:   } else if (ishdf5) {
1899: #if defined(PETSC_HAVE_HDF5)
1900:     PetscCall(DMPlexView_HDF5_Internal(dm, viewer));
1901: #else
1902:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1903: #endif
1904:   } else if (isvtk) {
1905:     PetscCall(DMPlexVTKWriteAll((PetscObject)dm, viewer));
1906:   } else if (isdraw) {
1907:     DM hdm;

1909:     PetscCall(DMPlexCreateHighOrderSurrogate_Internal(dm, &hdm));
1910:     PetscCall(DMPlexView_Draw(hdm, viewer));
1911:     PetscCall(DMDestroy(&hdm));
1912:   } else if (isglvis) {
1913:     PetscCall(DMPlexView_GLVis(dm, viewer));
1914: #if defined(PETSC_HAVE_EXODUSII)
1915:   } else if (isexodus) {
1916:     /*
1917:       exodusII requires that all sets be part of exactly one cell set.
1918:       If the dm does not have a "Cell Sets" label defined, we create one
1919:       with ID 1, containing all cells.
1920:       Note that if the Cell Sets label is defined but does not cover all cells,
1921:       we may still have a problem. This should probably be checked here or in the viewer;
1922:     */
1923:     PetscInt numCS;
1924:     PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS));
1925:     if (!numCS) {
1926:       PetscInt cStart, cEnd, c;
1927:       PetscCall(DMCreateLabel(dm, "Cell Sets"));
1928:       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1929:       for (c = cStart; c < cEnd; ++c) PetscCall(DMSetLabelValue(dm, "Cell Sets", c, 1));
1930:     }
1931:     PetscCall(DMView_PlexExodusII(dm, viewer));
1932: #endif
1933: #if defined(PETSC_HAVE_CGNS)
1934:   } else if (iscgns) {
1935:     PetscCall(DMView_PlexCGNS(dm, viewer));
1936: #endif
1937:   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlex writing", ((PetscObject)viewer)->type_name);
1938:   /* Optionally view the partition */
1939:   PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_partition_view", &flg));
1940:   if (flg) {
1941:     Vec ranks;
1942:     PetscCall(DMPlexCreateRankField(dm, &ranks));
1943:     PetscCall(VecView(ranks, viewer));
1944:     PetscCall(VecDestroy(&ranks));
1945:   }
1946:   /* Optionally view a label */
1947:   PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_label_view", name, sizeof(name), &flg));
1948:   if (flg) {
1949:     DMLabel label;
1950:     Vec     val;

1952:     PetscCall(DMGetLabel(dm, name, &label));
1953:     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Label %s provided to -dm_label_view does not exist in this DM", name);
1954:     PetscCall(DMPlexCreateLabelField(dm, label, &val));
1955:     PetscCall(VecView(val, viewer));
1956:     PetscCall(VecDestroy(&val));
1957:   }
1958:   PetscFunctionReturn(PETSC_SUCCESS);
1959: }

1961: /*@
1962:   DMPlexTopologyView - Saves a `DMPLEX` topology into a file

1964:   Collective

1966:   Input Parameters:
1967: + dm     - The `DM` whose topology is to be saved
1968: - viewer - The `PetscViewer` to save it in

1970:   Level: advanced

1972: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexCoordinatesView()`, `DMPlexLabelsView()`, `DMPlexTopologyLoad()`, `PetscViewer`
1973: @*/
1974: PetscErrorCode DMPlexTopologyView(DM dm, PetscViewer viewer)
1975: {
1976:   PetscBool ishdf5;

1978:   PetscFunctionBegin;
1981:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
1982:   PetscCall(PetscLogEventBegin(DMPLEX_TopologyView, viewer, 0, 0, 0));
1983:   if (ishdf5) {
1984: #if defined(PETSC_HAVE_HDF5)
1985:     PetscViewerFormat format;
1986:     PetscCall(PetscViewerGetFormat(viewer, &format));
1987:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1988:       IS globalPointNumbering;

1990:       PetscCall(DMPlexCreatePointNumbering(dm, &globalPointNumbering));
1991:       PetscCall(DMPlexTopologyView_HDF5_Internal(dm, globalPointNumbering, viewer));
1992:       PetscCall(ISDestroy(&globalPointNumbering));
1993:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
1994: #else
1995:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1996: #endif
1997:   }
1998:   PetscCall(PetscLogEventEnd(DMPLEX_TopologyView, viewer, 0, 0, 0));
1999:   PetscFunctionReturn(PETSC_SUCCESS);
2000: }

2002: /*@
2003:   DMPlexCoordinatesView - Saves `DMPLEX` coordinates into a file

2005:   Collective

2007:   Input Parameters:
2008: + dm     - The `DM` whose coordinates are to be saved
2009: - viewer - The `PetscViewer` for saving

2011:   Level: advanced

2013: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexTopologyView()`, `DMPlexLabelsView()`, `DMPlexCoordinatesLoad()`, `PetscViewer`
2014: @*/
2015: PetscErrorCode DMPlexCoordinatesView(DM dm, PetscViewer viewer)
2016: {
2017:   PetscBool ishdf5;

2019:   PetscFunctionBegin;
2022:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2023:   PetscCall(PetscLogEventBegin(DMPLEX_CoordinatesView, viewer, 0, 0, 0));
2024:   if (ishdf5) {
2025: #if defined(PETSC_HAVE_HDF5)
2026:     PetscViewerFormat format;
2027:     PetscCall(PetscViewerGetFormat(viewer, &format));
2028:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2029:       PetscCall(DMPlexCoordinatesView_HDF5_Internal(dm, viewer));
2030:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
2031: #else
2032:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2033: #endif
2034:   }
2035:   PetscCall(PetscLogEventEnd(DMPLEX_CoordinatesView, viewer, 0, 0, 0));
2036:   PetscFunctionReturn(PETSC_SUCCESS);
2037: }

2039: /*@
2040:   DMPlexLabelsView - Saves `DMPLEX` labels into a file

2042:   Collective

2044:   Input Parameters:
2045: + dm     - The `DM` whose labels are to be saved
2046: - viewer - The `PetscViewer` for saving

2048:   Level: advanced

2050: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexTopologyView()`, `DMPlexCoordinatesView()`, `DMPlexLabelsLoad()`, `PetscViewer`
2051: @*/
2052: PetscErrorCode DMPlexLabelsView(DM dm, PetscViewer viewer)
2053: {
2054:   PetscBool ishdf5;

2056:   PetscFunctionBegin;
2059:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2060:   PetscCall(PetscLogEventBegin(DMPLEX_LabelsView, viewer, 0, 0, 0));
2061:   if (ishdf5) {
2062: #if defined(PETSC_HAVE_HDF5)
2063:     IS                globalPointNumbering;
2064:     PetscViewerFormat format;

2066:     PetscCall(PetscViewerGetFormat(viewer, &format));
2067:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2068:       PetscCall(DMPlexCreatePointNumbering(dm, &globalPointNumbering));
2069:       PetscCall(DMPlexLabelsView_HDF5_Internal(dm, globalPointNumbering, viewer));
2070:       PetscCall(ISDestroy(&globalPointNumbering));
2071:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2072: #else
2073:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2074: #endif
2075:   }
2076:   PetscCall(PetscLogEventEnd(DMPLEX_LabelsView, viewer, 0, 0, 0));
2077:   PetscFunctionReturn(PETSC_SUCCESS);
2078: }

2080: /*@
2081:   DMPlexSectionView - Saves a section associated with a `DMPLEX`

2083:   Collective

2085:   Input Parameters:
2086: + dm        - The `DM` that contains the topology on which the section to be saved is defined
2087: . viewer    - The `PetscViewer` for saving
2088: - sectiondm - The `DM` that contains the section to be saved, can be `NULL`

2090:   Level: advanced

2092:   Notes:
2093:   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.

2095:   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.

2097: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMView()`, `DMPlexTopologyView()`, `DMPlexCoordinatesView()`, `DMPlexLabelsView()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`, `PetscSectionView()`, `DMPlexSectionLoad()`, `PetscViewer`
2098: @*/
2099: PetscErrorCode DMPlexSectionView(DM dm, PetscViewer viewer, DM sectiondm)
2100: {
2101:   PetscBool ishdf5;

2103:   PetscFunctionBegin;
2106:   if (!sectiondm) sectiondm = dm;
2108:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2109:   PetscCall(PetscLogEventBegin(DMPLEX_SectionView, viewer, 0, 0, 0));
2110:   if (ishdf5) {
2111: #if defined(PETSC_HAVE_HDF5)
2112:     PetscCall(DMPlexSectionView_HDF5_Internal(dm, viewer, sectiondm));
2113: #else
2114:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2115: #endif
2116:   }
2117:   PetscCall(PetscLogEventEnd(DMPLEX_SectionView, viewer, 0, 0, 0));
2118:   PetscFunctionReturn(PETSC_SUCCESS);
2119: }

2121: /*@
2122:   DMPlexGlobalVectorView - Saves a global vector

2124:   Collective

2126:   Input Parameters:
2127: + dm        - The `DM` that represents the topology
2128: . viewer    - The `PetscViewer` to save data with
2129: . sectiondm - The `DM` that contains the global section on which vec is defined, can be `NULL`
2130: - vec       - The global vector to be saved

2132:   Level: advanced

2134:   Notes:
2135:   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.

2137:   Calling sequence:
2138: .vb
2139:        DMCreate(PETSC_COMM_WORLD, &dm);
2140:        DMSetType(dm, DMPLEX);
2141:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2142:        DMClone(dm, &sectiondm);
2143:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2144:        PetscSectionCreate(PETSC_COMM_WORLD, &section);
2145:        DMPlexGetChart(sectiondm, &pStart, &pEnd);
2146:        PetscSectionSetChart(section, pStart, pEnd);
2147:        PetscSectionSetUp(section);
2148:        DMSetLocalSection(sectiondm, section);
2149:        PetscSectionDestroy(&section);
2150:        DMGetGlobalVector(sectiondm, &vec);
2151:        PetscObjectSetName((PetscObject)vec, "vec_name");
2152:        DMPlexTopologyView(dm, viewer);
2153:        DMPlexSectionView(dm, viewer, sectiondm);
2154:        DMPlexGlobalVectorView(dm, viewer, sectiondm, vec);
2155:        DMRestoreGlobalVector(sectiondm, &vec);
2156:        DMDestroy(&sectiondm);
2157:        DMDestroy(&dm);
2158: .ve

2160: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyView()`, `DMPlexSectionView()`, `DMPlexLocalVectorView()`, `DMPlexGlobalVectorLoad()`, `DMPlexLocalVectorLoad()`
2161: @*/
2162: PetscErrorCode DMPlexGlobalVectorView(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
2163: {
2164:   PetscBool ishdf5;

2166:   PetscFunctionBegin;
2169:   if (!sectiondm) sectiondm = dm;
2172:   /* Check consistency */
2173:   {
2174:     PetscSection section;
2175:     PetscBool    includesConstraints;
2176:     PetscInt     m, m1;

2178:     PetscCall(VecGetLocalSize(vec, &m1));
2179:     PetscCall(DMGetGlobalSection(sectiondm, &section));
2180:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2181:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2182:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2183:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global vector size (%" PetscInt_FMT ") != global section storage size (%" PetscInt_FMT ")", m1, m);
2184:   }
2185:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2186:   PetscCall(PetscLogEventBegin(DMPLEX_GlobalVectorView, viewer, 0, 0, 0));
2187:   if (ishdf5) {
2188: #if defined(PETSC_HAVE_HDF5)
2189:     PetscCall(DMPlexGlobalVectorView_HDF5_Internal(dm, viewer, sectiondm, vec));
2190: #else
2191:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2192: #endif
2193:   }
2194:   PetscCall(PetscLogEventEnd(DMPLEX_GlobalVectorView, viewer, 0, 0, 0));
2195:   PetscFunctionReturn(PETSC_SUCCESS);
2196: }

2198: /*@
2199:   DMPlexLocalVectorView - Saves a local vector

2201:   Collective

2203:   Input Parameters:
2204: + dm        - The `DM` that represents the topology
2205: . viewer    - The `PetscViewer` to save data with
2206: . sectiondm - The `DM` that contains the local section on which `vec` is defined, can be `NULL`
2207: - vec       - The local vector to be saved

2209:   Level: advanced

2211:   Note:
2212:   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.

2214:   Calling sequence:
2215: .vb
2216:        DMCreate(PETSC_COMM_WORLD, &dm);
2217:        DMSetType(dm, DMPLEX);
2218:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2219:        DMClone(dm, &sectiondm);
2220:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2221:        PetscSectionCreate(PETSC_COMM_WORLD, &section);
2222:        DMPlexGetChart(sectiondm, &pStart, &pEnd);
2223:        PetscSectionSetChart(section, pStart, pEnd);
2224:        PetscSectionSetUp(section);
2225:        DMSetLocalSection(sectiondm, section);
2226:        DMGetLocalVector(sectiondm, &vec);
2227:        PetscObjectSetName((PetscObject)vec, "vec_name");
2228:        DMPlexTopologyView(dm, viewer);
2229:        DMPlexSectionView(dm, viewer, sectiondm);
2230:        DMPlexLocalVectorView(dm, viewer, sectiondm, vec);
2231:        DMRestoreLocalVector(sectiondm, &vec);
2232:        DMDestroy(&sectiondm);
2233:        DMDestroy(&dm);
2234: .ve

2236: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyView()`, `DMPlexSectionView()`, `DMPlexGlobalVectorView()`, `DMPlexGlobalVectorLoad()`, `DMPlexLocalVectorLoad()`
2237: @*/
2238: PetscErrorCode DMPlexLocalVectorView(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
2239: {
2240:   PetscBool ishdf5;

2242:   PetscFunctionBegin;
2245:   if (!sectiondm) sectiondm = dm;
2248:   /* Check consistency */
2249:   {
2250:     PetscSection section;
2251:     PetscBool    includesConstraints;
2252:     PetscInt     m, m1;

2254:     PetscCall(VecGetLocalSize(vec, &m1));
2255:     PetscCall(DMGetLocalSection(sectiondm, &section));
2256:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2257:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2258:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2259:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local vector size (%" PetscInt_FMT ") != local section storage size (%" PetscInt_FMT ")", m1, m);
2260:   }
2261:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2262:   PetscCall(PetscLogEventBegin(DMPLEX_LocalVectorView, viewer, 0, 0, 0));
2263:   if (ishdf5) {
2264: #if defined(PETSC_HAVE_HDF5)
2265:     PetscCall(DMPlexLocalVectorView_HDF5_Internal(dm, viewer, sectiondm, vec));
2266: #else
2267:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2268: #endif
2269:   }
2270:   PetscCall(PetscLogEventEnd(DMPLEX_LocalVectorView, viewer, 0, 0, 0));
2271:   PetscFunctionReturn(PETSC_SUCCESS);
2272: }

2274: PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer)
2275: {
2276:   PetscBool ishdf5;

2278:   PetscFunctionBegin;
2281:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2282:   if (ishdf5) {
2283: #if defined(PETSC_HAVE_HDF5)
2284:     PetscViewerFormat format;
2285:     PetscCall(PetscViewerGetFormat(viewer, &format));
2286:     if (format == PETSC_VIEWER_HDF5_XDMF || format == PETSC_VIEWER_HDF5_VIZ) {
2287:       PetscCall(DMPlexLoad_HDF5_Xdmf_Internal(dm, viewer));
2288:     } else if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2289:       PetscCall(DMPlexLoad_HDF5_Internal(dm, viewer));
2290:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2291:     PetscFunctionReturn(PETSC_SUCCESS);
2292: #else
2293:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2294: #endif
2295:   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlex loading", ((PetscObject)viewer)->type_name);
2296: }

2298: /*@
2299:   DMPlexTopologyLoad - Loads a topology into a `DMPLEX`

2301:   Collective

2303:   Input Parameters:
2304: + dm     - The `DM` into which the topology is loaded
2305: - viewer - The `PetscViewer` for the saved topology

2307:   Output Parameter:
2308: . 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;
2309:   `NULL` if unneeded

2311:   Level: advanced

2313: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexCoordinatesLoad()`, `DMPlexLabelsLoad()`, `DMView()`, `PetscViewerHDF5Open()`, `PetscViewerPushFormat()`,
2314:           `PetscViewer`, `PetscSF`
2315: @*/
2316: PetscErrorCode DMPlexTopologyLoad(DM dm, PetscViewer viewer, PetscSF *globalToLocalPointSF)
2317: {
2318:   PetscBool ishdf5;

2320:   PetscFunctionBegin;
2323:   if (globalToLocalPointSF) PetscAssertPointer(globalToLocalPointSF, 3);
2324:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2325:   PetscCall(PetscLogEventBegin(DMPLEX_TopologyLoad, viewer, 0, 0, 0));
2326:   if (ishdf5) {
2327: #if defined(PETSC_HAVE_HDF5)
2328:     PetscViewerFormat format;
2329:     PetscCall(PetscViewerGetFormat(viewer, &format));
2330:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2331:       PetscCall(DMPlexTopologyLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF));
2332:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2333: #else
2334:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2335: #endif
2336:   }
2337:   PetscCall(PetscLogEventEnd(DMPLEX_TopologyLoad, viewer, 0, 0, 0));
2338:   PetscFunctionReturn(PETSC_SUCCESS);
2339: }

2341: /*@
2342:   DMPlexCoordinatesLoad - Loads coordinates into a `DMPLEX`

2344:   Collective

2346:   Input Parameters:
2347: + dm                   - The `DM` into which the coordinates are loaded
2348: . viewer               - The `PetscViewer` for the saved coordinates
2349: - globalToLocalPointSF - The `PetscSF` returned by `DMPlexTopologyLoad()` when loading dm from viewer

2351:   Level: advanced

2353: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexTopologyLoad()`, `DMPlexLabelsLoad()`, `DMView()`, `PetscViewerHDF5Open()`, `PetscViewerPushFormat()`,
2354:           `PetscSF`, `PetscViewer`
2355: @*/
2356: PetscErrorCode DMPlexCoordinatesLoad(DM dm, PetscViewer viewer, PetscSF globalToLocalPointSF)
2357: {
2358:   PetscBool ishdf5;

2360:   PetscFunctionBegin;
2364:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2365:   PetscCall(PetscLogEventBegin(DMPLEX_CoordinatesLoad, viewer, 0, 0, 0));
2366:   if (ishdf5) {
2367: #if defined(PETSC_HAVE_HDF5)
2368:     PetscViewerFormat format;
2369:     PetscCall(PetscViewerGetFormat(viewer, &format));
2370:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2371:       PetscCall(DMPlexCoordinatesLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF));
2372:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2373: #else
2374:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2375: #endif
2376:   }
2377:   PetscCall(PetscLogEventEnd(DMPLEX_CoordinatesLoad, viewer, 0, 0, 0));
2378:   PetscFunctionReturn(PETSC_SUCCESS);
2379: }

2381: /*@
2382:   DMPlexLabelsLoad - Loads labels into a `DMPLEX`

2384:   Collective

2386:   Input Parameters:
2387: + dm                   - The `DM` into which the labels are loaded
2388: . viewer               - The `PetscViewer` for the saved labels
2389: - globalToLocalPointSF - The `PetscSF` returned by `DMPlexTopologyLoad()` when loading `dm` from viewer

2391:   Level: advanced

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

2396: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexTopologyLoad()`, `DMPlexCoordinatesLoad()`, `DMView()`, `PetscViewerHDF5Open()`, `PetscViewerPushFormat()`,
2397:           `PetscSF`, `PetscViewer`
2398: @*/
2399: PetscErrorCode DMPlexLabelsLoad(DM dm, PetscViewer viewer, PetscSF globalToLocalPointSF)
2400: {
2401:   PetscBool ishdf5;

2403:   PetscFunctionBegin;
2407:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2408:   PetscCall(PetscLogEventBegin(DMPLEX_LabelsLoad, viewer, 0, 0, 0));
2409:   if (ishdf5) {
2410: #if defined(PETSC_HAVE_HDF5)
2411:     PetscViewerFormat format;

2413:     PetscCall(PetscViewerGetFormat(viewer, &format));
2414:     if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
2415:       PetscCall(DMPlexLabelsLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF));
2416:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
2417: #else
2418:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2419: #endif
2420:   }
2421:   PetscCall(PetscLogEventEnd(DMPLEX_LabelsLoad, viewer, 0, 0, 0));
2422:   PetscFunctionReturn(PETSC_SUCCESS);
2423: }

2425: /*@
2426:   DMPlexSectionLoad - Loads section into a `DMPLEX`

2428:   Collective

2430:   Input Parameters:
2431: + dm                   - The `DM` that represents the topology
2432: . viewer               - The `PetscViewer` that represents the on-disk section (sectionA)
2433: . sectiondm            - The `DM` into which the on-disk section (sectionA) is migrated, can be `NULL`
2434: - globalToLocalPointSF - The `PetscSF` returned by `DMPlexTopologyLoad(`) when loading dm from viewer

2436:   Output Parameters:
2437: + 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)
2438: - 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)

2440:   Level: advanced

2442:   Notes:
2443:   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.

2445:   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.

2447:   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.

2449:   Example using 2 processes:
2450: .vb
2451:   NX (number of points on dm): 4
2452:   sectionA                   : the on-disk section
2453:   vecA                       : a vector associated with sectionA
2454:   sectionB                   : sectiondm's local section constructed in this function
2455:   vecB (local)               : a vector associated with sectiondm's local section
2456:   vecB (global)              : a vector associated with sectiondm's global section

2458:                                      rank 0    rank 1
2459:   vecA (global)                  : [.0 .4 .1 | .2 .3]        <- to be loaded in DMPlexGlobalVectorLoad() or DMPlexLocalVectorLoad()
2460:   sectionA->atlasOff             :       0 2 | 1             <- loaded in PetscSectionLoad()
2461:   sectionA->atlasDof             :       1 3 | 1             <- loaded in PetscSectionLoad()
2462:   sectionA's global point numbers:       0 2 | 3             <- loaded in DMPlexSectionLoad()
2463:   [0, NX)                        :       0 1 | 2 3           <- conceptual partition used in globalToLocalPointSF
2464:   sectionB's global point numbers:     0 1 3 | 3 2           <- associated with [0, NX) by globalToLocalPointSF
2465:   sectionB->atlasDof             :     1 0 1 | 1 3
2466:   sectionB->atlasOff (no perm)   :     0 1 1 | 0 1
2467:   vecB (local)                   :   [.0 .4] | [.4 .1 .2 .3] <- to be constructed by calling DMPlexLocalVectorLoad() with localDofSF
2468:   vecB (global)                  :    [.0 .4 | .1 .2 .3]     <- to be constructed by calling DMPlexGlobalVectorLoad() with globalDofSF
2469: .ve
2470:   where "|" represents a partition of loaded data, and global point 3 is assumed to be owned by rank 0.

2472: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLoad()`, `DMPlexTopologyLoad()`, `DMPlexCoordinatesLoad()`, `DMPlexLabelsLoad()`, `DMPlexGlobalVectorLoad()`, `DMPlexLocalVectorLoad()`, `PetscSectionLoad()`, `DMPlexSectionView()`, `PetscSF`, `PetscViewer`
2473: @*/
2474: PetscErrorCode DMPlexSectionLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF globalToLocalPointSF, PetscSF *globalDofSF, PetscSF *localDofSF)
2475: {
2476:   PetscBool ishdf5;

2478:   PetscFunctionBegin;
2481:   if (!sectiondm) sectiondm = dm;
2484:   if (globalDofSF) PetscAssertPointer(globalDofSF, 5);
2485:   if (localDofSF) PetscAssertPointer(localDofSF, 6);
2486:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2487:   PetscCall(PetscLogEventBegin(DMPLEX_SectionLoad, viewer, 0, 0, 0));
2488:   if (ishdf5) {
2489: #if defined(PETSC_HAVE_HDF5)
2490:     PetscCall(DMPlexSectionLoad_HDF5_Internal(dm, viewer, sectiondm, globalToLocalPointSF, globalDofSF, localDofSF));
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_SectionLoad, viewer, 0, 0, 0));
2496:   PetscFunctionReturn(PETSC_SUCCESS);
2497: }

2499: /*@
2500:   DMPlexGlobalVectorLoad - Loads on-disk vector data into a global vector

2502:   Collective

2504:   Input Parameters:
2505: + dm        - The `DM` that represents the topology
2506: . viewer    - The `PetscViewer` that represents the on-disk vector data
2507: . sectiondm - The `DM` that contains the global section on which vec is defined, can be `NULL`
2508: . sf        - The `PetscSF` that migrates the on-disk vector data into vec
2509: - vec       - The global vector to set values of

2511:   Level: advanced

2513:   Notes:
2514:   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.

2516:   Calling sequence:
2517: .vb
2518:        DMCreate(PETSC_COMM_WORLD, &dm);
2519:        DMSetType(dm, DMPLEX);
2520:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2521:        DMPlexTopologyLoad(dm, viewer, &sfX);
2522:        DMClone(dm, &sectiondm);
2523:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2524:        DMPlexSectionLoad(dm, viewer, sectiondm, sfX, &gsf, NULL);
2525:        DMGetGlobalVector(sectiondm, &vec);
2526:        PetscObjectSetName((PetscObject)vec, "vec_name");
2527:        DMPlexGlobalVectorLoad(dm, viewer, sectiondm, gsf, vec);
2528:        DMRestoreGlobalVector(sectiondm, &vec);
2529:        PetscSFDestroy(&gsf);
2530:        PetscSFDestroy(&sfX);
2531:        DMDestroy(&sectiondm);
2532:        DMDestroy(&dm);
2533: .ve

2535: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyLoad()`, `DMPlexSectionLoad()`, `DMPlexLocalVectorLoad()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`,
2536:           `PetscSF`, `PetscViewer`
2537: @*/
2538: PetscErrorCode DMPlexGlobalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2539: {
2540:   PetscBool ishdf5;

2542:   PetscFunctionBegin;
2545:   if (!sectiondm) sectiondm = dm;
2549:   /* Check consistency */
2550:   {
2551:     PetscSection section;
2552:     PetscBool    includesConstraints;
2553:     PetscInt     m, m1;

2555:     PetscCall(VecGetLocalSize(vec, &m1));
2556:     PetscCall(DMGetGlobalSection(sectiondm, &section));
2557:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2558:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2559:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2560:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global vector size (%" PetscInt_FMT ") != global section storage size (%" PetscInt_FMT ")", m1, m);
2561:   }
2562:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2563:   PetscCall(PetscLogEventBegin(DMPLEX_GlobalVectorLoad, viewer, 0, 0, 0));
2564:   if (ishdf5) {
2565: #if defined(PETSC_HAVE_HDF5)
2566:     PetscCall(DMPlexVecLoad_HDF5_Internal(dm, viewer, sectiondm, sf, vec));
2567: #else
2568:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2569: #endif
2570:   }
2571:   PetscCall(PetscLogEventEnd(DMPLEX_GlobalVectorLoad, viewer, 0, 0, 0));
2572:   PetscFunctionReturn(PETSC_SUCCESS);
2573: }

2575: /*@
2576:   DMPlexLocalVectorLoad - Loads on-disk vector data into a local vector

2578:   Collective

2580:   Input Parameters:
2581: + dm        - The `DM` that represents the topology
2582: . viewer    - The `PetscViewer` that represents the on-disk vector data
2583: . sectiondm - The `DM` that contains the local section on which vec is defined, can be `NULL`
2584: . sf        - The `PetscSF` that migrates the on-disk vector data into vec
2585: - vec       - The local vector to set values of

2587:   Level: advanced

2589:   Notes:
2590:   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.

2592:   Calling sequence:
2593: .vb
2594:        DMCreate(PETSC_COMM_WORLD, &dm);
2595:        DMSetType(dm, DMPLEX);
2596:        PetscObjectSetName((PetscObject)dm, "topologydm_name");
2597:        DMPlexTopologyLoad(dm, viewer, &sfX);
2598:        DMClone(dm, &sectiondm);
2599:        PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2600:        DMPlexSectionLoad(dm, viewer, sectiondm, sfX, NULL, &lsf);
2601:        DMGetLocalVector(sectiondm, &vec);
2602:        PetscObjectSetName((PetscObject)vec, "vec_name");
2603:        DMPlexLocalVectorLoad(dm, viewer, sectiondm, lsf, vec);
2604:        DMRestoreLocalVector(sectiondm, &vec);
2605:        PetscSFDestroy(&lsf);
2606:        PetscSFDestroy(&sfX);
2607:        DMDestroy(&sectiondm);
2608:        DMDestroy(&dm);
2609: .ve

2611: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTopologyLoad()`, `DMPlexSectionLoad()`, `DMPlexGlobalVectorLoad()`, `DMPlexGlobalVectorView()`, `DMPlexLocalVectorView()`,
2612:           `PetscSF`, `PetscViewer`
2613: @*/
2614: PetscErrorCode DMPlexLocalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2615: {
2616:   PetscBool ishdf5;

2618:   PetscFunctionBegin;
2621:   if (!sectiondm) sectiondm = dm;
2625:   /* Check consistency */
2626:   {
2627:     PetscSection section;
2628:     PetscBool    includesConstraints;
2629:     PetscInt     m, m1;

2631:     PetscCall(VecGetLocalSize(vec, &m1));
2632:     PetscCall(DMGetLocalSection(sectiondm, &section));
2633:     PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints));
2634:     if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m));
2635:     else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m));
2636:     PetscCheck(m1 == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local vector size (%" PetscInt_FMT ") != local section storage size (%" PetscInt_FMT ")", m1, m);
2637:   }
2638:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
2639:   PetscCall(PetscLogEventBegin(DMPLEX_LocalVectorLoad, viewer, 0, 0, 0));
2640:   if (ishdf5) {
2641: #if defined(PETSC_HAVE_HDF5)
2642:     PetscCall(DMPlexVecLoad_HDF5_Internal(dm, viewer, sectiondm, sf, vec));
2643: #else
2644:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2645: #endif
2646:   }
2647:   PetscCall(PetscLogEventEnd(DMPLEX_LocalVectorLoad, viewer, 0, 0, 0));
2648:   PetscFunctionReturn(PETSC_SUCCESS);
2649: }

2651: PetscErrorCode DMDestroy_Plex(DM dm)
2652: {
2653:   DM_Plex *mesh = (DM_Plex *)dm->data;

2655:   PetscFunctionBegin;
2656:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
2657:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", NULL));
2658:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", NULL));
2659:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMInterpolateSolution_C", NULL));
2660:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertTimeDerivativeBoundaryValues_C", NULL));
2661:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", NULL));
2662:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeGetDefault_C", NULL));
2663:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexDistributeSetDefault_C", NULL));
2664:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", NULL));
2665:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderGetDefault_C", NULL));
2666:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexReorderSetDefault_C", NULL));
2667:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionGetDefault_C", NULL));
2668:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionSetDefault_C", NULL));
2669:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionGetType_C", NULL));
2670:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMReorderSectionSetType_C", NULL));
2671:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", NULL));
2672:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetOverlap_C", NULL));
2673:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetUseCeed_C", NULL));
2674:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexSetUseCeed_C", NULL));
2675:   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMGetIsoperiodicPointSF_C", NULL));
2676:   if (--mesh->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
2677:   PetscCall(PetscSectionDestroy(&mesh->coneSection));
2678:   PetscCall(PetscFree(mesh->cones));
2679:   PetscCall(PetscFree(mesh->coneOrientations));
2680:   PetscCall(PetscSectionDestroy(&mesh->supportSection));
2681:   PetscCall(PetscSectionDestroy(&mesh->subdomainSection));
2682:   PetscCall(PetscFree(mesh->supports));
2683:   PetscCall(PetscFree(mesh->cellTypes));
2684:   PetscCall(DMPlexTransformDestroy(&mesh->tr));
2685:   PetscCall(PetscFree(mesh->tetgenOpts));
2686:   PetscCall(PetscFree(mesh->triangleOpts));
2687:   PetscCall(PetscFree(mesh->transformType));
2688:   PetscCall(PetscFree(mesh->distributionName));
2689:   PetscCall(PetscPartitionerDestroy(&mesh->partitioner));
2690:   PetscCall(DMLabelDestroy(&mesh->subpointMap));
2691:   PetscCall(ISDestroy(&mesh->subpointIS));
2692:   PetscCall(ISDestroy(&mesh->globalVertexNumbers));
2693:   PetscCall(ISDestroy(&mesh->globalCellNumbers));
2694:   if (mesh->periodic.face_sfs) {
2695:     for (PetscInt i = 0; i < mesh->periodic.num_face_sfs; i++) PetscCall(PetscSFDestroy(&mesh->periodic.face_sfs[i]));
2696:     PetscCall(PetscFree(mesh->periodic.face_sfs));
2697:   }
2698:   PetscCall(PetscSFDestroy(&mesh->periodic.composed_sf));
2699:   if (mesh->periodic.periodic_points) {
2700:     for (PetscInt i = 0; i < mesh->periodic.num_face_sfs; i++) PetscCall(ISDestroy(&mesh->periodic.periodic_points[i]));
2701:     PetscCall(PetscFree(mesh->periodic.periodic_points));
2702:   }
2703:   if (mesh->periodic.transform) PetscCall(PetscFree(mesh->periodic.transform));
2704:   PetscCall(PetscSectionDestroy(&mesh->anchorSection));
2705:   PetscCall(ISDestroy(&mesh->anchorIS));
2706:   PetscCall(PetscSectionDestroy(&mesh->parentSection));
2707:   PetscCall(PetscFree(mesh->parents));
2708:   PetscCall(PetscFree(mesh->childIDs));
2709:   PetscCall(PetscSectionDestroy(&mesh->childSection));
2710:   PetscCall(PetscFree(mesh->children));
2711:   PetscCall(DMDestroy(&mesh->referenceTree));
2712:   PetscCall(PetscGridHashDestroy(&mesh->lbox));
2713:   PetscCall(PetscFree(mesh->neighbors));
2714:   if (mesh->metricCtx) PetscCall(PetscFree(mesh->metricCtx));
2715:   if (mesh->nonempty_comm != MPI_COMM_NULL && mesh->nonempty_comm != MPI_COMM_SELF) PetscCallMPI(MPI_Comm_free(&mesh->nonempty_comm));
2716:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
2717:   PetscCall(PetscFree(mesh));
2718:   PetscFunctionReturn(PETSC_SUCCESS);
2719: }

2721: PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J)
2722: {
2723:   PetscSection           sectionGlobal, sectionLocal;
2724:   PetscInt               bs = -1, mbs;
2725:   PetscInt               localSize, localStart = 0;
2726:   PetscBool              isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isMatIS;
2727:   MatType                mtype;
2728:   ISLocalToGlobalMapping ltog;

2730:   PetscFunctionBegin;
2731:   PetscCall(MatInitializePackage());
2732:   mtype = dm->mattype;
2733:   PetscCall(DMGetLocalSection(dm, &sectionLocal));
2734:   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
2735:   /* PetscCall(PetscSectionGetStorageSize(sectionGlobal, &localSize)); */
2736:   PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize));
2737:   PetscCallMPI(MPI_Exscan(&localSize, &localStart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
2738:   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
2739:   PetscCall(MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE));
2740:   PetscCall(MatSetType(*J, mtype));
2741:   PetscCall(MatSetFromOptions(*J));
2742:   PetscCall(MatGetBlockSize(*J, &mbs));
2743:   if (mbs > 1) bs = mbs;
2744:   PetscCall(PetscStrcmp(mtype, MATSHELL, &isShell));
2745:   PetscCall(PetscStrcmp(mtype, MATBAIJ, &isBlock));
2746:   PetscCall(PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock));
2747:   PetscCall(PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock));
2748:   PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
2749:   PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
2750:   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
2751:   PetscCall(PetscStrcmp(mtype, MATIS, &isMatIS));
2752:   if (!isShell) {
2753:     // There are three states with pblocks, since block starts can have no dofs:
2754:     // UNKNOWN) New Block:   An open block has been signalled by pblocks[p] == 1
2755:     // TRUE)    Block Start: The first entry in a block has been added
2756:     // FALSE)   Block Add:   An additional block entry has been added, since pblocks[p] == 0
2757:     PetscBT         blst;
2758:     PetscBool3      bstate     = PETSC_BOOL3_UNKNOWN;
2759:     PetscBool       fillMatrix = (PetscBool)(!dm->prealloc_only && !isMatIS);
2760:     const PetscInt *perm       = NULL;
2761:     PetscInt       *dnz, *onz, *dnzu, *onzu, bsLocal[2], bsMinMax[2], *pblocks;
2762:     PetscInt        pStart, pEnd, dof, cdof, num_fields;

2764:     PetscCall(DMGetLocalToGlobalMapping(dm, &ltog));
2765:     PetscCall(PetscSectionGetBlockStarts(sectionLocal, &blst));
2766:     if (sectionLocal->perm) PetscCall(ISGetIndices(sectionLocal->perm, &perm));

2768:     PetscCall(PetscCalloc1(localSize, &pblocks));
2769:     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
2770:     PetscCall(PetscSectionGetNumFields(sectionGlobal, &num_fields));
2771:     // We need to process in the permuted order to get block sizes right
2772:     for (PetscInt point = pStart; point < pEnd; ++point) {
2773:       const PetscInt p = perm ? perm[point] : point;

2775:       switch (dm->blocking_type) {
2776:       case DM_BLOCKING_TOPOLOGICAL_POINT: { // One block per topological point
2777:         PetscInt bdof, offset;

2779:         PetscCall(PetscSectionGetDof(sectionGlobal, p, &dof));
2780:         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &offset));
2781:         PetscCall(PetscSectionGetConstraintDof(sectionGlobal, p, &cdof));
2782:         if (blst && PetscBTLookup(blst, p)) bstate = PETSC_BOOL3_UNKNOWN;
2783:         if (dof > 0) {
2784:           // State change
2785:           if (bstate == PETSC_BOOL3_UNKNOWN) bstate = PETSC_BOOL3_TRUE;
2786:           else if (bstate == PETSC_BOOL3_TRUE && blst && !PetscBTLookup(blst, p)) bstate = PETSC_BOOL3_FALSE;

2788:           for (PetscInt i = 0; i < dof - cdof; ++i) pblocks[offset - localStart + i] = dof - cdof;
2789:           // Signal block concatenation
2790:           if (bstate == PETSC_BOOL3_FALSE && dof - cdof) pblocks[offset - localStart] = -(dof - cdof);
2791:         }
2792:         dof  = dof < 0 ? -(dof + 1) : dof;
2793:         bdof = cdof && (dof - cdof) ? 1 : dof;
2794:         if (dof) {
2795:           if (bs < 0) {
2796:             bs = bdof;
2797:           } else if (bs != bdof) {
2798:             bs = 1;
2799:           }
2800:         }
2801:       } break;
2802:       case DM_BLOCKING_FIELD_NODE: {
2803:         for (PetscInt field = 0; field < num_fields; field++) {
2804:           PetscInt num_comp, bdof, offset;
2805:           PetscCall(PetscSectionGetFieldComponents(sectionGlobal, field, &num_comp));
2806:           PetscCall(PetscSectionGetFieldDof(sectionGlobal, p, field, &dof));
2807:           if (dof < 0) continue;
2808:           PetscCall(PetscSectionGetFieldOffset(sectionGlobal, p, field, &offset));
2809:           PetscCall(PetscSectionGetFieldConstraintDof(sectionGlobal, p, field, &cdof));
2810:           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);
2811:           PetscInt num_nodes = dof / num_comp;
2812:           for (PetscInt i = 0; i < dof - cdof; i++) pblocks[offset - localStart + i] = (dof - cdof) / num_nodes;
2813:           // Handle possibly constant block size (unlikely)
2814:           bdof = cdof && (dof - cdof) ? 1 : dof;
2815:           if (dof) {
2816:             if (bs < 0) {
2817:               bs = bdof;
2818:             } else if (bs != bdof) {
2819:               bs = 1;
2820:             }
2821:           }
2822:         }
2823:       } break;
2824:       }
2825:     }
2826:     if (sectionLocal->perm) PetscCall(ISRestoreIndices(sectionLocal->perm, &perm));
2827:     /* Must have same blocksize on all procs (some might have no points) */
2828:     bsLocal[0] = bs < 0 ? PETSC_INT_MAX : bs;
2829:     bsLocal[1] = bs;
2830:     PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax));
2831:     if (bsMinMax[0] != bsMinMax[1]) bs = 1;
2832:     else bs = bsMinMax[0];
2833:     bs = PetscMax(1, bs);
2834:     PetscCall(MatSetLocalToGlobalMapping(*J, ltog, ltog));
2835:     if (dm->prealloc_skip) { // User will likely use MatSetPreallocationCOO(), but still set structural parameters
2836:       PetscCall(MatSetBlockSize(*J, bs));
2837:       PetscCall(MatSetUp(*J));
2838:     } else {
2839:       PetscCall(PetscCalloc4(localSize / bs, &dnz, localSize / bs, &onz, localSize / bs, &dnzu, localSize / bs, &onzu));
2840:       PetscCall(DMPlexPreallocateOperator(dm, bs, dnz, onz, dnzu, onzu, *J, fillMatrix));
2841:       PetscCall(PetscFree4(dnz, onz, dnzu, onzu));
2842:     }
2843:     if (pblocks) { // Consolidate blocks
2844:       PetscInt nblocks = 0;
2845:       pblocks[0]       = PetscAbs(pblocks[0]);
2846:       for (PetscInt i = 0; i < localSize; i += PetscMax(1, pblocks[i])) {
2847:         if (pblocks[i] == 0) continue;
2848:         // Negative block size indicates the blocks should be concatenated
2849:         if (pblocks[i] < 0) {
2850:           pblocks[i] = -pblocks[i];
2851:           pblocks[nblocks - 1] += pblocks[i];
2852:         } else {
2853:           pblocks[nblocks++] = pblocks[i]; // nblocks always <= i
2854:         }
2855:         for (PetscInt j = 1; j < pblocks[i]; j++)
2856:           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);
2857:       }
2858:       PetscCall(MatSetVariableBlockSizes(*J, nblocks, pblocks));
2859:     }
2860:     PetscCall(PetscFree(pblocks));
2861:   }
2862:   PetscCall(MatSetDM(*J, dm));
2863:   PetscFunctionReturn(PETSC_SUCCESS);
2864: }

2866: /*@
2867:   DMPlexGetSubdomainSection - Returns the section associated with the subdomain

2869:   Not Collective

2871:   Input Parameter:
2872: . dm - The `DMPLEX`

2874:   Output Parameter:
2875: . subsection - The subdomain section

2877:   Level: developer

2879: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSection`
2880: @*/
2881: PetscErrorCode DMPlexGetSubdomainSection(DM dm, PetscSection *subsection)
2882: {
2883:   DM_Plex *mesh = (DM_Plex *)dm->data;

2885:   PetscFunctionBegin;
2887:   if (!mesh->subdomainSection) {
2888:     PetscSection section;
2889:     PetscSF      sf;

2891:     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sf));
2892:     PetscCall(DMGetLocalSection(dm, &section));
2893:     PetscCall(PetscSectionCreateGlobalSection(section, sf, PETSC_TRUE, PETSC_FALSE, PETSC_TRUE, &mesh->subdomainSection));
2894:     PetscCall(PetscSFDestroy(&sf));
2895:   }
2896:   *subsection = mesh->subdomainSection;
2897:   PetscFunctionReturn(PETSC_SUCCESS);
2898: }

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

2903:   Not Collective

2905:   Input Parameter:
2906: . dm - The `DMPLEX`

2908:   Output Parameters:
2909: + pStart - The first mesh point
2910: - pEnd   - The upper bound for mesh points

2912:   Level: beginner

2914: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetChart()`
2915: @*/
2916: PetscErrorCode DMPlexGetChart(DM dm, PetscInt *pStart, PetscInt *pEnd)
2917: {
2918:   DM_Plex *mesh = (DM_Plex *)dm->data;

2920:   PetscFunctionBegin;
2922:   if (mesh->tr) PetscCall(DMPlexTransformGetChart(mesh->tr, pStart, pEnd));
2923:   else PetscCall(PetscSectionGetChart(mesh->coneSection, pStart, pEnd));
2924:   PetscFunctionReturn(PETSC_SUCCESS);
2925: }

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

2930:   Not Collective

2932:   Input Parameters:
2933: + dm     - The `DMPLEX`
2934: . pStart - The first mesh point
2935: - pEnd   - The upper bound for mesh points

2937:   Level: beginner

2939: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetChart()`
2940: @*/
2941: PetscErrorCode DMPlexSetChart(DM dm, PetscInt pStart, PetscInt pEnd)
2942: {
2943:   DM_Plex *mesh = (DM_Plex *)dm->data;

2945:   PetscFunctionBegin;
2947:   PetscCall(PetscSectionSetChart(mesh->coneSection, pStart, pEnd));
2948:   PetscCall(PetscSectionSetChart(mesh->supportSection, pStart, pEnd));
2949:   PetscCall(PetscFree(mesh->cellTypes));
2950:   PetscFunctionReturn(PETSC_SUCCESS);
2951: }

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

2956:   Not Collective

2958:   Input Parameters:
2959: + dm - The `DMPLEX`
2960: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

2962:   Output Parameter:
2963: . size - The cone size for point `p`

2965:   Level: beginner

2967: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`
2968: @*/
2969: PetscErrorCode DMPlexGetConeSize(DM dm, PetscInt p, PetscInt *size)
2970: {
2971:   DM_Plex *mesh = (DM_Plex *)dm->data;

2973:   PetscFunctionBegin;
2975:   PetscAssertPointer(size, 3);
2976:   if (mesh->tr) PetscCall(DMPlexTransformGetConeSize(mesh->tr, p, size));
2977:   else PetscCall(PetscSectionGetDof(mesh->coneSection, p, size));
2978:   PetscFunctionReturn(PETSC_SUCCESS);
2979: }

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

2984:   Not Collective

2986:   Input Parameters:
2987: + dm   - The `DMPLEX`
2988: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
2989: - size - The cone size for point `p`

2991:   Level: beginner

2993:   Note:
2994:   This should be called after `DMPlexSetChart()`.

2996: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetCone()`, `DMPlexCreate()`, `DMPlexGetConeSize()`, `DMPlexSetChart()`
2997: @*/
2998: PetscErrorCode DMPlexSetConeSize(DM dm, PetscInt p, PetscInt size)
2999: {
3000:   DM_Plex *mesh = (DM_Plex *)dm->data;

3002:   PetscFunctionBegin;
3004:   PetscCheck(!mesh->tr, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot call DMPlexSetConeSize() on a mesh with a transform defined.");
3005:   PetscCall(PetscSectionSetDof(mesh->coneSection, p, size));
3006:   PetscFunctionReturn(PETSC_SUCCESS);
3007: }

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

3012:   Not Collective

3014:   Input Parameters:
3015: + dm - The `DMPLEX`
3016: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

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

3021:   Level: beginner

3023:   Fortran Notes:
3024:   `cone` must be declared with
3025: .vb
3026:   PetscInt, pointer :: cone(:)
3027: .ve

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

3032: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSize()`, `DMPlexSetCone()`, `DMPlexGetConeTuple()`, `DMPlexSetChart()`, `DMPlexRestoreCone()`
3033: @*/
3034: PetscErrorCode DMPlexGetCone(DM dm, PetscInt p, const PetscInt *cone[])
3035: {
3036:   DM_Plex *mesh = (DM_Plex *)dm->data;
3037:   PetscInt off;

3039:   PetscFunctionBegin;
3041:   PetscAssertPointer(cone, 3);
3042:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3043:   *cone = PetscSafePointerPlusOffset(mesh->cones, off);
3044:   PetscFunctionReturn(PETSC_SUCCESS);
3045: }

3047: /*@
3048:   DMPlexGetConeTuple - Return the points on the in-edges of several points in the DAG

3050:   Not Collective

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

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

3060:   Level: intermediate

3062: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeRecursive()`, `DMPlexSetChart()`, `PetscSection`, `IS`
3063: @*/
3064: PetscErrorCode DMPlexGetConeTuple(DM dm, IS p, PetscSection *pConesSection, IS *pCones)
3065: {
3066:   PetscSection cs, newcs;
3067:   PetscInt    *cones;
3068:   PetscInt    *newarr = NULL;
3069:   PetscInt     n;

3071:   PetscFunctionBegin;
3072:   PetscCall(DMPlexGetCones(dm, &cones));
3073:   PetscCall(DMPlexGetConeSection(dm, &cs));
3074:   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_INT, cones, p, &newcs, pCones ? ((void **)&newarr) : NULL));
3075:   if (pConesSection) *pConesSection = newcs;
3076:   if (pCones) {
3077:     PetscCall(PetscSectionGetStorageSize(newcs, &n));
3078:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)p), n, newarr, PETSC_OWN_POINTER, pCones));
3079:   }
3080:   PetscFunctionReturn(PETSC_SUCCESS);
3081: }

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

3086:   Not Collective

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

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

3095:   Level: advanced

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

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

3102: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeTuple()`, `DMPlexGetConeRecursive()`, `DMPlexRestoreConeRecursive()`,
3103:           `DMPlexGetDepth()`, `IS`
3104: @*/
3105: PetscErrorCode DMPlexGetConeRecursiveVertices(DM dm, IS points, IS *expandedPoints)
3106: {
3107:   IS      *expandedPointsAll;
3108:   PetscInt depth;

3110:   PetscFunctionBegin;
3113:   PetscAssertPointer(expandedPoints, 3);
3114:   PetscCall(DMPlexGetConeRecursive(dm, points, &depth, &expandedPointsAll, NULL));
3115:   *expandedPoints = expandedPointsAll[0];
3116:   PetscCall(PetscObjectReference((PetscObject)expandedPointsAll[0]));
3117:   PetscCall(DMPlexRestoreConeRecursive(dm, points, &depth, &expandedPointsAll, NULL));
3118:   PetscFunctionReturn(PETSC_SUCCESS);
3119: }

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

3125:   Not Collective

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

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

3136:   Level: advanced

3138:   Notes:
3139:   Like `DMPlexGetConeTuple()` but recursive.

3141:   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.
3142:   For example, for d=0 it contains only vertices, for d=1 it can contain vertices and edges, etc.

3144:   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\:
3145:   (1) DAG points in `expandedPoints`[d+1] with `depth` d+1 to their cone points in `expandedPoints`[d];
3146:   (2) DAG points in `expandedPoints`[d+1] with `depth` in [0,d] to the same points in `expandedPoints`[d].

3148: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeTuple()`, `DMPlexRestoreConeRecursive()`, `DMPlexGetConeRecursiveVertices()`,
3149:           `DMPlexGetDepth()`, `PetscSection`, `IS`
3150: @*/
3151: PetscErrorCode DMPlexGetConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
3152: {
3153:   const PetscInt *arr0 = NULL, *cone = NULL;
3154:   PetscInt       *arr = NULL, *newarr = NULL;
3155:   PetscInt        d, depth_, i, n, newn, cn, co, start, end;
3156:   IS             *expandedPoints_;
3157:   PetscSection   *sections_;

3159:   PetscFunctionBegin;
3162:   if (depth) PetscAssertPointer(depth, 3);
3163:   if (expandedPoints) PetscAssertPointer(expandedPoints, 4);
3164:   if (sections) PetscAssertPointer(sections, 5);
3165:   PetscCall(ISGetLocalSize(points, &n));
3166:   PetscCall(ISGetIndices(points, &arr0));
3167:   PetscCall(DMPlexGetDepth(dm, &depth_));
3168:   PetscCall(PetscCalloc1(depth_, &expandedPoints_));
3169:   PetscCall(PetscCalloc1(depth_, &sections_));
3170:   arr = (PetscInt *)arr0; /* this is ok because first generation of arr is not modified */
3171:   for (d = depth_ - 1; d >= 0; d--) {
3172:     PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &sections_[d]));
3173:     PetscCall(PetscSectionSetChart(sections_[d], 0, n));
3174:     for (i = 0; i < n; i++) {
3175:       PetscCall(DMPlexGetDepthStratum(dm, d + 1, &start, &end));
3176:       if (arr[i] >= start && arr[i] < end) {
3177:         PetscCall(DMPlexGetConeSize(dm, arr[i], &cn));
3178:         PetscCall(PetscSectionSetDof(sections_[d], i, cn));
3179:       } else {
3180:         PetscCall(PetscSectionSetDof(sections_[d], i, 1));
3181:       }
3182:     }
3183:     PetscCall(PetscSectionSetUp(sections_[d]));
3184:     PetscCall(PetscSectionGetStorageSize(sections_[d], &newn));
3185:     PetscCall(PetscMalloc1(newn, &newarr));
3186:     for (i = 0; i < n; i++) {
3187:       PetscCall(PetscSectionGetDof(sections_[d], i, &cn));
3188:       PetscCall(PetscSectionGetOffset(sections_[d], i, &co));
3189:       if (cn > 1) {
3190:         PetscCall(DMPlexGetCone(dm, arr[i], &cone));
3191:         PetscCall(PetscMemcpy(&newarr[co], cone, cn * sizeof(PetscInt)));
3192:       } else {
3193:         newarr[co] = arr[i];
3194:       }
3195:     }
3196:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newarr, PETSC_OWN_POINTER, &expandedPoints_[d]));
3197:     arr = newarr;
3198:     n   = newn;
3199:   }
3200:   PetscCall(ISRestoreIndices(points, &arr0));
3201:   *depth = depth_;
3202:   if (expandedPoints) *expandedPoints = expandedPoints_;
3203:   else {
3204:     for (d = 0; d < depth_; d++) PetscCall(ISDestroy(&expandedPoints_[d]));
3205:     PetscCall(PetscFree(expandedPoints_));
3206:   }
3207:   if (sections) *sections = sections_;
3208:   else {
3209:     for (d = 0; d < depth_; d++) PetscCall(PetscSectionDestroy(&sections_[d]));
3210:     PetscCall(PetscFree(sections_));
3211:   }
3212:   PetscFunctionReturn(PETSC_SUCCESS);
3213: }

3215: /*@
3216:   DMPlexRestoreConeRecursive - Deallocates arrays created by `DMPlexGetConeRecursive()`

3218:   Not Collective

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

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

3229:   Level: advanced

3231:   Note:
3232:   See `DMPlexGetConeRecursive()`

3234: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexGetConeTuple()`, `DMPlexGetConeRecursive()`, `DMPlexGetConeRecursiveVertices()`,
3235:           `DMPlexGetDepth()`, `IS`, `PetscSection`
3236: @*/
3237: PetscErrorCode DMPlexRestoreConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
3238: {
3239:   PetscInt d, depth_;

3241:   PetscFunctionBegin;
3242:   PetscCall(DMPlexGetDepth(dm, &depth_));
3243:   PetscCheck(!depth || *depth == depth_, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "depth changed since last call to DMPlexGetConeRecursive");
3244:   if (depth) *depth = 0;
3245:   if (expandedPoints) {
3246:     for (d = 0; d < depth_; d++) PetscCall(ISDestroy(&(*expandedPoints)[d]));
3247:     PetscCall(PetscFree(*expandedPoints));
3248:   }
3249:   if (sections) {
3250:     for (d = 0; d < depth_; d++) PetscCall(PetscSectionDestroy(&(*sections)[d]));
3251:     PetscCall(PetscFree(*sections));
3252:   }
3253:   PetscFunctionReturn(PETSC_SUCCESS);
3254: }

3256: /*@
3257:   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

3259:   Not Collective

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

3266:   Level: beginner

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

3271: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`, `DMPlexSetSupport()`, `DMPlexSetSupportSize()`
3272: @*/
3273: PetscErrorCode DMPlexSetCone(DM dm, PetscInt p, const PetscInt cone[])
3274: {
3275:   DM_Plex *mesh = (DM_Plex *)dm->data;
3276:   PetscInt dof, off, c;

3278:   PetscFunctionBegin;
3280:   PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3281:   if (dof) PetscAssertPointer(cone, 3);
3282:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3283:   if (PetscDefined(USE_DEBUG)) {
3284:     PetscInt pStart, pEnd;
3285:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3286:     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);
3287:     for (c = 0; c < dof; ++c) {
3288:       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);
3289:       mesh->cones[off + c] = cone[c];
3290:     }
3291:   } else {
3292:     for (c = 0; c < dof; ++c) mesh->cones[off + c] = cone[c];
3293:   }
3294:   PetscFunctionReturn(PETSC_SUCCESS);
3295: }

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

3300:   Not Collective

3302:   Input Parameters:
3303: + dm - The `DMPLEX`
3304: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

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

3310:   Level: beginner

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

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

3322: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetConeSize()`, `DMPolytopeTypeComposeOrientation()`, `DMPolytopeTypeComposeOrientationInv()`,
3323:           `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetCone()`, `DMPlexSetChart()`
3324: @*/
3325: PetscErrorCode DMPlexGetConeOrientation(DM dm, PetscInt p, const PetscInt *coneOrientation[])
3326: {
3327:   DM_Plex *mesh = (DM_Plex *)dm->data;
3328:   PetscInt off;

3330:   PetscFunctionBegin;
3332:   if (PetscDefined(USE_DEBUG)) {
3333:     PetscInt dof;
3334:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3335:     if (dof) PetscAssertPointer(coneOrientation, 3);
3336:   }
3337:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));

3339:   *coneOrientation = &mesh->coneOrientations[off];
3340:   PetscFunctionReturn(PETSC_SUCCESS);
3341: }

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

3346:   Not Collective

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

3353:   Level: beginner

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

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

3360: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetConeOrientation()`, `DMPlexSetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3361: @*/
3362: PetscErrorCode DMPlexSetConeOrientation(DM dm, PetscInt p, const PetscInt coneOrientation[])
3363: {
3364:   DM_Plex *mesh = (DM_Plex *)dm->data;
3365:   PetscInt pStart, pEnd;
3366:   PetscInt dof, off, c;

3368:   PetscFunctionBegin;
3370:   PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3371:   if (dof) PetscAssertPointer(coneOrientation, 3);
3372:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3373:   if (PetscDefined(USE_DEBUG)) {
3374:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3375:     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);
3376:     for (c = 0; c < dof; ++c) {
3377:       PetscInt cdof, o = coneOrientation[c];

3379:       PetscCall(PetscSectionGetDof(mesh->coneSection, mesh->cones[off + c], &cdof));
3380:       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);
3381:       mesh->coneOrientations[off + c] = o;
3382:     }
3383:   } else {
3384:     for (c = 0; c < dof; ++c) mesh->coneOrientations[off + c] = coneOrientation[c];
3385:   }
3386:   PetscFunctionReturn(PETSC_SUCCESS);
3387: }

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

3392:   Not Collective

3394:   Input Parameters:
3395: + dm        - The `DMPLEX`
3396: . p         - The point, which must lie in the chart set with `DMPlexSetChart()`
3397: . conePos   - The local index in the cone where the point should be put
3398: - conePoint - The mesh point to insert

3400:   Level: beginner

3402: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3403: @*/
3404: PetscErrorCode DMPlexInsertCone(DM dm, PetscInt p, PetscInt conePos, PetscInt conePoint)
3405: {
3406:   DM_Plex *mesh = (DM_Plex *)dm->data;
3407:   PetscInt pStart, pEnd;
3408:   PetscInt dof, off;

3410:   PetscFunctionBegin;
3412:   if (PetscDefined(USE_DEBUG)) {
3413:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3414:     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);
3415:     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);
3416:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3417:     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);
3418:   }
3419:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3420:   mesh->cones[off + conePos] = conePoint;
3421:   PetscFunctionReturn(PETSC_SUCCESS);
3422: }

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

3427:   Not Collective

3429:   Input Parameters:
3430: + dm              - The `DMPLEX`
3431: . p               - The point, which must lie in the chart set with `DMPlexSetChart()`
3432: . conePos         - The local index in the cone where the point should be put
3433: - coneOrientation - The point orientation to insert

3435:   Level: beginner

3437:   Note:
3438:   The meaning of coneOrientation values is detailed in `DMPlexGetConeOrientation()`.

3440: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3441: @*/
3442: PetscErrorCode DMPlexInsertConeOrientation(DM dm, PetscInt p, PetscInt conePos, PetscInt coneOrientation)
3443: {
3444:   DM_Plex *mesh = (DM_Plex *)dm->data;
3445:   PetscInt pStart, pEnd;
3446:   PetscInt dof, off;

3448:   PetscFunctionBegin;
3450:   if (PetscDefined(USE_DEBUG)) {
3451:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
3452:     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);
3453:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3454:     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);
3455:   }
3456:   PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3457:   mesh->coneOrientations[off + conePos] = coneOrientation;
3458:   PetscFunctionReturn(PETSC_SUCCESS);
3459: }

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

3464:   Not collective

3466:   Input Parameters:
3467: + dm - The DMPlex
3468: - p  - The point, which must lie in the chart set with DMPlexSetChart()

3470:   Output Parameters:
3471: + cone - An array of points which are on the in-edges for point `p`
3472: - ornt - An array of orientations which are on the in-edges for point `p`. An orientation is an
3473:          integer giving the prescription for cone traversal.

3475:   Level: beginner

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

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

3485:   Fortran Notes:
3486:   `cone` and `ornt` must be declared with
3487: .vb
3488:   PetscInt, pointer :: cone(:)
3489:   PetscInt, pointer :: ornt(:)
3490: .ve

3492: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreOrientedCone()`, `DMPlexGetConeSize()`, `DMPlexGetCone()`, `DMPlexGetChart()`
3493: @*/
3494: PetscErrorCode DMPlexGetOrientedCone(DM dm, PetscInt p, const PetscInt *cone[], const PetscInt *ornt[])
3495: {
3496:   DM_Plex *mesh = (DM_Plex *)dm->data;

3498:   PetscFunctionBegin;
3500:   if (mesh->tr) {
3501:     PetscCall(DMPlexTransformGetCone(mesh->tr, p, cone, ornt));
3502:   } else {
3503:     PetscInt off;
3504:     if (PetscDefined(USE_DEBUG)) {
3505:       PetscInt dof;
3506:       PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
3507:       if (dof) {
3508:         if (cone) PetscAssertPointer(cone, 3);
3509:         if (ornt) PetscAssertPointer(ornt, 4);
3510:       }
3511:     }
3512:     PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
3513:     if (cone) *cone = PetscSafePointerPlusOffset(mesh->cones, off);
3514:     if (ornt) *ornt = PetscSafePointerPlusOffset(mesh->coneOrientations, off);
3515:   }
3516:   PetscFunctionReturn(PETSC_SUCCESS);
3517: }

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

3522:   Not Collective

3524:   Input Parameters:
3525: + dm   - The DMPlex
3526: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
3527: . cone - An array of points which are on the in-edges for point p
3528: - ornt - An array of orientations which are on the in-edges for point `p`. An orientation is an
3529:          integer giving the prescription for cone traversal.

3531:   Level: beginner

3533: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetOrientedCone()`, `DMPlexGetConeSize()`, `DMPlexGetCone()`, `DMPlexGetChart()`
3534: @*/
3535: PetscErrorCode DMPlexRestoreOrientedCone(DM dm, PetscInt p, const PetscInt *cone[], const PetscInt *ornt[])
3536: {
3537:   DM_Plex *mesh = (DM_Plex *)dm->data;

3539:   PetscFunctionBegin;
3541:   if (mesh->tr) PetscCall(DMPlexTransformRestoreCone(mesh->tr, p, cone, ornt));
3542:   PetscFunctionReturn(PETSC_SUCCESS);
3543: }

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

3548:   Not Collective

3550:   Input Parameters:
3551: + dm - The `DMPLEX`
3552: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

3554:   Output Parameter:
3555: . size - The support size for point `p`

3557:   Level: beginner

3559: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`, `DMPlexGetConeSize()`
3560: @*/
3561: PetscErrorCode DMPlexGetSupportSize(DM dm, PetscInt p, PetscInt *size)
3562: {
3563:   DM_Plex *mesh = (DM_Plex *)dm->data;

3565:   PetscFunctionBegin;
3567:   PetscAssertPointer(size, 3);
3568:   PetscCall(PetscSectionGetDof(mesh->supportSection, p, size));
3569:   PetscFunctionReturn(PETSC_SUCCESS);
3570: }

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

3575:   Not Collective

3577:   Input Parameters:
3578: + dm   - The `DMPLEX`
3579: . p    - The point, which must lie in the chart set with `DMPlexSetChart()`
3580: - size - The support size for point `p`

3582:   Level: beginner

3584:   Note:
3585:   This should be called after `DMPlexSetChart()`.

3587: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetSupportSize()`, `DMPlexSetChart()`
3588: @*/
3589: PetscErrorCode DMPlexSetSupportSize(DM dm, PetscInt p, PetscInt size)
3590: {
3591:   DM_Plex *mesh = (DM_Plex *)dm->data;

3593:   PetscFunctionBegin;
3595:   PetscCall(PetscSectionSetDof(mesh->supportSection, p, size));
3596:   PetscFunctionReturn(PETSC_SUCCESS);
3597: }

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

3602:   Not Collective

3604:   Input Parameters:
3605: + dm - The `DMPLEX`
3606: - p  - The point, which must lie in the chart set with `DMPlexSetChart()`

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

3611:   Level: beginner

3613:   Fortran Notes:
3614:   `support` must be declared with
3615: .vb
3616:   PetscInt, pointer :: support(:)
3617: .ve

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

3622: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetSupportSize()`, `DMPlexSetSupport()`, `DMPlexGetCone()`, `DMPlexSetChart()`
3623: @*/
3624: PetscErrorCode DMPlexGetSupport(DM dm, PetscInt p, const PetscInt *support[])
3625: {
3626:   DM_Plex *mesh = (DM_Plex *)dm->data;
3627:   PetscInt off;

3629:   PetscFunctionBegin;
3631:   PetscAssertPointer(support, 3);
3632:   PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
3633:   *support = PetscSafePointerPlusOffset(mesh->supports, off);
3634:   PetscFunctionReturn(PETSC_SUCCESS);
3635: }

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

3640:   Not Collective

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

3647:   Level: beginner

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

3652: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetCone()`, `DMPlexSetConeSize()`, `DMPlexCreate()`, `DMPlexGetSupport()`, `DMPlexSetChart()`, `DMPlexSetSupportSize()`, `DMSetUp()`
3653: @*/
3654: PetscErrorCode DMPlexSetSupport(DM dm, PetscInt p, const PetscInt support[])
3655: {
3656:   DM_Plex *mesh = (DM_Plex *)dm->data;
3657:   PetscInt pStart, pEnd;
3658:   PetscInt dof, off, c;

3660:   PetscFunctionBegin;
3662:   PetscCall(PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd));
3663:   PetscCall(PetscSectionGetDof(mesh->supportSection, p, &dof));
3664:   if (dof) PetscAssertPointer(support, 3);
3665:   PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
3666:   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);
3667:   for (c = 0; c < dof; ++c) {
3668:     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);
3669:     mesh->supports[off + c] = support[c];
3670:   }
3671:   PetscFunctionReturn(PETSC_SUCCESS);
3672: }

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

3677:   Not Collective

3679:   Input Parameters:
3680: + dm           - The `DMPLEX`
3681: . p            - The point, which must lie in the chart set with `DMPlexSetChart()`
3682: . supportPos   - The local index in the cone where the point should be put
3683: - supportPoint - The mesh point to insert

3685:   Level: beginner

3687: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexGetCone()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMSetUp()`
3688: @*/
3689: PetscErrorCode DMPlexInsertSupport(DM dm, PetscInt p, PetscInt supportPos, PetscInt supportPoint)
3690: {
3691:   DM_Plex *mesh = (DM_Plex *)dm->data;
3692:   PetscInt pStart, pEnd;
3693:   PetscInt dof, off;

3695:   PetscFunctionBegin;
3697:   PetscCall(PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd));
3698:   PetscCall(PetscSectionGetDof(mesh->supportSection, p, &dof));
3699:   PetscCall(PetscSectionGetOffset(mesh->supportSection, p, &off));
3700:   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);
3701:   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);
3702:   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);
3703:   mesh->supports[off + supportPos] = supportPoint;
3704:   PetscFunctionReturn(PETSC_SUCCESS);
3705: }

3707: /* Converts an orientation o in the current numbering to the previous scheme used in Plex */
3708: PetscInt DMPolytopeConvertNewOrientation_Internal(DMPolytopeType ct, PetscInt o)
3709: {
3710:   switch (ct) {
3711:   case DM_POLYTOPE_SEGMENT:
3712:     if (o == -1) return -2;
3713:     break;
3714:   case DM_POLYTOPE_TRIANGLE:
3715:     if (o == -3) return -1;
3716:     if (o == -2) return -3;
3717:     if (o == -1) return -2;
3718:     break;
3719:   case DM_POLYTOPE_QUADRILATERAL:
3720:     if (o == -4) return -2;
3721:     if (o == -3) return -1;
3722:     if (o == -2) return -4;
3723:     if (o == -1) return -3;
3724:     break;
3725:   default:
3726:     return o;
3727:   }
3728:   return o;
3729: }

3731: /* Converts an orientation o in the previous scheme used in Plex to the current numbering */
3732: PetscInt DMPolytopeConvertOldOrientation_Internal(DMPolytopeType ct, PetscInt o)
3733: {
3734:   switch (ct) {
3735:   case DM_POLYTOPE_SEGMENT:
3736:     if ((o == -2) || (o == 1)) return -1;
3737:     if (o == -1) return 0;
3738:     break;
3739:   case DM_POLYTOPE_TRIANGLE:
3740:     if (o == -3) return -2;
3741:     if (o == -2) return -1;
3742:     if (o == -1) return -3;
3743:     break;
3744:   case DM_POLYTOPE_QUADRILATERAL:
3745:     if (o == -4) return -2;
3746:     if (o == -3) return -1;
3747:     if (o == -2) return -4;
3748:     if (o == -1) return -3;
3749:     break;
3750:   default:
3751:     return o;
3752:   }
3753:   return o;
3754: }

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

3761:   PetscFunctionBegin;
3762:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
3763:   for (p = pStart; p < pEnd; ++p) {
3764:     const PetscInt *cone, *ornt;
3765:     PetscInt        coneSize, c;

3767:     PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
3768:     PetscCall(DMPlexGetCone(dm, p, &cone));
3769:     PetscCall(DMPlexGetConeOrientation(dm, p, &ornt));
3770:     for (c = 0; c < coneSize; ++c) {
3771:       DMPolytopeType ct;
3772:       const PetscInt o = ornt[c];

3774:       PetscCall(DMPlexGetCellType(dm, cone[c], &ct));
3775:       switch (ct) {
3776:       case DM_POLYTOPE_SEGMENT:
3777:         if ((o == -2) || (o == 1)) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -1));
3778:         if (o == -1) PetscCall(DMPlexInsertConeOrientation(dm, p, c, 0));
3779:         break;
3780:       case DM_POLYTOPE_TRIANGLE:
3781:         if (o == -3) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -2));
3782:         if (o == -2) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -1));
3783:         if (o == -1) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -3));
3784:         break;
3785:       case DM_POLYTOPE_QUADRILATERAL:
3786:         if (o == -4) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -2));
3787:         if (o == -3) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -1));
3788:         if (o == -2) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -4));
3789:         if (o == -1) PetscCall(DMPlexInsertConeOrientation(dm, p, c, -3));
3790:         break;
3791:       default:
3792:         break;
3793:       }
3794:     }
3795:   }
3796:   PetscFunctionReturn(PETSC_SUCCESS);
3797: }

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

3803:   PetscFunctionBeginHot;
3804:   if (PetscDefined(USE_DEBUG) || mesh->tr) {
3805:     if (useCone) {
3806:       PetscCall(DMPlexGetConeSize(dm, p, size));
3807:       PetscCall(DMPlexGetOrientedCone(dm, p, arr, ornt));
3808:     } else {
3809:       PetscCall(DMPlexGetSupportSize(dm, p, size));
3810:       PetscCall(DMPlexGetSupport(dm, p, arr));
3811:     }
3812:   } else {
3813:     if (useCone) {
3814:       const PetscSection s   = mesh->coneSection;
3815:       const PetscInt     ps  = p - s->pStart;
3816:       const PetscInt     off = s->atlasOff[ps];

3818:       *size = s->atlasDof[ps];
3819:       *arr  = mesh->cones + off;
3820:       *ornt = mesh->coneOrientations + off;
3821:     } else {
3822:       const PetscSection s   = mesh->supportSection;
3823:       const PetscInt     ps  = p - s->pStart;
3824:       const PetscInt     off = s->atlasOff[ps];

3826:       *size = s->atlasDof[ps];
3827:       *arr  = mesh->supports + off;
3828:     }
3829:   }
3830:   PetscFunctionReturn(PETSC_SUCCESS);
3831: }

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

3837:   PetscFunctionBeginHot;
3838:   if (PetscDefined(USE_DEBUG) || mesh->tr) {
3839:     if (useCone) PetscCall(DMPlexRestoreOrientedCone(dm, p, arr, ornt));
3840:   }
3841:   PetscFunctionReturn(PETSC_SUCCESS);
3842: }

3844: static PetscErrorCode DMPlexGetTransitiveClosure_Depth1_Private(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
3845: {
3846:   DMPolytopeType  ct = DM_POLYTOPE_UNKNOWN;
3847:   PetscInt       *closure;
3848:   const PetscInt *tmp = NULL, *tmpO = NULL;
3849:   PetscInt        off = 0, tmpSize, t;

3851:   PetscFunctionBeginHot;
3852:   if (ornt) {
3853:     PetscCall(DMPlexGetCellType(dm, p, &ct));
3854:     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;
3855:   }
3856:   if (*points) {
3857:     closure = *points;
3858:   } else {
3859:     PetscInt maxConeSize, maxSupportSize;
3860:     PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
3861:     PetscCall(DMGetWorkArray(dm, 2 * (PetscMax(maxConeSize, maxSupportSize) + 1), MPIU_INT, &closure));
3862:   }
3863:   PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, p, useCone, &tmpSize, &tmp, &tmpO));
3864:   if (ct == DM_POLYTOPE_UNKNOWN) {
3865:     closure[off++] = p;
3866:     closure[off++] = 0;
3867:     for (t = 0; t < tmpSize; ++t) {
3868:       closure[off++] = tmp[t];
3869:       closure[off++] = tmpO ? tmpO[t] : 0;
3870:     }
3871:   } else {
3872:     const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, ornt);

3874:     /* We assume that cells with a valid type have faces with a valid type */
3875:     closure[off++] = p;
3876:     closure[off++] = ornt;
3877:     for (t = 0; t < tmpSize; ++t) {
3878:       DMPolytopeType ft;

3880:       PetscCall(DMPlexGetCellType(dm, tmp[t], &ft));
3881:       closure[off++] = tmp[arr[t]];
3882:       closure[off++] = tmpO ? DMPolytopeTypeComposeOrientation(ft, ornt, tmpO[t]) : 0;
3883:     }
3884:   }
3885:   PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, p, useCone, &tmpSize, &tmp, &tmpO));
3886:   if (numPoints) *numPoints = tmpSize + 1;
3887:   if (points) *points = closure;
3888:   PetscFunctionReturn(PETSC_SUCCESS);
3889: }

3891: /* We need a special tensor version because we want to allow duplicate points in the endcaps for hybrid cells */
3892: static PetscErrorCode DMPlexTransitiveClosure_Tensor_Internal(DM dm, PetscInt point, DMPolytopeType ct, PetscInt o, PetscBool useCone, PetscInt *numPoints, PetscInt **points)
3893: {
3894:   const PetscInt *arr = DMPolytopeTypeGetArrangement(ct, o);
3895:   const PetscInt *cone, *ornt;
3896:   PetscInt       *pts, *closure = NULL;
3897:   DMPolytopeType  ft;
3898:   PetscInt        maxConeSize, maxSupportSize, coneSeries, supportSeries, maxSize;
3899:   PetscInt        dim, coneSize, c, d, clSize, cl;

3901:   PetscFunctionBeginHot;
3902:   PetscCall(DMGetDimension(dm, &dim));
3903:   PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, point, PETSC_TRUE, &coneSize, &cone, &ornt));
3904:   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
3905:   coneSeries    = (maxConeSize > 1) ? ((PetscPowInt(maxConeSize, dim + 1) - 1) / (maxConeSize - 1)) : dim + 1;
3906:   supportSeries = (maxSupportSize > 1) ? ((PetscPowInt(maxSupportSize, dim + 1) - 1) / (maxSupportSize - 1)) : dim + 1;
3907:   maxSize       = PetscMax(coneSeries, supportSeries);
3908:   if (*points) {
3909:     pts = *points;
3910:   } else PetscCall(DMGetWorkArray(dm, 2 * maxSize, MPIU_INT, &pts));
3911:   c        = 0;
3912:   pts[c++] = point;
3913:   pts[c++] = o;
3914:   PetscCall(DMPlexGetCellType(dm, cone[arr[0 * 2 + 0]], &ft));
3915:   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, cone[arr[0 * 2 + 0]], DMPolytopeTypeComposeOrientation(ft, arr[0 * 2 + 1], ornt[0]), useCone, &clSize, &closure));
3916:   for (cl = 0; cl < clSize * 2; cl += 2) {
3917:     pts[c++] = closure[cl];
3918:     pts[c++] = closure[cl + 1];
3919:   }
3920:   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, cone[arr[1 * 2 + 0]], DMPolytopeTypeComposeOrientation(ft, arr[1 * 2 + 1], ornt[1]), useCone, &clSize, &closure));
3921:   for (cl = 0; cl < clSize * 2; cl += 2) {
3922:     pts[c++] = closure[cl];
3923:     pts[c++] = closure[cl + 1];
3924:   }
3925:   PetscCall(DMPlexRestoreTransitiveClosure(dm, cone[0], useCone, &clSize, &closure));
3926:   for (d = 2; d < coneSize; ++d) {
3927:     PetscCall(DMPlexGetCellType(dm, cone[arr[d * 2 + 0]], &ft));
3928:     pts[c++] = cone[arr[d * 2 + 0]];
3929:     pts[c++] = DMPolytopeTypeComposeOrientation(ft, arr[d * 2 + 1], ornt[d]);
3930:   }
3931:   PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, point, PETSC_TRUE, &coneSize, &cone, &ornt));
3932:   if (dim >= 3) {
3933:     for (d = 2; d < coneSize; ++d) {
3934:       const PetscInt  fpoint = cone[arr[d * 2 + 0]];
3935:       const PetscInt *fcone, *fornt;
3936:       PetscInt        fconeSize, fc, i;

3938:       PetscCall(DMPlexGetCellType(dm, fpoint, &ft));
3939:       const PetscInt *farr = DMPolytopeTypeGetArrangement(ft, DMPolytopeTypeComposeOrientation(ft, arr[d * 2 + 1], ornt[d]));
3940:       PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, fpoint, PETSC_TRUE, &fconeSize, &fcone, &fornt));
3941:       for (fc = 0; fc < fconeSize; ++fc) {
3942:         const PetscInt cp = fcone[farr[fc * 2 + 0]];
3943:         const PetscInt co = farr[fc * 2 + 1];

3945:         for (i = 0; i < c; i += 2)
3946:           if (pts[i] == cp) break;
3947:         if (i == c) {
3948:           PetscCall(DMPlexGetCellType(dm, cp, &ft));
3949:           pts[c++] = cp;
3950:           pts[c++] = DMPolytopeTypeComposeOrientation(ft, co, fornt[farr[fc * 2 + 0]]);
3951:         }
3952:       }
3953:       PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, fpoint, PETSC_TRUE, &fconeSize, &fcone, &fornt));
3954:     }
3955:   }
3956:   *numPoints = c / 2;
3957:   *points    = pts;
3958:   PetscFunctionReturn(PETSC_SUCCESS);
3959: }

3961: PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
3962: {
3963:   DMPolytopeType ct;
3964:   PetscInt      *closure, *fifo;
3965:   PetscInt       closureSize = 0, fifoStart = 0, fifoSize = 0;
3966:   PetscInt       maxConeSize, maxSupportSize, coneSeries, supportSeries;
3967:   PetscInt       depth, maxSize;

3969:   PetscFunctionBeginHot;
3970:   PetscCall(DMPlexGetDepth(dm, &depth));
3971:   if (depth == 1) {
3972:     PetscCall(DMPlexGetTransitiveClosure_Depth1_Private(dm, p, ornt, useCone, numPoints, points));
3973:     PetscFunctionReturn(PETSC_SUCCESS);
3974:   }
3975:   PetscCall(DMPlexGetCellType(dm, p, &ct));
3976:   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;
3977:   if (DMPolytopeTypeIsHybrid(ct) && ct != DM_POLYTOPE_POINT_PRISM_TENSOR) {
3978:     PetscCall(DMPlexTransitiveClosure_Tensor_Internal(dm, p, ct, ornt, useCone, numPoints, points));
3979:     PetscFunctionReturn(PETSC_SUCCESS);
3980:   }
3981:   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
3982:   coneSeries    = (maxConeSize > 1) ? ((PetscPowInt(maxConeSize, depth + 1) - 1) / (maxConeSize - 1)) : depth + 1;
3983:   supportSeries = (maxSupportSize > 1) ? ((PetscPowInt(maxSupportSize, depth + 1) - 1) / (maxSupportSize - 1)) : depth + 1;
3984:   maxSize       = PetscMax(coneSeries, supportSeries);
3985:   PetscCall(DMGetWorkArray(dm, 3 * maxSize, MPIU_INT, &fifo));
3986:   if (*points) {
3987:     closure = *points;
3988:   } else PetscCall(DMGetWorkArray(dm, 2 * maxSize, MPIU_INT, &closure));
3989:   closure[closureSize++] = p;
3990:   closure[closureSize++] = ornt;
3991:   fifo[fifoSize++]       = p;
3992:   fifo[fifoSize++]       = ornt;
3993:   fifo[fifoSize++]       = ct;
3994:   /* Should kick out early when depth is reached, rather than checking all vertices for empty cones */
3995:   while (fifoSize - fifoStart) {
3996:     const PetscInt       q    = fifo[fifoStart++];
3997:     const PetscInt       o    = fifo[fifoStart++];
3998:     const DMPolytopeType qt   = (DMPolytopeType)fifo[fifoStart++];
3999:     const PetscInt      *qarr = DMPolytopeTypeGetArrangement(qt, o);
4000:     const PetscInt      *tmp, *tmpO = NULL;
4001:     PetscInt             tmpSize, t;

4003:     if (PetscDefined(USE_DEBUG)) {
4004:       PetscInt nO = DMPolytopeTypeGetNumArrangements(qt) / 2;
4005:       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);
4006:     }
4007:     PetscCall(DMPlexGetTransitiveClosure_Hot_Private(dm, q, useCone, &tmpSize, &tmp, &tmpO));
4008:     for (t = 0; t < tmpSize; ++t) {
4009:       const PetscInt ip = useCone && qarr ? qarr[t * 2] : t;
4010:       const PetscInt io = useCone && qarr ? qarr[t * 2 + 1] : 0;
4011:       const PetscInt cp = tmp[ip];
4012:       PetscCall(DMPlexGetCellType(dm, cp, &ct));
4013:       const PetscInt co = tmpO ? DMPolytopeTypeComposeOrientation(ct, io, tmpO[ip]) : 0;
4014:       PetscInt       c;

4016:       /* Check for duplicate */
4017:       for (c = 0; c < closureSize; c += 2) {
4018:         if (closure[c] == cp) break;
4019:       }
4020:       if (c == closureSize) {
4021:         closure[closureSize++] = cp;
4022:         closure[closureSize++] = co;
4023:         fifo[fifoSize++]       = cp;
4024:         fifo[fifoSize++]       = co;
4025:         fifo[fifoSize++]       = ct;
4026:       }
4027:     }
4028:     PetscCall(DMPlexRestoreTransitiveClosure_Hot_Private(dm, q, useCone, &tmpSize, &tmp, &tmpO));
4029:   }
4030:   PetscCall(DMRestoreWorkArray(dm, 3 * maxSize, MPIU_INT, &fifo));
4031:   if (numPoints) *numPoints = closureSize / 2;
4032:   if (points) *points = closure;
4033:   PetscFunctionReturn(PETSC_SUCCESS);
4034: }

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

4039:   Not Collective

4041:   Input Parameters:
4042: + dm      - The `DMPLEX`
4043: . p       - The mesh point
4044: - useCone - `PETSC_TRUE` for the closure, otherwise return the star

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

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

4054:   Level: beginner

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

4059:   Fortran Notes:
4060:   `points` must be declared with
4061: .vb
4062:   PetscInt, pointer :: points(:)
4063: .ve
4064:   and is always allocated by the function.

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

4068: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreTransitiveClosure()`, `DMPlexCreate()`, `DMPlexSetCone()`, `DMPlexSetChart()`, `DMPlexGetCone()`
4069: @*/
4070: PetscErrorCode DMPlexGetTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
4071: {
4072:   PetscFunctionBeginHot;
4074:   if (numPoints) PetscAssertPointer(numPoints, 4);
4075:   if (points) PetscAssertPointer(points, 5);
4076:   if (PetscDefined(USE_DEBUG)) {
4077:     PetscInt pStart, pEnd;
4078:     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4079:     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);
4080:   }
4081:   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, p, 0, useCone, numPoints, points));
4082:   PetscFunctionReturn(PETSC_SUCCESS);
4083: }

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

4088:   Not Collective

4090:   Input Parameters:
4091: + dm        - The `DMPLEX`
4092: . p         - The mesh point
4093: . useCone   - `PETSC_TRUE` for the closure, otherwise return the star
4094: . numPoints - The number of points in the closure, so points[] is of size 2*`numPoints`
4095: - points    - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...]

4097:   Level: beginner

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

4102: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetTransitiveClosure()`, `DMPlexCreate()`, `DMPlexSetCone()`, `DMPlexSetChart()`, `DMPlexGetCone()`
4103: @*/
4104: PetscErrorCode DMPlexRestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
4105: {
4106:   PetscFunctionBeginHot;
4108:   if (numPoints) *numPoints = 0;
4109:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, points));
4110:   PetscFunctionReturn(PETSC_SUCCESS);
4111: }

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

4116:   Not Collective

4118:   Input Parameter:
4119: . dm - The `DMPLEX`

4121:   Output Parameters:
4122: + maxConeSize    - The maximum number of in-edges
4123: - maxSupportSize - The maximum number of out-edges

4125:   Level: beginner

4127: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`
4128: @*/
4129: PetscErrorCode DMPlexGetMaxSizes(DM dm, PetscInt *maxConeSize, PetscInt *maxSupportSize)
4130: {
4131:   DM_Plex *mesh = (DM_Plex *)dm->data;

4133:   PetscFunctionBegin;
4135:   if (maxConeSize) PetscCall(PetscSectionGetMaxDof(mesh->coneSection, maxConeSize));
4136:   if (maxSupportSize) PetscCall(PetscSectionGetMaxDof(mesh->supportSection, maxSupportSize));
4137:   PetscFunctionReturn(PETSC_SUCCESS);
4138: }

4140: PetscErrorCode DMSetUp_Plex(DM dm)
4141: {
4142:   DM_Plex *mesh = (DM_Plex *)dm->data;
4143:   PetscInt size, maxSupportSize;

4145:   PetscFunctionBegin;
4147:   PetscCall(PetscSectionSetUp(mesh->coneSection));
4148:   PetscCall(PetscSectionGetStorageSize(mesh->coneSection, &size));
4149:   PetscCall(PetscMalloc1(size, &mesh->cones));
4150:   PetscCall(PetscCalloc1(size, &mesh->coneOrientations));
4151:   PetscCall(PetscSectionGetMaxDof(mesh->supportSection, &maxSupportSize));
4152:   if (maxSupportSize) {
4153:     PetscCall(PetscSectionSetUp(mesh->supportSection));
4154:     PetscCall(PetscSectionGetStorageSize(mesh->supportSection, &size));
4155:     PetscCall(PetscMalloc1(size, &mesh->supports));
4156:   }
4157:   PetscFunctionReturn(PETSC_SUCCESS);
4158: }

4160: PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
4161: {
4162:   PetscFunctionBegin;
4163:   if (subdm) PetscCall(DMClone(dm, subdm));
4164:   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, NULL, NULL, is, subdm));
4165:   if (subdm) (*subdm)->useNatural = dm->useNatural;
4166:   if (dm->useNatural && dm->sfMigration) {
4167:     PetscSF sfNatural;

4169:     (*subdm)->sfMigration = dm->sfMigration;
4170:     PetscCall(PetscObjectReference((PetscObject)dm->sfMigration));
4171:     PetscCall(DMPlexCreateGlobalToNaturalSF(*subdm, NULL, (*subdm)->sfMigration, &sfNatural));
4172:     (*subdm)->sfNatural = sfNatural;
4173:   }
4174:   PetscFunctionReturn(PETSC_SUCCESS);
4175: }

4177: PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm)
4178: {
4179:   PetscInt i = 0;

4181:   PetscFunctionBegin;
4182:   PetscCall(DMClone(dms[0], superdm));
4183:   PetscCall(DMCreateSectionSuperDM(dms, len, is, superdm));
4184:   (*superdm)->useNatural = PETSC_FALSE;
4185:   for (i = 0; i < len; i++) {
4186:     if (dms[i]->useNatural && dms[i]->sfMigration) {
4187:       PetscSF sfNatural;

4189:       (*superdm)->sfMigration = dms[i]->sfMigration;
4190:       PetscCall(PetscObjectReference((PetscObject)dms[i]->sfMigration));
4191:       (*superdm)->useNatural = PETSC_TRUE;
4192:       PetscCall(DMPlexCreateGlobalToNaturalSF(*superdm, NULL, (*superdm)->sfMigration, &sfNatural));
4193:       (*superdm)->sfNatural = sfNatural;
4194:       break;
4195:     }
4196:   }
4197:   PetscFunctionReturn(PETSC_SUCCESS);
4198: }

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

4203:   Not Collective

4205:   Input Parameter:
4206: . dm - The `DMPLEX`

4208:   Level: beginner

4210:   Note:
4211:   This should be called after all calls to `DMPlexSetCone()`

4213: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSetChart()`, `DMPlexSetConeSize()`, `DMPlexSetCone()`
4214: @*/
4215: PetscErrorCode DMPlexSymmetrize(DM dm)
4216: {
4217:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4218:   PetscInt *offsets;
4219:   PetscInt  supportSize;
4220:   PetscInt  pStart, pEnd, p;

4222:   PetscFunctionBegin;
4224:   PetscCheck(!mesh->supports, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Supports were already setup in this DMPlex");
4225:   PetscCall(PetscLogEventBegin(DMPLEX_Symmetrize, dm, 0, 0, 0));
4226:   /* Calculate support sizes */
4227:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4228:   for (p = pStart; p < pEnd; ++p) {
4229:     PetscInt dof, off, c;

4231:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
4232:     PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
4233:     for (c = off; c < off + dof; ++c) PetscCall(PetscSectionAddDof(mesh->supportSection, mesh->cones[c], 1));
4234:   }
4235:   PetscCall(PetscSectionSetUp(mesh->supportSection));
4236:   /* Calculate supports */
4237:   PetscCall(PetscSectionGetStorageSize(mesh->supportSection, &supportSize));
4238:   PetscCall(PetscMalloc1(supportSize, &mesh->supports));
4239:   PetscCall(PetscCalloc1(pEnd - pStart, &offsets));
4240:   for (p = pStart; p < pEnd; ++p) {
4241:     PetscInt dof, off, c;

4243:     PetscCall(PetscSectionGetDof(mesh->coneSection, p, &dof));
4244:     PetscCall(PetscSectionGetOffset(mesh->coneSection, p, &off));
4245:     for (c = off; c < off + dof; ++c) {
4246:       const PetscInt q = mesh->cones[c];
4247:       PetscInt       offS;

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

4251:       mesh->supports[offS + offsets[q]] = p;
4252:       ++offsets[q];
4253:     }
4254:   }
4255:   PetscCall(PetscFree(offsets));
4256:   PetscCall(PetscLogEventEnd(DMPLEX_Symmetrize, dm, 0, 0, 0));
4257:   PetscFunctionReturn(PETSC_SUCCESS);
4258: }

4260: static PetscErrorCode DMPlexCreateDepthStratum(DM dm, DMLabel label, PetscInt depth, PetscInt pStart, PetscInt pEnd)
4261: {
4262:   IS stratumIS;

4264:   PetscFunctionBegin;
4265:   if (pStart >= pEnd) PetscFunctionReturn(PETSC_SUCCESS);
4266:   if (PetscDefined(USE_DEBUG)) {
4267:     PetscInt  qStart, qEnd, numLevels, level;
4268:     PetscBool overlap = PETSC_FALSE;
4269:     PetscCall(DMLabelGetNumValues(label, &numLevels));
4270:     for (level = 0; level < numLevels; level++) {
4271:       PetscCall(DMLabelGetStratumBounds(label, level, &qStart, &qEnd));
4272:       if ((pStart >= qStart && pStart < qEnd) || (pEnd > qStart && pEnd <= qEnd)) {
4273:         overlap = PETSC_TRUE;
4274:         break;
4275:       }
4276:     }
4277:     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);
4278:   }
4279:   PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &stratumIS));
4280:   PetscCall(DMLabelSetStratumIS(label, depth, stratumIS));
4281:   PetscCall(ISDestroy(&stratumIS));
4282:   PetscFunctionReturn(PETSC_SUCCESS);
4283: }

4285: static PetscErrorCode DMPlexStratify_CellType_Private(DM dm, DMLabel label)
4286: {
4287:   PetscInt *pMin, *pMax;
4288:   PetscInt  pStart, pEnd;
4289:   PetscInt  dmin = PETSC_INT_MAX, dmax = PETSC_INT_MIN;

4291:   PetscFunctionBegin;
4292:   {
4293:     DMLabel label2;

4295:     PetscCall(DMPlexGetCellTypeLabel(dm, &label2));
4296:     PetscCall(PetscObjectViewFromOptions((PetscObject)label2, NULL, "-ct_view"));
4297:   }
4298:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4299:   for (PetscInt p = pStart; p < pEnd; ++p) {
4300:     DMPolytopeType ct;

4302:     PetscCall(DMPlexGetCellType(dm, p, &ct));
4303:     dmin = PetscMin(DMPolytopeTypeGetDim(ct), dmin);
4304:     dmax = PetscMax(DMPolytopeTypeGetDim(ct), dmax);
4305:   }
4306:   PetscCall(PetscMalloc2(dmax + 1, &pMin, dmax + 1, &pMax));
4307:   for (PetscInt d = dmin; d <= dmax; ++d) {
4308:     pMin[d] = PETSC_INT_MAX;
4309:     pMax[d] = PETSC_INT_MIN;
4310:   }
4311:   for (PetscInt p = pStart; p < pEnd; ++p) {
4312:     DMPolytopeType ct;
4313:     PetscInt       d;

4315:     PetscCall(DMPlexGetCellType(dm, p, &ct));
4316:     d       = DMPolytopeTypeGetDim(ct);
4317:     pMin[d] = PetscMin(p, pMin[d]);
4318:     pMax[d] = PetscMax(p, pMax[d]);
4319:   }
4320:   for (PetscInt d = dmin; d <= dmax; ++d) {
4321:     if (pMin[d] > pMax[d]) continue;
4322:     PetscCall(DMPlexCreateDepthStratum(dm, label, d, pMin[d], pMax[d] + 1));
4323:   }
4324:   PetscCall(PetscFree2(pMin, pMax));
4325:   PetscFunctionReturn(PETSC_SUCCESS);
4326: }

4328: static PetscErrorCode DMPlexStratify_Topological_Private(DM dm, DMLabel label)
4329: {
4330:   PetscInt pStart, pEnd;
4331:   PetscInt numRoots = 0, numLeaves = 0;

4333:   PetscFunctionBegin;
4334:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4335:   {
4336:     /* Initialize roots and count leaves */
4337:     PetscInt sMin = PETSC_INT_MAX;
4338:     PetscInt sMax = PETSC_INT_MIN;
4339:     PetscInt coneSize, supportSize;

4341:     for (PetscInt p = pStart; p < pEnd; ++p) {
4342:       PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
4343:       PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
4344:       if (!coneSize && supportSize) {
4345:         sMin = PetscMin(p, sMin);
4346:         sMax = PetscMax(p, sMax);
4347:         ++numRoots;
4348:       } else if (!supportSize && coneSize) {
4349:         ++numLeaves;
4350:       } else if (!supportSize && !coneSize) {
4351:         /* Isolated points */
4352:         sMin = PetscMin(p, sMin);
4353:         sMax = PetscMax(p, sMax);
4354:       }
4355:     }
4356:     PetscCall(DMPlexCreateDepthStratum(dm, label, 0, sMin, sMax + 1));
4357:   }

4359:   if (numRoots + numLeaves == (pEnd - pStart)) {
4360:     PetscInt sMin = PETSC_INT_MAX;
4361:     PetscInt sMax = PETSC_INT_MIN;
4362:     PetscInt coneSize, supportSize;

4364:     for (PetscInt p = pStart; p < pEnd; ++p) {
4365:       PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
4366:       PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
4367:       if (!supportSize && coneSize) {
4368:         sMin = PetscMin(p, sMin);
4369:         sMax = PetscMax(p, sMax);
4370:       }
4371:     }
4372:     PetscCall(DMPlexCreateDepthStratum(dm, label, 1, sMin, sMax + 1));
4373:   } else {
4374:     PetscInt level = 0;
4375:     PetscInt qStart, qEnd;

4377:     PetscCall(DMLabelGetStratumBounds(label, level, &qStart, &qEnd));
4378:     while (qEnd > qStart) {
4379:       PetscInt sMin = PETSC_INT_MAX;
4380:       PetscInt sMax = PETSC_INT_MIN;

4382:       for (PetscInt q = qStart; q < qEnd; ++q) {
4383:         const PetscInt *support;
4384:         PetscInt        supportSize;

4386:         PetscCall(DMPlexGetSupportSize(dm, q, &supportSize));
4387:         PetscCall(DMPlexGetSupport(dm, q, &support));
4388:         for (PetscInt s = 0; s < supportSize; ++s) {
4389:           sMin = PetscMin(support[s], sMin);
4390:           sMax = PetscMax(support[s], sMax);
4391:         }
4392:       }
4393:       PetscCall(DMLabelGetNumValues(label, &level));
4394:       PetscCall(DMPlexCreateDepthStratum(dm, label, level, sMin, sMax + 1));
4395:       PetscCall(DMLabelGetStratumBounds(label, level, &qStart, &qEnd));
4396:     }
4397:   }
4398:   PetscFunctionReturn(PETSC_SUCCESS);
4399: }

4401: /*@
4402:   DMPlexStratify - Computes the strata for all points in the `DMPLEX`

4404:   Collective

4406:   Input Parameter:
4407: . dm - The `DMPLEX`

4409:   Level: beginner

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

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

4423:   The depth of a point is calculated by executing a breadth-first search (BFS) on the DAG. This could produce surprising results
4424:   if run on a partially interpolated mesh, meaning one that had some edges and faces, but not others. For example, suppose that
4425:   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
4426:   to interpolate only that one (e0), so that
4427: .vb
4428:   cone(c0) = {e0, v2}
4429:   cone(e0) = {v0, v1}
4430: .ve
4431:   If `DMPlexStratify()` is run on this mesh, it will give depths
4432: .vb
4433:    depth 0 = {v0, v1, v2}
4434:    depth 1 = {e0, c0}
4435: .ve
4436:   where the triangle has been given depth 1, instead of 2, because it is reachable from vertex v2.

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

4440: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSymmetrize()`, `DMPlexComputeCellTypes()`
4441: @*/
4442: PetscErrorCode DMPlexStratify(DM dm)
4443: {
4444:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4445:   DMLabel   label;
4446:   PetscBool flg = PETSC_FALSE;

4448:   PetscFunctionBegin;
4450:   PetscCall(PetscLogEventBegin(DMPLEX_Stratify, dm, 0, 0, 0));

4452:   // Create depth label
4453:   PetscCall(DMRemoveLabel(dm, "depth", NULL));
4454:   PetscCall(DMCreateLabel(dm, "depth"));
4455:   PetscCall(DMPlexGetDepthLabel(dm, &label));

4457:   PetscCall(PetscOptionsGetBool(NULL, dm->hdr.prefix, "-dm_plex_stratify_celltype", &flg, NULL));
4458:   if (flg) PetscCall(DMPlexStratify_CellType_Private(dm, label));
4459:   else PetscCall(DMPlexStratify_Topological_Private(dm, label));

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

4464:     PetscCall(DMLabelGetNumValues(label, &numValues));
4465:     PetscCallMPI(MPIU_Allreduce(&numValues, &maxValues, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
4466:     for (v = numValues; v < maxValues; v++) PetscCall(DMLabelAddStratum(label, v));
4467:   }
4468:   PetscCall(PetscObjectStateGet((PetscObject)label, &mesh->depthState));
4469:   PetscCall(PetscLogEventEnd(DMPLEX_Stratify, dm, 0, 0, 0));
4470:   PetscFunctionReturn(PETSC_SUCCESS);
4471: }

4473: PetscErrorCode DMPlexComputeCellType_Internal(DM dm, PetscInt p, PetscInt pdepth, DMPolytopeType *pt)
4474: {
4475:   DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
4476:   PetscInt       dim, depth, pheight, coneSize;

4478:   PetscFunctionBeginHot;
4479:   PetscCall(DMGetDimension(dm, &dim));
4480:   PetscCall(DMPlexGetDepth(dm, &depth));
4481:   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
4482:   pheight = depth - pdepth;
4483:   if (depth <= 1) {
4484:     switch (pdepth) {
4485:     case 0:
4486:       ct = DM_POLYTOPE_POINT;
4487:       break;
4488:     case 1:
4489:       switch (coneSize) {
4490:       case 2:
4491:         ct = DM_POLYTOPE_SEGMENT;
4492:         break;
4493:       case 3:
4494:         ct = DM_POLYTOPE_TRIANGLE;
4495:         break;
4496:       case 4:
4497:         switch (dim) {
4498:         case 2:
4499:           ct = DM_POLYTOPE_QUADRILATERAL;
4500:           break;
4501:         case 3:
4502:           ct = DM_POLYTOPE_TETRAHEDRON;
4503:           break;
4504:         default:
4505:           break;
4506:         }
4507:         break;
4508:       case 5:
4509:         ct = DM_POLYTOPE_PYRAMID;
4510:         break;
4511:       case 6:
4512:         ct = DM_POLYTOPE_TRI_PRISM_TENSOR;
4513:         break;
4514:       case 8:
4515:         ct = DM_POLYTOPE_HEXAHEDRON;
4516:         break;
4517:       default:
4518:         break;
4519:       }
4520:     }
4521:   } else {
4522:     if (pdepth == 0) {
4523:       ct = DM_POLYTOPE_POINT;
4524:     } else if (pheight == 0) {
4525:       switch (dim) {
4526:       case 1:
4527:         switch (coneSize) {
4528:         case 2:
4529:           ct = DM_POLYTOPE_SEGMENT;
4530:           break;
4531:         default:
4532:           break;
4533:         }
4534:         break;
4535:       case 2:
4536:         switch (coneSize) {
4537:         case 3:
4538:           ct = DM_POLYTOPE_TRIANGLE;
4539:           break;
4540:         case 4:
4541:           ct = DM_POLYTOPE_QUADRILATERAL;
4542:           break;
4543:         default:
4544:           break;
4545:         }
4546:         break;
4547:       case 3:
4548:         switch (coneSize) {
4549:         case 4:
4550:           ct = DM_POLYTOPE_TETRAHEDRON;
4551:           break;
4552:         case 5: {
4553:           const PetscInt *cone;
4554:           PetscInt        faceConeSize;

4556:           PetscCall(DMPlexGetCone(dm, p, &cone));
4557:           PetscCall(DMPlexGetConeSize(dm, cone[0], &faceConeSize));
4558:           switch (faceConeSize) {
4559:           case 3:
4560:             ct = DM_POLYTOPE_TRI_PRISM_TENSOR;
4561:             break;
4562:           case 4:
4563:             ct = DM_POLYTOPE_PYRAMID;
4564:             break;
4565:           }
4566:         } break;
4567:         case 6:
4568:           ct = DM_POLYTOPE_HEXAHEDRON;
4569:           break;
4570:         default:
4571:           break;
4572:         }
4573:         break;
4574:       default:
4575:         break;
4576:       }
4577:     } else if (pheight > 0) {
4578:       switch (coneSize) {
4579:       case 2:
4580:         ct = DM_POLYTOPE_SEGMENT;
4581:         break;
4582:       case 3:
4583:         ct = DM_POLYTOPE_TRIANGLE;
4584:         break;
4585:       case 4:
4586:         ct = DM_POLYTOPE_QUADRILATERAL;
4587:         break;
4588:       default:
4589:         break;
4590:       }
4591:     }
4592:   }
4593:   *pt = ct;
4594:   PetscFunctionReturn(PETSC_SUCCESS);
4595: }

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

4600:   Collective

4602:   Input Parameter:
4603: . dm - The `DMPLEX`

4605:   Level: developer

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

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

4614: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreate()`, `DMPlexSymmetrize()`, `DMPlexStratify()`, `DMGetLabel()`, `DMCreateLabel()`
4615: @*/
4616: PetscErrorCode DMPlexComputeCellTypes(DM dm)
4617: {
4618:   DM_Plex *mesh;
4619:   DMLabel  ctLabel;
4620:   PetscInt pStart, pEnd, p;

4622:   PetscFunctionBegin;
4624:   mesh = (DM_Plex *)dm->data;
4625:   PetscCall(DMCreateLabel(dm, "celltype"));
4626:   PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel));
4627:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4628:   PetscCall(PetscFree(mesh->cellTypes));
4629:   PetscCall(PetscMalloc1(pEnd - pStart, &mesh->cellTypes));
4630:   for (p = pStart; p < pEnd; ++p) {
4631:     DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
4632:     PetscInt       pdepth;

4634:     PetscCall(DMPlexGetPointDepth(dm, p, &pdepth));
4635:     PetscCall(DMPlexComputeCellType_Internal(dm, p, pdepth, &ct));
4636:     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]);
4637:     PetscCall(DMLabelSetValue(ctLabel, p, ct));
4638:     mesh->cellTypes[p - pStart].value_as_uint8 = (uint8_t)ct;
4639:   }
4640:   PetscCall(PetscObjectStateGet((PetscObject)ctLabel, &mesh->celltypeState));
4641:   PetscCall(PetscObjectViewFromOptions((PetscObject)ctLabel, NULL, "-dm_plex_celltypes_view"));
4642:   PetscFunctionReturn(PETSC_SUCCESS);
4643: }

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

4648:   Not Collective

4650:   Input Parameters:
4651: + dm        - The `DMPLEX` object
4652: . numPoints - The number of input points for the join
4653: - points    - The input points

4655:   Output Parameters:
4656: + numCoveredPoints - The number of points in the join
4657: - coveredPoints    - The points in the join

4659:   Level: intermediate

4661:   Note:
4662:   Currently, this is restricted to a single level join

4664:   Fortran Notes:
4665:   `converedPoints` must be declared with
4666: .vb
4667:   PetscInt, pointer :: coveredPints(:)
4668: .ve

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

4672: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreJoin()`, `DMPlexGetMeet()`
4673: @*/
4674: PetscErrorCode DMPlexGetJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4675: {
4676:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4677:   PetscInt *join[2];
4678:   PetscInt  joinSize, i = 0;
4679:   PetscInt  dof, off, p, c, m;
4680:   PetscInt  maxSupportSize;

4682:   PetscFunctionBegin;
4684:   PetscAssertPointer(points, 3);
4685:   PetscAssertPointer(numCoveredPoints, 4);
4686:   PetscAssertPointer(coveredPoints, 5);
4687:   PetscCall(PetscSectionGetMaxDof(mesh->supportSection, &maxSupportSize));
4688:   PetscCall(DMGetWorkArray(dm, maxSupportSize, MPIU_INT, &join[0]));
4689:   PetscCall(DMGetWorkArray(dm, maxSupportSize, MPIU_INT, &join[1]));
4690:   /* Copy in support of first point */
4691:   PetscCall(PetscSectionGetDof(mesh->supportSection, points[0], &dof));
4692:   PetscCall(PetscSectionGetOffset(mesh->supportSection, points[0], &off));
4693:   for (joinSize = 0; joinSize < dof; ++joinSize) join[i][joinSize] = mesh->supports[off + joinSize];
4694:   /* Check each successive support */
4695:   for (p = 1; p < numPoints; ++p) {
4696:     PetscInt newJoinSize = 0;

4698:     PetscCall(PetscSectionGetDof(mesh->supportSection, points[p], &dof));
4699:     PetscCall(PetscSectionGetOffset(mesh->supportSection, points[p], &off));
4700:     for (c = 0; c < dof; ++c) {
4701:       const PetscInt point = mesh->supports[off + c];

4703:       for (m = 0; m < joinSize; ++m) {
4704:         if (point == join[i][m]) {
4705:           join[1 - i][newJoinSize++] = point;
4706:           break;
4707:         }
4708:       }
4709:     }
4710:     joinSize = newJoinSize;
4711:     i        = 1 - i;
4712:   }
4713:   *numCoveredPoints = joinSize;
4714:   *coveredPoints    = join[i];
4715:   PetscCall(DMRestoreWorkArray(dm, maxSupportSize, MPIU_INT, &join[1 - i]));
4716:   PetscFunctionReturn(PETSC_SUCCESS);
4717: }

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

4722:   Not Collective

4724:   Input Parameters:
4725: + dm        - The `DMPLEX` object
4726: . numPoints - The number of input points for the join
4727: - points    - The input points

4729:   Output Parameters:
4730: + numCoveredPoints - The number of points in the join
4731: - coveredPoints    - The points in the join

4733:   Level: intermediate

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

4738: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetJoin()`, `DMPlexGetFullJoin()`, `DMPlexGetMeet()`
4739: @*/
4740: PetscErrorCode DMPlexRestoreJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4741: {
4742:   PetscFunctionBegin;
4744:   if (points) PetscAssertPointer(points, 3);
4745:   if (numCoveredPoints) PetscAssertPointer(numCoveredPoints, 4);
4746:   PetscAssertPointer(coveredPoints, 5);
4747:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, (void *)coveredPoints));
4748:   if (numCoveredPoints) *numCoveredPoints = 0;
4749:   PetscFunctionReturn(PETSC_SUCCESS);
4750: }

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

4755:   Not Collective

4757:   Input Parameters:
4758: + dm        - The `DMPLEX` object
4759: . numPoints - The number of input points for the join
4760: - points    - The input points, its length is `numPoints`

4762:   Output Parameters:
4763: + numCoveredPoints - The number of points in the join
4764: - coveredPoints    - The points in the join, its length is `numCoveredPoints`

4766:   Level: intermediate

4768:   Fortran Notes:
4769:   `points` and `converedPoints` must be declared with
4770: .vb
4771:   PetscInt, pointer :: points(:)
4772:   PetscInt, pointer :: coveredPints(:)
4773: .ve

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

4777: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetJoin()`, `DMPlexRestoreJoin()`, `DMPlexGetMeet()`
4778: @*/
4779: PetscErrorCode DMPlexGetFullJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4780: {
4781:   PetscInt *offsets, **closures;
4782:   PetscInt *join[2];
4783:   PetscInt  depth = 0, maxSize, joinSize = 0, i = 0;
4784:   PetscInt  p, d, c, m, ms;

4786:   PetscFunctionBegin;
4788:   PetscAssertPointer(points, 3);
4789:   PetscAssertPointer(numCoveredPoints, 4);
4790:   PetscAssertPointer(coveredPoints, 5);

4792:   PetscCall(DMPlexGetDepth(dm, &depth));
4793:   PetscCall(PetscCalloc1(numPoints, &closures));
4794:   PetscCall(DMGetWorkArray(dm, numPoints * (depth + 2), MPIU_INT, &offsets));
4795:   PetscCall(DMPlexGetMaxSizes(dm, NULL, &ms));
4796:   maxSize = (ms > 1) ? ((PetscPowInt(ms, depth + 1) - 1) / (ms - 1)) : depth + 1;
4797:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &join[0]));
4798:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &join[1]));

4800:   for (p = 0; p < numPoints; ++p) {
4801:     PetscInt closureSize;

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

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

4809:       PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
4810:       for (i = offsets[p * (depth + 2) + d]; i < closureSize; ++i) {
4811:         if ((pStart > closures[p][i * 2]) || (pEnd <= closures[p][i * 2])) {
4812:           offsets[p * (depth + 2) + d + 1] = i;
4813:           break;
4814:         }
4815:       }
4816:       if (i == closureSize) offsets[p * (depth + 2) + d + 1] = i;
4817:     }
4818:     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);
4819:   }
4820:   for (d = 0; d < depth + 1; ++d) {
4821:     PetscInt dof;

4823:     /* Copy in support of first point */
4824:     dof = offsets[d + 1] - offsets[d];
4825:     for (joinSize = 0; joinSize < dof; ++joinSize) join[i][joinSize] = closures[0][(offsets[d] + joinSize) * 2];
4826:     /* Check each successive cone */
4827:     for (p = 1; p < numPoints && joinSize; ++p) {
4828:       PetscInt newJoinSize = 0;

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

4834:         for (m = 0; m < joinSize; ++m) {
4835:           if (point == join[i][m]) {
4836:             join[1 - i][newJoinSize++] = point;
4837:             break;
4838:           }
4839:         }
4840:       }
4841:       joinSize = newJoinSize;
4842:       i        = 1 - i;
4843:     }
4844:     if (joinSize) break;
4845:   }
4846:   *numCoveredPoints = joinSize;
4847:   *coveredPoints    = join[i];
4848:   for (p = 0; p < numPoints; ++p) PetscCall(DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, NULL, &closures[p]));
4849:   PetscCall(PetscFree(closures));
4850:   PetscCall(DMRestoreWorkArray(dm, numPoints * (depth + 2), MPIU_INT, &offsets));
4851:   PetscCall(DMRestoreWorkArray(dm, ms, MPIU_INT, &join[1 - i]));
4852:   PetscFunctionReturn(PETSC_SUCCESS);
4853: }

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

4858:   Not Collective

4860:   Input Parameters:
4861: + dm        - The `DMPLEX` object
4862: . numPoints - The number of input points for the meet
4863: - points    - The input points, of length `numPoints`

4865:   Output Parameters:
4866: + numCoveringPoints - The number of points in the meet
4867: - coveringPoints    - The points in the meet, of length `numCoveringPoints`

4869:   Level: intermediate

4871:   Note:
4872:   Currently, this is restricted to a single level meet

4874:   Fortran Notes:
4875:   `coveringPoints` must be declared with
4876: .vb
4877:   PetscInt, pointer :: coveringPoints(:)
4878: .ve

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

4882: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexRestoreMeet()`, `DMPlexGetJoin()`
4883: @*/
4884: PetscErrorCode DMPlexGetMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveringPoints, const PetscInt *coveringPoints[])
4885: {
4886:   DM_Plex  *mesh = (DM_Plex *)dm->data;
4887:   PetscInt *meet[2];
4888:   PetscInt  meetSize, i = 0;
4889:   PetscInt  dof, off, p, c, m;
4890:   PetscInt  maxConeSize;

4892:   PetscFunctionBegin;
4894:   PetscAssertPointer(points, 3);
4895:   PetscAssertPointer(numCoveringPoints, 4);
4896:   PetscAssertPointer(coveringPoints, 5);
4897:   PetscCall(PetscSectionGetMaxDof(mesh->coneSection, &maxConeSize));
4898:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &meet[0]));
4899:   PetscCall(DMGetWorkArray(dm, maxConeSize, MPIU_INT, &meet[1]));
4900:   /* Copy in cone of first point */
4901:   PetscCall(PetscSectionGetDof(mesh->coneSection, points[0], &dof));
4902:   PetscCall(PetscSectionGetOffset(mesh->coneSection, points[0], &off));
4903:   for (meetSize = 0; meetSize < dof; ++meetSize) meet[i][meetSize] = mesh->cones[off + meetSize];
4904:   /* Check each successive cone */
4905:   for (p = 1; p < numPoints; ++p) {
4906:     PetscInt newMeetSize = 0;

4908:     PetscCall(PetscSectionGetDof(mesh->coneSection, points[p], &dof));
4909:     PetscCall(PetscSectionGetOffset(mesh->coneSection, points[p], &off));
4910:     for (c = 0; c < dof; ++c) {
4911:       const PetscInt point = mesh->cones[off + c];

4913:       for (m = 0; m < meetSize; ++m) {
4914:         if (point == meet[i][m]) {
4915:           meet[1 - i][newMeetSize++] = point;
4916:           break;
4917:         }
4918:       }
4919:     }
4920:     meetSize = newMeetSize;
4921:     i        = 1 - i;
4922:   }
4923:   *numCoveringPoints = meetSize;
4924:   *coveringPoints    = meet[i];
4925:   PetscCall(DMRestoreWorkArray(dm, maxConeSize, MPIU_INT, &meet[1 - i]));
4926:   PetscFunctionReturn(PETSC_SUCCESS);
4927: }

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

4932:   Not Collective

4934:   Input Parameters:
4935: + dm        - The `DMPLEX` object
4936: . numPoints - The number of input points for the meet
4937: - points    - The input points

4939:   Output Parameters:
4940: + numCoveredPoints - The number of points in the meet
4941: - coveredPoints    - The points in the meet

4943:   Level: intermediate

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

4948: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetMeet()`, `DMPlexGetFullMeet()`, `DMPlexGetJoin()`
4949: @*/
4950: PetscErrorCode DMPlexRestoreMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4951: {
4952:   PetscFunctionBegin;
4954:   if (points) PetscAssertPointer(points, 3);
4955:   if (numCoveredPoints) PetscAssertPointer(numCoveredPoints, 4);
4956:   PetscAssertPointer(coveredPoints, 5);
4957:   PetscCall(DMRestoreWorkArray(dm, 0, MPIU_INT, (void *)coveredPoints));
4958:   if (numCoveredPoints) *numCoveredPoints = 0;
4959:   PetscFunctionReturn(PETSC_SUCCESS);
4960: }

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

4965:   Not Collective

4967:   Input Parameters:
4968: + dm        - The `DMPLEX` object
4969: . numPoints - The number of input points for the meet
4970: - points    - The input points, of length  `numPoints`

4972:   Output Parameters:
4973: + numCoveredPoints - The number of points in the meet
4974: - coveredPoints    - The points in the meet, of length  `numCoveredPoints`

4976:   Level: intermediate

4978:   Fortran Notes:
4979:   `points` and `coveredPoints` must be declared with
4980: .vb
4981:   PetscInt, pointer :: points(:)
4982:   PetscInt, pointer :: coveredPoints(:)
4983: .ve

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

4987: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetMeet()`, `DMPlexRestoreMeet()`, `DMPlexGetJoin()`
4988: @*/
4989: PetscErrorCode DMPlexGetFullMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt *coveredPoints[])
4990: {
4991:   PetscInt *offsets, **closures;
4992:   PetscInt *meet[2];
4993:   PetscInt  height = 0, maxSize, meetSize = 0, i = 0;
4994:   PetscInt  p, h, c, m, mc;

4996:   PetscFunctionBegin;
4998:   PetscAssertPointer(points, 3);
4999:   PetscAssertPointer(numCoveredPoints, 4);
5000:   PetscAssertPointer(coveredPoints, 5);

5002:   PetscCall(DMPlexGetDepth(dm, &height));
5003:   PetscCall(PetscMalloc1(numPoints, &closures));
5004:   PetscCall(DMGetWorkArray(dm, numPoints * (height + 2), MPIU_INT, &offsets));
5005:   PetscCall(DMPlexGetMaxSizes(dm, &mc, NULL));
5006:   maxSize = (mc > 1) ? ((PetscPowInt(mc, height + 1) - 1) / (mc - 1)) : height + 1;
5007:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &meet[0]));
5008:   PetscCall(DMGetWorkArray(dm, maxSize, MPIU_INT, &meet[1]));

5010:   for (p = 0; p < numPoints; ++p) {
5011:     PetscInt closureSize;

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

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

5019:       PetscCall(DMPlexGetHeightStratum(dm, h, &pStart, &pEnd));
5020:       for (i = offsets[p * (height + 2) + h]; i < closureSize; ++i) {
5021:         if ((pStart > closures[p][i * 2]) || (pEnd <= closures[p][i * 2])) {
5022:           offsets[p * (height + 2) + h + 1] = i;
5023:           break;
5024:         }
5025:       }
5026:       if (i == closureSize) offsets[p * (height + 2) + h + 1] = i;
5027:     }
5028:     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);
5029:   }
5030:   for (h = 0; h < height + 1; ++h) {
5031:     PetscInt dof;

5033:     /* Copy in cone of first point */
5034:     dof = offsets[h + 1] - offsets[h];
5035:     for (meetSize = 0; meetSize < dof; ++meetSize) meet[i][meetSize] = closures[0][(offsets[h] + meetSize) * 2];
5036:     /* Check each successive cone */
5037:     for (p = 1; p < numPoints && meetSize; ++p) {
5038:       PetscInt newMeetSize = 0;

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

5044:         for (m = 0; m < meetSize; ++m) {
5045:           if (point == meet[i][m]) {
5046:             meet[1 - i][newMeetSize++] = point;
5047:             break;
5048:           }
5049:         }
5050:       }
5051:       meetSize = newMeetSize;
5052:       i        = 1 - i;
5053:     }
5054:     if (meetSize) break;
5055:   }
5056:   *numCoveredPoints = meetSize;
5057:   *coveredPoints    = meet[i];
5058:   for (p = 0; p < numPoints; ++p) PetscCall(DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, NULL, &closures[p]));
5059:   PetscCall(PetscFree(closures));
5060:   PetscCall(DMRestoreWorkArray(dm, numPoints * (height + 2), MPIU_INT, &offsets));
5061:   PetscCall(DMRestoreWorkArray(dm, mc, MPIU_INT, &meet[1 - i]));
5062:   PetscFunctionReturn(PETSC_SUCCESS);
5063: }

5065: /*@
5066:   DMPlexEqual - Determine if two `DM` have the same topology

5068:   Not Collective

5070:   Input Parameters:
5071: + dmA - A `DMPLEX` object
5072: - dmB - A `DMPLEX` object

5074:   Output Parameter:
5075: . equal - `PETSC_TRUE` if the topologies are identical

5077:   Level: intermediate

5079:   Note:
5080:   We are not solving graph isomorphism, so we do not permute.

5082: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCone()`
5083: @*/
5084: PetscErrorCode DMPlexEqual(DM dmA, DM dmB, PetscBool *equal)
5085: {
5086:   PetscInt depth, depthB, pStart, pEnd, pStartB, pEndB, p;

5088:   PetscFunctionBegin;
5091:   PetscAssertPointer(equal, 3);

5093:   *equal = PETSC_FALSE;
5094:   PetscCall(DMPlexGetDepth(dmA, &depth));
5095:   PetscCall(DMPlexGetDepth(dmB, &depthB));
5096:   if (depth != depthB) PetscFunctionReturn(PETSC_SUCCESS);
5097:   PetscCall(DMPlexGetChart(dmA, &pStart, &pEnd));
5098:   PetscCall(DMPlexGetChart(dmB, &pStartB, &pEndB));
5099:   if ((pStart != pStartB) || (pEnd != pEndB)) PetscFunctionReturn(PETSC_SUCCESS);
5100:   for (p = pStart; p < pEnd; ++p) {
5101:     const PetscInt *cone, *coneB, *ornt, *orntB, *support, *supportB;
5102:     PetscInt        coneSize, coneSizeB, c, supportSize, supportSizeB, s;

5104:     PetscCall(DMPlexGetConeSize(dmA, p, &coneSize));
5105:     PetscCall(DMPlexGetCone(dmA, p, &cone));
5106:     PetscCall(DMPlexGetConeOrientation(dmA, p, &ornt));
5107:     PetscCall(DMPlexGetConeSize(dmB, p, &coneSizeB));
5108:     PetscCall(DMPlexGetCone(dmB, p, &coneB));
5109:     PetscCall(DMPlexGetConeOrientation(dmB, p, &orntB));
5110:     if (coneSize != coneSizeB) PetscFunctionReturn(PETSC_SUCCESS);
5111:     for (c = 0; c < coneSize; ++c) {
5112:       if (cone[c] != coneB[c]) PetscFunctionReturn(PETSC_SUCCESS);
5113:       if (ornt[c] != orntB[c]) PetscFunctionReturn(PETSC_SUCCESS);
5114:     }
5115:     PetscCall(DMPlexGetSupportSize(dmA, p, &supportSize));
5116:     PetscCall(DMPlexGetSupport(dmA, p, &support));
5117:     PetscCall(DMPlexGetSupportSize(dmB, p, &supportSizeB));
5118:     PetscCall(DMPlexGetSupport(dmB, p, &supportB));
5119:     if (supportSize != supportSizeB) PetscFunctionReturn(PETSC_SUCCESS);
5120:     for (s = 0; s < supportSize; ++s) {
5121:       if (support[s] != supportB[s]) PetscFunctionReturn(PETSC_SUCCESS);
5122:     }
5123:   }
5124:   *equal = PETSC_TRUE;
5125:   PetscFunctionReturn(PETSC_SUCCESS);
5126: }

5128: /*@
5129:   DMPlexGetNumFaceVertices - Returns the number of vertices on a face

5131:   Not Collective

5133:   Input Parameters:
5134: + dm         - The `DMPLEX`
5135: . cellDim    - The cell dimension
5136: - numCorners - The number of vertices on a cell

5138:   Output Parameter:
5139: . numFaceVertices - The number of vertices on a face

5141:   Level: developer

5143:   Note:
5144:   Of course this can only work for a restricted set of symmetric shapes

5146: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCone()`
5147: @*/
5148: PetscErrorCode DMPlexGetNumFaceVertices(DM dm, PetscInt cellDim, PetscInt numCorners, PetscInt *numFaceVertices)
5149: {
5150:   MPI_Comm comm;

5152:   PetscFunctionBegin;
5153:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5154:   PetscAssertPointer(numFaceVertices, 4);
5155:   switch (cellDim) {
5156:   case 0:
5157:     *numFaceVertices = 0;
5158:     break;
5159:   case 1:
5160:     *numFaceVertices = 1;
5161:     break;
5162:   case 2:
5163:     switch (numCorners) {
5164:     case 3:                 /* triangle */
5165:       *numFaceVertices = 2; /* Edge has 2 vertices */
5166:       break;
5167:     case 4:                 /* quadrilateral */
5168:       *numFaceVertices = 2; /* Edge has 2 vertices */
5169:       break;
5170:     case 6:                 /* quadratic triangle, tri and quad cohesive Lagrange cells */
5171:       *numFaceVertices = 3; /* Edge has 3 vertices */
5172:       break;
5173:     case 9:                 /* quadratic quadrilateral, quadratic quad cohesive Lagrange cells */
5174:       *numFaceVertices = 3; /* Edge has 3 vertices */
5175:       break;
5176:     default:
5177:       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %" PetscInt_FMT " for dimension %" PetscInt_FMT, numCorners, cellDim);
5178:     }
5179:     break;
5180:   case 3:
5181:     switch (numCorners) {
5182:     case 4:                 /* tetradehdron */
5183:       *numFaceVertices = 3; /* Face has 3 vertices */
5184:       break;
5185:     case 6:                 /* tet cohesive cells */
5186:       *numFaceVertices = 4; /* Face has 4 vertices */
5187:       break;
5188:     case 8:                 /* hexahedron */
5189:       *numFaceVertices = 4; /* Face has 4 vertices */
5190:       break;
5191:     case 9:                 /* tet cohesive Lagrange cells */
5192:       *numFaceVertices = 6; /* Face has 6 vertices */
5193:       break;
5194:     case 10:                /* quadratic tetrahedron */
5195:       *numFaceVertices = 6; /* Face has 6 vertices */
5196:       break;
5197:     case 12:                /* hex cohesive Lagrange cells */
5198:       *numFaceVertices = 6; /* Face has 6 vertices */
5199:       break;
5200:     case 18:                /* quadratic tet cohesive Lagrange cells */
5201:       *numFaceVertices = 6; /* Face has 6 vertices */
5202:       break;
5203:     case 27:                /* quadratic hexahedron, quadratic hex cohesive Lagrange cells */
5204:       *numFaceVertices = 9; /* Face has 9 vertices */
5205:       break;
5206:     default:
5207:       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %" PetscInt_FMT " for dimension %" PetscInt_FMT, numCorners, cellDim);
5208:     }
5209:     break;
5210:   default:
5211:     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid cell dimension %" PetscInt_FMT, cellDim);
5212:   }
5213:   PetscFunctionReturn(PETSC_SUCCESS);
5214: }

5216: /*@
5217:   DMPlexGetDepthLabel - Get the `DMLabel` recording the depth of each point

5219:   Not Collective

5221:   Input Parameter:
5222: . dm - The `DMPLEX` object

5224:   Output Parameter:
5225: . depthLabel - The `DMLabel` recording point depth

5227:   Level: developer

5229: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetDepth()`, `DMPlexGetHeightStratum()`, `DMPlexGetDepthStratum()`, `DMPlexGetPointDepth()`,
5230: @*/
5231: PetscErrorCode DMPlexGetDepthLabel(DM dm, DMLabel *depthLabel)
5232: {
5233:   PetscFunctionBegin;
5235:   PetscAssertPointer(depthLabel, 2);
5236:   *depthLabel = dm->depthLabel;
5237:   PetscFunctionReturn(PETSC_SUCCESS);
5238: }

5240: /*@
5241:   DMPlexGetDepth - Get the depth of the DAG representing this mesh

5243:   Not Collective

5245:   Input Parameter:
5246: . dm - The `DMPLEX` object

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

5251:   Level: developer

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

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

5258:   An empty mesh gives -1.

5260: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetDepthLabel()`, `DMPlexGetDepthStratum()`, `DMPlexGetPointDepth()`, `DMPlexSymmetrize()`
5261: @*/
5262: PetscErrorCode DMPlexGetDepth(DM dm, PetscInt *depth)
5263: {
5264:   DM_Plex *mesh = (DM_Plex *)dm->data;
5265:   DMLabel  label;
5266:   PetscInt d = -1;

5268:   PetscFunctionBegin;
5270:   PetscAssertPointer(depth, 2);
5271:   if (mesh->tr) {
5272:     PetscCall(DMPlexTransformGetDepth(mesh->tr, depth));
5273:   } else {
5274:     PetscCall(DMPlexGetDepthLabel(dm, &label));
5275:     // Allow missing depths
5276:     if (label) PetscCall(DMLabelGetValueBounds(label, NULL, &d));
5277:     *depth = d;
5278:   }
5279:   PetscFunctionReturn(PETSC_SUCCESS);
5280: }

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

5285:   Not Collective

5287:   Input Parameters:
5288: + dm    - The `DMPLEX` object
5289: - depth - The requested depth

5291:   Output Parameters:
5292: + start - The first point at this `depth`
5293: - end   - One beyond the last point at this `depth`

5295:   Level: developer

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

5302: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetHeightStratum()`, `DMPlexGetCellTypeStratum()`, `DMPlexGetDepth()`, `DMPlexGetDepthLabel()`, `DMPlexGetPointDepth()`, `DMPlexSymmetrize()`, `DMPlexInterpolate()`
5303: @*/
5304: PetscErrorCode DMPlexGetDepthStratum(DM dm, PetscInt depth, PetscInt *start, PetscInt *end)
5305: {
5306:   DM_Plex *mesh = (DM_Plex *)dm->data;
5307:   DMLabel  label;
5308:   PetscInt pStart, pEnd;

5310:   PetscFunctionBegin;
5312:   if (start) {
5313:     PetscAssertPointer(start, 3);
5314:     *start = 0;
5315:   }
5316:   if (end) {
5317:     PetscAssertPointer(end, 4);
5318:     *end = 0;
5319:   }
5320:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5321:   if (pStart == pEnd) PetscFunctionReturn(PETSC_SUCCESS);
5322:   if (depth < 0) {
5323:     if (start) *start = pStart;
5324:     if (end) *end = pEnd;
5325:     PetscFunctionReturn(PETSC_SUCCESS);
5326:   }
5327:   if (mesh->tr) {
5328:     PetscCall(DMPlexTransformGetDepthStratum(mesh->tr, depth, start, end));
5329:   } else {
5330:     PetscCall(DMPlexGetDepthLabel(dm, &label));
5331:     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "No label named depth was found");
5332:     PetscCall(DMLabelGetStratumBounds(label, depth, start, end));
5333:   }
5334:   PetscFunctionReturn(PETSC_SUCCESS);
5335: }

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

5340:   Not Collective

5342:   Input Parameters:
5343: + dm     - The `DMPLEX` object
5344: - height - The requested height

5346:   Output Parameters:
5347: + start - The first point at this `height`
5348: - end   - One beyond the last point at this `height`

5350:   Level: developer

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

5357: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetDepthStratum()`, `DMPlexGetCellTypeStratum()`, `DMPlexGetDepth()`, `DMPlexGetPointHeight()`
5358: @*/
5359: PetscErrorCode DMPlexGetHeightStratum(DM dm, PetscInt height, PetscInt *start, PetscInt *end)
5360: {
5361:   DMLabel  label;
5362:   PetscInt depth, pStart, pEnd;

5364:   PetscFunctionBegin;
5366:   if (start) {
5367:     PetscAssertPointer(start, 3);
5368:     *start = 0;
5369:   }
5370:   if (end) {
5371:     PetscAssertPointer(end, 4);
5372:     *end = 0;
5373:   }
5374:   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5375:   if (pStart == pEnd) PetscFunctionReturn(PETSC_SUCCESS);
5376:   if (height < 0) {
5377:     if (start) *start = pStart;
5378:     if (end) *end = pEnd;
5379:     PetscFunctionReturn(PETSC_SUCCESS);
5380:   }
5381:   PetscCall(DMPlexGetDepthLabel(dm, &label));
5382:   if (label) PetscCall(DMLabelGetNumValues(label, &depth));
5383:   else PetscCall(DMGetDimension(dm, &depth));
5384:   PetscCheck(depth >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Depth not yet computed");
5385:   PetscCall(DMPlexGetDepthStratum(dm, depth - 1 - height, start, end));
5386:   PetscFunctionReturn(PETSC_SUCCESS);
5387: }

5389: /*@
5390:   DMPlexGetPointDepth - Get the `depth` of a given point

5392:   Not Collective

5394:   Input Parameters:
5395: + dm    - The `DMPLEX` object
5396: - point - The point

5398:   Output Parameter:
5399: . depth - The depth of the `point`

5401:   Level: intermediate

5403: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellType()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`, `DMPlexGetPointHeight()`
5404: @*/
5405: PetscErrorCode DMPlexGetPointDepth(DM dm, PetscInt point, PetscInt *depth)
5406: {
5407:   PetscFunctionBegin;
5409:   PetscAssertPointer(depth, 3);
5410:   PetscCall(DMLabelGetValue(dm->depthLabel, point, depth));
5411:   PetscFunctionReturn(PETSC_SUCCESS);
5412: }

5414: /*@
5415:   DMPlexGetPointHeight - Get the `height` of a given point

5417:   Not Collective

5419:   Input Parameters:
5420: + dm    - The `DMPLEX` object
5421: - point - The point

5423:   Output Parameter:
5424: . height - The height of the `point`

5426:   Level: intermediate

5428: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellType()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`, `DMPlexGetPointDepth()`
5429: @*/
5430: PetscErrorCode DMPlexGetPointHeight(DM dm, PetscInt point, PetscInt *height)
5431: {
5432:   PetscInt n, pDepth;

5434:   PetscFunctionBegin;
5436:   PetscAssertPointer(height, 3);
5437:   PetscCall(DMLabelGetNumValues(dm->depthLabel, &n));
5438:   PetscCall(DMLabelGetValue(dm->depthLabel, point, &pDepth));
5439:   *height = n - 1 - pDepth; /* DAG depth is n-1 */
5440:   PetscFunctionReturn(PETSC_SUCCESS);
5441: }

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

5446:   Not Collective

5448:   Input Parameter:
5449: . dm - The `DMPLEX` object

5451:   Output Parameter:
5452: . celltypeLabel - The `DMLabel` recording cell polytope type

5454:   Level: developer

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

5460: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellType()`, `DMPlexGetDepthLabel()`, `DMCreateLabel()`
5461: @*/
5462: PetscErrorCode DMPlexGetCellTypeLabel(DM dm, DMLabel *celltypeLabel)
5463: {
5464:   PetscFunctionBegin;
5466:   PetscAssertPointer(celltypeLabel, 2);
5467:   if (!dm->celltypeLabel) PetscCall(DMPlexComputeCellTypes(dm));
5468:   *celltypeLabel = dm->celltypeLabel;
5469:   PetscFunctionReturn(PETSC_SUCCESS);
5470: }

5472: /*@
5473:   DMPlexGetCellType - Get the polytope type of a given cell

5475:   Not Collective

5477:   Input Parameters:
5478: + dm   - The `DMPLEX` object
5479: - cell - The cell

5481:   Output Parameter:
5482: . celltype - The polytope type of the cell

5484:   Level: intermediate

5486: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPolytopeType`, `DMPlexGetCellTypeLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`
5487: @*/
5488: PetscErrorCode DMPlexGetCellType(DM dm, PetscInt cell, DMPolytopeType *celltype)
5489: {
5490:   DM_Plex *mesh = (DM_Plex *)dm->data;
5491:   DMLabel  label;
5492:   PetscInt ct;

5494:   PetscFunctionBegin;
5496:   PetscAssertPointer(celltype, 3);
5497:   if (mesh->tr) {
5498:     PetscCall(DMPlexTransformGetCellType(mesh->tr, cell, celltype));
5499:   } else {
5500:     PetscInt pStart, pEnd;

5502:     PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, NULL));
5503:     if (!mesh->cellTypes) { /* XXX remove? optimize? */
5504:       PetscCall(PetscSectionGetChart(mesh->coneSection, NULL, &pEnd));
5505:       PetscCall(PetscMalloc1(pEnd - pStart, &mesh->cellTypes));
5506:       PetscCall(DMPlexGetCellTypeLabel(dm, &label));
5507:       for (PetscInt p = pStart; p < pEnd; p++) {
5508:         PetscCall(DMLabelGetValue(label, p, &ct));
5509:         mesh->cellTypes[p - pStart].value_as_uint8 = (uint8_t)ct;
5510:       }
5511:     }
5512:     *celltype = (DMPolytopeType)mesh->cellTypes[cell - pStart].value_as_uint8;
5513:     if (PetscDefined(USE_DEBUG)) {
5514:       PetscCall(DMPlexGetCellTypeLabel(dm, &label));
5515:       PetscCall(DMLabelGetValue(label, cell, &ct));
5516:       PetscCheck(ct >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %" PetscInt_FMT " has not been assigned a cell type", cell);
5517:       PetscCheck(ct == (PetscInt)*celltype, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid cellType for %" PetscInt_FMT ": %d != %" PetscInt_FMT, cell, (int)*celltype, ct);
5518:     }
5519:   }
5520:   PetscFunctionReturn(PETSC_SUCCESS);
5521: }

5523: /*@
5524:   DMPlexSetCellType - Set the polytope type of a given cell

5526:   Not Collective

5528:   Input Parameters:
5529: + dm       - The `DMPLEX` object
5530: . cell     - The cell
5531: - celltype - The polytope type of the cell

5533:   Level: advanced

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

5541: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellTypeLabel()`, `DMPlexGetDepthLabel()`, `DMPlexGetDepth()`, `DMPlexComputeCellTypes()`, `DMCreateLabel()`
5542: @*/
5543: PetscErrorCode DMPlexSetCellType(DM dm, PetscInt cell, DMPolytopeType celltype)
5544: {
5545:   DM_Plex *mesh = (DM_Plex *)dm->data;
5546:   DMLabel  label;
5547:   PetscInt pStart, pEnd;

5549:   PetscFunctionBegin;
5551:   PetscCall(PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd));
5552:   PetscCall(DMPlexGetCellTypeLabel(dm, &label));
5553:   PetscCall(DMLabelSetValue(label, cell, celltype));
5554:   if (!mesh->cellTypes) PetscCall(PetscMalloc1(pEnd - pStart, &mesh->cellTypes));
5555:   mesh->cellTypes[cell - pStart].value_as_uint8 = (uint8_t)celltype;
5556:   PetscFunctionReturn(PETSC_SUCCESS);
5557: }

5559: PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm)
5560: {
5561:   PetscSection section;
5562:   PetscInt     maxHeight;
5563:   const char  *prefix;

5565:   PetscFunctionBegin;
5566:   PetscCall(DMClone(dm, cdm));
5567:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
5568:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix));
5569:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_"));
5570:   PetscCall(DMPlexGetMaxProjectionHeight(dm, &maxHeight));
5571:   PetscCall(DMPlexSetMaxProjectionHeight(*cdm, maxHeight));
5572:   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
5573:   PetscCall(DMSetLocalSection(*cdm, section));
5574:   PetscCall(PetscSectionDestroy(&section));

5576:   PetscCall(DMSetNumFields(*cdm, 1));
5577:   PetscCall(DMCreateDS(*cdm));
5578:   (*cdm)->cloneOpts = PETSC_TRUE;
5579:   if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(*cdm));
5580:   PetscFunctionReturn(PETSC_SUCCESS);
5581: }

5583: PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field)
5584: {
5585:   Vec coordsLocal, cellCoordsLocal;
5586:   DM  coordsDM, cellCoordsDM;

5588:   PetscFunctionBegin;
5589:   *field = NULL;
5590:   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
5591:   PetscCall(DMGetCoordinateDM(dm, &coordsDM));
5592:   PetscCall(DMGetCellCoordinatesLocal(dm, &cellCoordsLocal));
5593:   PetscCall(DMGetCellCoordinateDM(dm, &cellCoordsDM));
5594:   if (coordsLocal && coordsDM) {
5595:     if (cellCoordsLocal && cellCoordsDM) PetscCall(DMFieldCreateDSWithDG(coordsDM, cellCoordsDM, 0, coordsLocal, cellCoordsLocal, field));
5596:     else PetscCall(DMFieldCreateDS(coordsDM, 0, coordsLocal, field));
5597:   }
5598:   PetscFunctionReturn(PETSC_SUCCESS);
5599: }

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

5604:   Not Collective

5606:   Input Parameter:
5607: . dm - The `DMPLEX` object

5609:   Output Parameter:
5610: . section - The `PetscSection` object

5612:   Level: developer

5614: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetSupportSection()`, `DMPlexGetCones()`, `DMPlexGetConeOrientations()`, `PetscSection`
5615: @*/
5616: PetscErrorCode DMPlexGetConeSection(DM dm, PetscSection *section)
5617: {
5618:   DM_Plex *mesh = (DM_Plex *)dm->data;

5620:   PetscFunctionBegin;
5622:   if (section) *section = mesh->coneSection;
5623:   PetscFunctionReturn(PETSC_SUCCESS);
5624: }

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

5629:   Not Collective

5631:   Input Parameter:
5632: . dm - The `DMPLEX` object

5634:   Output Parameter:
5635: . section - The `PetscSection` object

5637:   Level: developer

5639: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSection()`, `PetscSection`
5640: @*/
5641: PetscErrorCode DMPlexGetSupportSection(DM dm, PetscSection *section)
5642: {
5643:   DM_Plex *mesh = (DM_Plex *)dm->data;

5645:   PetscFunctionBegin;
5647:   if (section) *section = mesh->supportSection;
5648:   PetscFunctionReturn(PETSC_SUCCESS);
5649: }

5651: /*@C
5652:   DMPlexGetCones - Return cone data

5654:   Not Collective

5656:   Input Parameter:
5657: . dm - The `DMPLEX` object

5659:   Output Parameter:
5660: . cones - The cone for each point

5662:   Level: developer

5664: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSection()`
5665: @*/
5666: PetscErrorCode DMPlexGetCones(DM dm, PetscInt *cones[])
5667: {
5668:   DM_Plex *mesh = (DM_Plex *)dm->data;

5670:   PetscFunctionBegin;
5672:   if (cones) *cones = mesh->cones;
5673:   PetscFunctionReturn(PETSC_SUCCESS);
5674: }

5676: /*@C
5677:   DMPlexGetConeOrientations - Return cone orientation data

5679:   Not Collective

5681:   Input Parameter:
5682: . dm - The `DMPLEX` object

5684:   Output Parameter:
5685: . coneOrientations - The array of cone orientations for all points

5687:   Level: developer

5689:   Notes:
5690:   The `PetscSection` returned by `DMPlexGetConeSection()` partitions coneOrientations into cone orientations of particular points
5691:   as returned by `DMPlexGetConeOrientation()`.

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

5695: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetConeSection()`, `DMPlexGetConeOrientation()`, `PetscSection`
5696: @*/
5697: PetscErrorCode DMPlexGetConeOrientations(DM dm, PetscInt *coneOrientations[])
5698: {
5699:   DM_Plex *mesh = (DM_Plex *)dm->data;

5701:   PetscFunctionBegin;
5703:   if (coneOrientations) *coneOrientations = mesh->coneOrientations;
5704:   PetscFunctionReturn(PETSC_SUCCESS);
5705: }

5707: /* FEM Support */

5709: PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
5710: {
5711:   PetscInt depth;

5713:   PetscFunctionBegin;
5714:   PetscCall(DMPlexGetDepth(plex, &depth));
5715:   PetscCall(DMGetStratumIS(plex, "dim", depth, cellIS));
5716:   if (!*cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, cellIS));
5717:   PetscFunctionReturn(PETSC_SUCCESS);
5718: }

5720: PetscErrorCode DMPlexGetAllFaces_Internal(DM plex, IS *faceIS)
5721: {
5722:   PetscInt depth;

5724:   PetscFunctionBegin;
5725:   PetscCall(DMPlexGetDepth(plex, &depth));
5726:   PetscCall(DMGetStratumIS(plex, "dim", depth - 1, faceIS));
5727:   if (!*faceIS) PetscCall(DMGetStratumIS(plex, "depth", depth - 1, faceIS));
5728:   PetscFunctionReturn(PETSC_SUCCESS);
5729: }

5731: /*
5732:  Returns number of components and tensor degree for the field.  For interpolated meshes, line should be a point
5733:  representing a line in the section.
5734: */
5735: static PetscErrorCode PetscSectionFieldGetTensorDegree_Private(DM dm, PetscSection section, PetscInt field, PetscInt line, PetscInt *Nc, PetscInt *k, PetscBool *continuous, PetscBool *tensor)
5736: {
5737:   PetscObject  obj;
5738:   PetscClassId id;
5739:   PetscFE      fe = NULL;

5741:   PetscFunctionBeginHot;
5742:   PetscCall(PetscSectionGetFieldComponents(section, field, Nc));
5743:   PetscCall(DMGetField(dm, field, NULL, &obj));
5744:   PetscCall(PetscObjectGetClassId(obj, &id));
5745:   if (id == PETSCFE_CLASSID) fe = (PetscFE)obj;

5747:   if (!fe) {
5748:     /* Assume the full interpolated mesh is in the chart; lines in particular */
5749:     /* An order k SEM disc has k-1 dofs on an edge */
5750:     PetscCall(PetscSectionGetFieldDof(section, line, field, k));
5751:     *k = *k / *Nc + 1;
5752:   } else {
5753:     PetscInt       dual_space_size, dim;
5754:     PetscDualSpace dsp;

5756:     PetscCall(DMGetDimension(dm, &dim));
5757:     PetscCall(PetscFEGetDualSpace(fe, &dsp));
5758:     PetscCall(PetscDualSpaceGetDimension(dsp, &dual_space_size));
5759:     *k = (PetscInt)PetscCeilReal(PetscPowReal(dual_space_size / *Nc, 1.0 / dim)) - 1;
5760:     PetscCall(PetscDualSpaceLagrangeGetContinuity(dsp, continuous));
5761:     PetscCall(PetscDualSpaceLagrangeGetTensor(dsp, tensor));
5762:   }
5763:   PetscFunctionReturn(PETSC_SUCCESS);
5764: }

5766: static PetscErrorCode GetFieldSize_Private(PetscInt dim, PetscInt k, PetscBool tensor, PetscInt *dof)
5767: {
5768:   PetscFunctionBeginHot;
5769:   if (tensor) {
5770:     *dof = PetscPowInt(k + 1, dim);
5771:   } else {
5772:     switch (dim) {
5773:     case 1:
5774:       *dof = k + 1;
5775:       break;
5776:     case 2:
5777:       *dof = ((k + 1) * (k + 2)) / 2;
5778:       break;
5779:     case 3:
5780:       *dof = ((k + 1) * (k + 2) * (k + 3)) / 6;
5781:       break;
5782:     default:
5783:       *dof = 0;
5784:     }
5785:   }
5786:   PetscFunctionReturn(PETSC_SUCCESS);
5787: }

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

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

5799:   Example:
5800:   A typical interpolated single-quad mesh might order points as
5801: .vb
5802:   [c0, v1, v2, v3, v4, e5, e6, e7, e8]

5804:   v4 -- e6 -- v3
5805:   |           |
5806:   e7    c0    e8
5807:   |           |
5808:   v1 -- e5 -- v2
5809: .ve

5811:   (There is no significance to the ordering described here.)  The default section for a Q3 quad might typically assign
5812:   dofs in the order of points, e.g.,
5813: .vb
5814:     c0 -> [0,1,2,3]
5815:     v1 -> [4]
5816:     ...
5817:     e5 -> [8, 9]
5818: .ve

5820:   which corresponds to the dofs
5821: .vb
5822:     6   10  11  7
5823:     13  2   3   15
5824:     12  0   1   14
5825:     4   8   9   5
5826: .ve

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

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

5838:   Level: developer

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

5844:   This is required to run with libCEED.

5846: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMGetLocalSection()`, `PetscSectionSetClosurePermutation()`, `DMSetGlobalSection()`
5847: @*/
5848: PetscErrorCode DMPlexSetClosurePermutationTensor(DM dm, PetscInt point, PetscSection section)
5849: {
5850:   DMLabel   label;
5851:   PetscInt  dim, depth = -1, eStart = -1, Nf;
5852:   PetscBool continuous = PETSC_TRUE, tensor = PETSC_TRUE;

5854:   PetscFunctionBegin;
5855:   PetscCall(DMGetDimension(dm, &dim));
5856:   if (dim < 1) PetscFunctionReturn(PETSC_SUCCESS);
5857:   if (point < 0) {
5858:     PetscInt sStart, sEnd;

5860:     PetscCall(DMPlexGetDepthStratum(dm, 1, &sStart, &sEnd));
5861:     point = sEnd - sStart ? sStart : point;
5862:   }
5863:   PetscCall(DMPlexGetDepthLabel(dm, &label));
5864:   if (point >= 0) PetscCall(DMLabelGetValue(label, point, &depth));
5865:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
5866:   if (depth == 1) {
5867:     eStart = point;
5868:   } else if (depth == dim) {
5869:     const PetscInt *cone;

5871:     PetscCall(DMPlexGetCone(dm, point, &cone));
5872:     if (dim == 2) eStart = cone[0];
5873:     else if (dim == 3) {
5874:       const PetscInt *cone2;
5875:       PetscCall(DMPlexGetCone(dm, cone[0], &cone2));
5876:       eStart = cone2[0];
5877:     } 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);
5878:   } 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);

5880:   PetscCall(PetscSectionGetNumFields(section, &Nf));
5881:   for (PetscInt d = 1; d <= dim; d++) {
5882:     PetscInt  k, f, Nc, c, i, j, size = 0, offset = 0, foffset = 0;
5883:     PetscInt *perm;

5885:     for (f = 0; f < Nf; ++f) {
5886:       PetscInt dof;

5888:       PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
5889:       PetscCheck(dim == 1 || tensor || !continuous, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Continuous field %" PetscInt_FMT " must have a tensor product discretization", f);
5890:       if (!continuous && d < dim) continue;
5891:       PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
5892:       size += dof * Nc;
5893:     }
5894:     PetscCall(PetscMalloc1(size, &perm));
5895:     for (f = 0; f < Nf; ++f) {
5896:       switch (d) {
5897:       case 1:
5898:         PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
5899:         if (!continuous && d < dim) continue;
5900:         /*
5901:          Original ordering is [ edge of length k-1; vtx0; vtx1 ]
5902:          We want              [ vtx0; edge of length k-1; vtx1 ]
5903:          */
5904:         if (continuous) {
5905:           for (c = 0; c < Nc; c++, offset++) perm[offset] = (k - 1) * Nc + c + foffset;
5906:           for (i = 0; i < k - 1; i++)
5907:             for (c = 0; c < Nc; c++, offset++) perm[offset] = i * Nc + c + foffset;
5908:           for (c = 0; c < Nc; c++, offset++) perm[offset] = k * Nc + c + foffset;
5909:           foffset = offset;
5910:         } else {
5911:           PetscInt dof;

5913:           PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
5914:           for (i = 0; i < dof * Nc; ++i, ++offset) perm[offset] = i + foffset;
5915:           foffset = offset;
5916:         }
5917:         break;
5918:       case 2:
5919:         /* The original quad closure is oriented clockwise, {f, e_b, e_r, e_t, e_l, v_lb, v_rb, v_tr, v_tl} */
5920:         PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
5921:         if (!continuous && d < dim) continue;
5922:         /* The SEM order is

5924:          v_lb, {e_b}, v_rb,
5925:          e^{(k-1)-i}_l, {f^{i*(k-1)}}, e^i_r,
5926:          v_lt, reverse {e_t}, v_rt
5927:          */
5928:         if (continuous) {
5929:           const PetscInt of   = 0;
5930:           const PetscInt oeb  = of + PetscSqr(k - 1);
5931:           const PetscInt oer  = oeb + (k - 1);
5932:           const PetscInt oet  = oer + (k - 1);
5933:           const PetscInt oel  = oet + (k - 1);
5934:           const PetscInt ovlb = oel + (k - 1);
5935:           const PetscInt ovrb = ovlb + 1;
5936:           const PetscInt ovrt = ovrb + 1;
5937:           const PetscInt ovlt = ovrt + 1;
5938:           PetscInt       o;

5940:           /* bottom */
5941:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlb * Nc + c + foffset;
5942:           for (o = oeb; o < oer; ++o)
5943:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
5944:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrb * Nc + c + foffset;
5945:           /* middle */
5946:           for (i = 0; i < k - 1; ++i) {
5947:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oel + (k - 2) - i) * Nc + c + foffset;
5948:             for (o = of + (k - 1) * i; o < of + (k - 1) * (i + 1); ++o)
5949:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
5950:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oer + i) * Nc + c + foffset;
5951:           }
5952:           /* top */
5953:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlt * Nc + c + foffset;
5954:           for (o = oel - 1; o >= oet; --o)
5955:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
5956:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrt * Nc + c + foffset;
5957:           foffset = offset;
5958:         } else {
5959:           PetscInt dof;

5961:           PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
5962:           for (i = 0; i < dof * Nc; ++i, ++offset) perm[offset] = i + foffset;
5963:           foffset = offset;
5964:         }
5965:         break;
5966:       case 3:
5967:         /* The original hex closure is

5969:          {c,
5970:          f_b, f_t, f_f, f_b, f_r, f_l,
5971:          e_bl, e_bb, e_br, e_bf,  e_tf, e_tr, e_tb, e_tl,  e_rf, e_lf, e_lb, e_rb,
5972:          v_blf, v_blb, v_brb, v_brf, v_tlf, v_trf, v_trb, v_tlb}
5973:          */
5974:         PetscCall(PetscSectionFieldGetTensorDegree_Private(dm, section, f, eStart, &Nc, &k, &continuous, &tensor));
5975:         if (!continuous && d < dim) continue;
5976:         /* The SEM order is
5977:          Bottom Slice
5978:          v_blf, {e^{(k-1)-n}_bf}, v_brf,
5979:          e^{i}_bl, f^{n*(k-1)+(k-1)-i}_b, e^{(k-1)-i}_br,
5980:          v_blb, {e_bb}, v_brb,

5982:          Middle Slice (j)
5983:          {e^{(k-1)-j}_lf}, {f^{j*(k-1)+n}_f}, e^j_rf,
5984:          f^{i*(k-1)+j}_l, {c^{(j*(k-1) + i)*(k-1)+n}_t}, f^{j*(k-1)+i}_r,
5985:          e^j_lb, {f^{j*(k-1)+(k-1)-n}_b}, e^{(k-1)-j}_rb,

5987:          Top Slice
5988:          v_tlf, {e_tf}, v_trf,
5989:          e^{(k-1)-i}_tl, {f^{i*(k-1)}_t}, e^{i}_tr,
5990:          v_tlb, {e^{(k-1)-n}_tb}, v_trb,
5991:          */
5992:         if (continuous) {
5993:           const PetscInt oc    = 0;
5994:           const PetscInt ofb   = oc + PetscSqr(k - 1) * (k - 1);
5995:           const PetscInt oft   = ofb + PetscSqr(k - 1);
5996:           const PetscInt off   = oft + PetscSqr(k - 1);
5997:           const PetscInt ofk   = off + PetscSqr(k - 1);
5998:           const PetscInt ofr   = ofk + PetscSqr(k - 1);
5999:           const PetscInt ofl   = ofr + PetscSqr(k - 1);
6000:           const PetscInt oebl  = ofl + PetscSqr(k - 1);
6001:           const PetscInt oebb  = oebl + (k - 1);
6002:           const PetscInt oebr  = oebb + (k - 1);
6003:           const PetscInt oebf  = oebr + (k - 1);
6004:           const PetscInt oetf  = oebf + (k - 1);
6005:           const PetscInt oetr  = oetf + (k - 1);
6006:           const PetscInt oetb  = oetr + (k - 1);
6007:           const PetscInt oetl  = oetb + (k - 1);
6008:           const PetscInt oerf  = oetl + (k - 1);
6009:           const PetscInt oelf  = oerf + (k - 1);
6010:           const PetscInt oelb  = oelf + (k - 1);
6011:           const PetscInt oerb  = oelb + (k - 1);
6012:           const PetscInt ovblf = oerb + (k - 1);
6013:           const PetscInt ovblb = ovblf + 1;
6014:           const PetscInt ovbrb = ovblb + 1;
6015:           const PetscInt ovbrf = ovbrb + 1;
6016:           const PetscInt ovtlf = ovbrf + 1;
6017:           const PetscInt ovtrf = ovtlf + 1;
6018:           const PetscInt ovtrb = ovtrf + 1;
6019:           const PetscInt ovtlb = ovtrb + 1;
6020:           PetscInt       o, n;

6022:           /* Bottom Slice */
6023:           /*   bottom */
6024:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblf * Nc + c + foffset;
6025:           for (o = oetf - 1; o >= oebf; --o)
6026:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6027:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrf * Nc + c + foffset;
6028:           /*   middle */
6029:           for (i = 0; i < k - 1; ++i) {
6030:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebl + i) * Nc + c + foffset;
6031:             for (n = 0; n < k - 1; ++n) {
6032:               o = ofb + n * (k - 1) + i;
6033:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6034:             }
6035:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebr + (k - 2) - i) * Nc + c + foffset;
6036:           }
6037:           /*   top */
6038:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblb * Nc + c + foffset;
6039:           for (o = oebb; o < oebr; ++o)
6040:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6041:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrb * Nc + c + foffset;

6043:           /* Middle Slice */
6044:           for (j = 0; j < k - 1; ++j) {
6045:             /*   bottom */
6046:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelf + (k - 2) - j) * Nc + c + foffset;
6047:             for (o = off + j * (k - 1); o < off + (j + 1) * (k - 1); ++o)
6048:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6049:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerf + j) * Nc + c + foffset;
6050:             /*   middle */
6051:             for (i = 0; i < k - 1; ++i) {
6052:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofl + i * (k - 1) + j) * Nc + c + foffset;
6053:               for (n = 0; n < k - 1; ++n)
6054:                 for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oc + (j * (k - 1) + i) * (k - 1) + n) * Nc + c + foffset;
6055:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofr + j * (k - 1) + i) * Nc + c + foffset;
6056:             }
6057:             /*   top */
6058:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelb + j) * Nc + c + foffset;
6059:             for (o = ofk + j * (k - 1) + (k - 2); o >= ofk + j * (k - 1); --o)
6060:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6061:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerb + (k - 2) - j) * Nc + c + foffset;
6062:           }

6064:           /* Top Slice */
6065:           /*   bottom */
6066:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlf * Nc + c + foffset;
6067:           for (o = oetf; o < oetr; ++o)
6068:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6069:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrf * Nc + c + foffset;
6070:           /*   middle */
6071:           for (i = 0; i < k - 1; ++i) {
6072:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetl + (k - 2) - i) * Nc + c + foffset;
6073:             for (n = 0; n < k - 1; ++n)
6074:               for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oft + i * (k - 1) + n) * Nc + c + foffset;
6075:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetr + i) * Nc + c + foffset;
6076:           }
6077:           /*   top */
6078:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlb * Nc + c + foffset;
6079:           for (o = oetl - 1; o >= oetb; --o)
6080:             for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o * Nc + c + foffset;
6081:           for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrb * Nc + c + foffset;

6083:           foffset = offset;
6084:         } else {
6085:           PetscInt dof;

6087:           PetscCall(GetFieldSize_Private(d, k, tensor, &dof));
6088:           for (i = 0; i < dof * Nc; ++i, ++offset) perm[offset] = i + foffset;
6089:           foffset = offset;
6090:         }
6091:         break;
6092:       default:
6093:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No spectral ordering for dimension %" PetscInt_FMT, d);
6094:       }
6095:     }
6096:     PetscCheck(offset == size, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Number of permutation entries %" PetscInt_FMT " != %" PetscInt_FMT, offset, size);
6097:     /* Check permutation */
6098:     {
6099:       PetscInt *check;

6101:       PetscCall(PetscMalloc1(size, &check));
6102:       for (i = 0; i < size; ++i) {
6103:         check[i] = -1;
6104:         PetscCheck(perm[i] >= 0 && perm[i] < size, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid permutation index p[%" PetscInt_FMT "] = %" PetscInt_FMT, i, perm[i]);
6105:       }
6106:       for (i = 0; i < size; ++i) check[perm[i]] = i;
6107:       for (i = 0; i < size; ++i) PetscCheck(check[i] >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing permutation index %" PetscInt_FMT, i);
6108:       PetscCall(PetscFree(check));
6109:     }
6110:     PetscCall(PetscSectionSetClosurePermutation_Internal(section, (PetscObject)dm, d, size, PETSC_OWN_POINTER, perm));
6111:     if (d == dim) { // Add permutation for localized (in case this is a coordinate DM)
6112:       PetscInt *loc_perm;
6113:       PetscCall(PetscMalloc1(size * 2, &loc_perm));
6114:       for (PetscInt i = 0; i < size; i++) {
6115:         loc_perm[i]        = perm[i];
6116:         loc_perm[size + i] = size + perm[i];
6117:       }
6118:       PetscCall(PetscSectionSetClosurePermutation_Internal(section, (PetscObject)dm, d, size * 2, PETSC_OWN_POINTER, loc_perm));
6119:     }
6120:   }
6121:   PetscFunctionReturn(PETSC_SUCCESS);
6122: }

6124: PetscErrorCode DMPlexGetPointDualSpaceFEM(DM dm, PetscInt point, PetscInt field, PetscDualSpace *dspace)
6125: {
6126:   PetscDS  prob;
6127:   PetscInt depth, Nf, h;
6128:   DMLabel  label;

6130:   PetscFunctionBeginHot;
6131:   PetscCall(DMGetDS(dm, &prob));
6132:   Nf      = prob->Nf;
6133:   label   = dm->depthLabel;
6134:   *dspace = NULL;
6135:   if (field < Nf) {
6136:     PetscObject disc = prob->disc[field];

6138:     if (disc->classid == PETSCFE_CLASSID) {
6139:       PetscDualSpace dsp;

6141:       PetscCall(PetscFEGetDualSpace((PetscFE)disc, &dsp));
6142:       PetscCall(DMLabelGetNumValues(label, &depth));
6143:       PetscCall(DMLabelGetValue(label, point, &h));
6144:       h = depth - 1 - h;
6145:       if (h) {
6146:         PetscCall(PetscDualSpaceGetHeightSubspace(dsp, h, dspace));
6147:       } else {
6148:         *dspace = dsp;
6149:       }
6150:     }
6151:   }
6152:   PetscFunctionReturn(PETSC_SUCCESS);
6153: }

6155: static inline PetscErrorCode DMPlexVecGetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
6156: {
6157:   PetscScalar       *array;
6158:   const PetscScalar *vArray;
6159:   const PetscInt    *cone, *coneO;
6160:   PetscInt           pStart, pEnd, p, numPoints, size = 0, offset = 0;

6162:   PetscFunctionBeginHot;
6163:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6164:   PetscCall(DMPlexGetConeSize(dm, point, &numPoints));
6165:   PetscCall(DMPlexGetCone(dm, point, &cone));
6166:   PetscCall(DMPlexGetConeOrientation(dm, point, &coneO));
6167:   if (!values || !*values) {
6168:     if ((point >= pStart) && (point < pEnd)) {
6169:       PetscInt dof;

6171:       PetscCall(PetscSectionGetDof(section, point, &dof));
6172:       size += dof;
6173:     }
6174:     for (p = 0; p < numPoints; ++p) {
6175:       const PetscInt cp = cone[p];
6176:       PetscInt       dof;

6178:       if ((cp < pStart) || (cp >= pEnd)) continue;
6179:       PetscCall(PetscSectionGetDof(section, cp, &dof));
6180:       size += dof;
6181:     }
6182:     if (!values) {
6183:       if (csize) *csize = size;
6184:       PetscFunctionReturn(PETSC_SUCCESS);
6185:     }
6186:     PetscCall(DMGetWorkArray(dm, size, MPIU_SCALAR, &array));
6187:   } else {
6188:     array = *values;
6189:   }
6190:   size = 0;
6191:   PetscCall(VecGetArrayRead(v, &vArray));
6192:   if ((point >= pStart) && (point < pEnd)) {
6193:     PetscInt           dof, off, d;
6194:     const PetscScalar *varr;

6196:     PetscCall(PetscSectionGetDof(section, point, &dof));
6197:     PetscCall(PetscSectionGetOffset(section, point, &off));
6198:     varr = PetscSafePointerPlusOffset(vArray, off);
6199:     for (d = 0; d < dof; ++d, ++offset) array[offset] = varr[d];
6200:     size += dof;
6201:   }
6202:   for (p = 0; p < numPoints; ++p) {
6203:     const PetscInt     cp = cone[p];
6204:     PetscInt           o  = coneO[p];
6205:     PetscInt           dof, off, d;
6206:     const PetscScalar *varr;

6208:     if ((cp < pStart) || (cp >= pEnd)) continue;
6209:     PetscCall(PetscSectionGetDof(section, cp, &dof));
6210:     PetscCall(PetscSectionGetOffset(section, cp, &off));
6211:     varr = PetscSafePointerPlusOffset(vArray, off);
6212:     if (o >= 0) {
6213:       for (d = 0; d < dof; ++d, ++offset) array[offset] = varr[d];
6214:     } else {
6215:       for (d = dof - 1; d >= 0; --d, ++offset) array[offset] = varr[d];
6216:     }
6217:     size += dof;
6218:   }
6219:   PetscCall(VecRestoreArrayRead(v, &vArray));
6220:   if (!*values) {
6221:     if (csize) *csize = size;
6222:     *values = array;
6223:   } else {
6224:     PetscCheck(size <= *csize, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %" PetscInt_FMT " < actual size %" PetscInt_FMT, *csize, size);
6225:     *csize = size;
6226:   }
6227:   PetscFunctionReturn(PETSC_SUCCESS);
6228: }

6230: /* Compress out points not in the section */
6231: static inline PetscErrorCode CompressPoints_Private(PetscSection section, PetscInt *numPoints, PetscInt points[])
6232: {
6233:   const PetscInt np = *numPoints;
6234:   PetscInt       pStart, pEnd, p, q;

6236:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6237:   for (p = 0, q = 0; p < np; ++p) {
6238:     const PetscInt r = points[p * 2];
6239:     if ((r >= pStart) && (r < pEnd)) {
6240:       points[q * 2]     = r;
6241:       points[q * 2 + 1] = points[p * 2 + 1];
6242:       ++q;
6243:     }
6244:   }
6245:   *numPoints = q;
6246:   return PETSC_SUCCESS;
6247: }

6249: /* Compressed closure does not apply closure permutation */
6250: PetscErrorCode DMPlexGetCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt ornt, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
6251: {
6252:   const PetscInt *cla = NULL;
6253:   PetscInt        np, *pts = NULL;

6255:   PetscFunctionBeginHot;
6256:   PetscCall(PetscSectionGetClosureIndex(section, (PetscObject)dm, clSec, clPoints));
6257:   if (!ornt && *clPoints) {
6258:     PetscInt dof, off;

6260:     PetscCall(PetscSectionGetDof(*clSec, point, &dof));
6261:     PetscCall(PetscSectionGetOffset(*clSec, point, &off));
6262:     PetscCall(ISGetIndices(*clPoints, &cla));
6263:     np  = dof / 2;
6264:     pts = PetscSafePointerPlusOffset((PetscInt *)cla, off);
6265:   } else {
6266:     PetscCall(DMPlexGetTransitiveClosure_Internal(dm, point, ornt, PETSC_TRUE, &np, &pts));
6267:     PetscCall(CompressPoints_Private(section, &np, pts));
6268:   }
6269:   *numPoints = np;
6270:   *points    = pts;
6271:   *clp       = cla;
6272:   PetscFunctionReturn(PETSC_SUCCESS);
6273: }

6275: PetscErrorCode DMPlexRestoreCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
6276: {
6277:   PetscFunctionBeginHot;
6278:   if (!*clPoints) {
6279:     PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, numPoints, points));
6280:   } else {
6281:     PetscCall(ISRestoreIndices(*clPoints, clp));
6282:   }
6283:   *numPoints = 0;
6284:   *points    = NULL;
6285:   *clSec     = NULL;
6286:   *clPoints  = NULL;
6287:   *clp       = NULL;
6288:   PetscFunctionReturn(PETSC_SUCCESS);
6289: }

6291: static inline PetscErrorCode DMPlexVecGetClosure_Static(DM dm, PetscSection section, PetscInt numPoints, const PetscInt points[], const PetscInt clperm[], const PetscScalar vArray[], PetscInt *size, PetscScalar array[])
6292: {
6293:   PetscInt            offset = 0, p;
6294:   const PetscInt    **perms  = NULL;
6295:   const PetscScalar **flips  = NULL;

6297:   PetscFunctionBeginHot;
6298:   *size = 0;
6299:   PetscCall(PetscSectionGetPointSyms(section, numPoints, points, &perms, &flips));
6300:   for (p = 0; p < numPoints; p++) {
6301:     const PetscInt     point = points[2 * p];
6302:     const PetscInt    *perm  = perms ? perms[p] : NULL;
6303:     const PetscScalar *flip  = flips ? flips[p] : NULL;
6304:     PetscInt           dof, off, d;
6305:     const PetscScalar *varr;

6307:     PetscCall(PetscSectionGetDof(section, point, &dof));
6308:     PetscCall(PetscSectionGetOffset(section, point, &off));
6309:     varr = PetscSafePointerPlusOffset(vArray, off);
6310:     if (clperm) {
6311:       if (perm) {
6312:         for (d = 0; d < dof; d++) array[clperm[offset + perm[d]]] = varr[d];
6313:       } else {
6314:         for (d = 0; d < dof; d++) array[clperm[offset + d]] = varr[d];
6315:       }
6316:       if (flip) {
6317:         for (d = 0; d < dof; d++) array[clperm[offset + d]] *= flip[d];
6318:       }
6319:     } else {
6320:       if (perm) {
6321:         for (d = 0; d < dof; d++) array[offset + perm[d]] = varr[d];
6322:       } else {
6323:         for (d = 0; d < dof; d++) array[offset + d] = varr[d];
6324:       }
6325:       if (flip) {
6326:         for (d = 0; d < dof; d++) array[offset + d] *= flip[d];
6327:       }
6328:     }
6329:     offset += dof;
6330:   }
6331:   PetscCall(PetscSectionRestorePointSyms(section, numPoints, points, &perms, &flips));
6332:   *size = offset;
6333:   PetscFunctionReturn(PETSC_SUCCESS);
6334: }

6336: 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[])
6337: {
6338:   PetscInt offset = 0, f;

6340:   PetscFunctionBeginHot;
6341:   *size = 0;
6342:   for (f = 0; f < numFields; ++f) {
6343:     PetscInt            p;
6344:     const PetscInt    **perms = NULL;
6345:     const PetscScalar **flips = NULL;

6347:     PetscCall(PetscSectionGetFieldPointSyms(section, f, numPoints, points, &perms, &flips));
6348:     for (p = 0; p < numPoints; p++) {
6349:       const PetscInt     point = points[2 * p];
6350:       PetscInt           fdof, foff, b;
6351:       const PetscScalar *varr;
6352:       const PetscInt    *perm = perms ? perms[p] : NULL;
6353:       const PetscScalar *flip = flips ? flips[p] : NULL;

6355:       PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
6356:       PetscCall(PetscSectionGetFieldOffset(section, point, f, &foff));
6357:       varr = &vArray[foff];
6358:       if (clperm) {
6359:         if (perm) {
6360:           for (b = 0; b < fdof; b++) array[clperm[offset + perm[b]]] = varr[b];
6361:         } else {
6362:           for (b = 0; b < fdof; b++) array[clperm[offset + b]] = varr[b];
6363:         }
6364:         if (flip) {
6365:           for (b = 0; b < fdof; b++) array[clperm[offset + b]] *= flip[b];
6366:         }
6367:       } else {
6368:         if (perm) {
6369:           for (b = 0; b < fdof; b++) array[offset + perm[b]] = varr[b];
6370:         } else {
6371:           for (b = 0; b < fdof; b++) array[offset + b] = varr[b];
6372:         }
6373:         if (flip) {
6374:           for (b = 0; b < fdof; b++) array[offset + b] *= flip[b];
6375:         }
6376:       }
6377:       offset += fdof;
6378:     }
6379:     PetscCall(PetscSectionRestoreFieldPointSyms(section, f, numPoints, points, &perms, &flips));
6380:   }
6381:   *size = offset;
6382:   PetscFunctionReturn(PETSC_SUCCESS);
6383: }

6385: PetscErrorCode DMPlexVecGetOrientedClosure_Internal(DM dm, PetscSection section, PetscBool useClPerm, Vec v, PetscInt point, PetscInt ornt, PetscInt *csize, PetscScalar *values[])
6386: {
6387:   PetscSection    clSection;
6388:   IS              clPoints;
6389:   PetscInt       *points = NULL;
6390:   const PetscInt *clp, *perm = NULL;
6391:   PetscInt        depth, numFields, numPoints, asize;

6393:   PetscFunctionBeginHot;
6395:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
6398:   PetscCall(DMPlexGetDepth(dm, &depth));
6399:   PetscCall(PetscSectionGetNumFields(section, &numFields));
6400:   if (depth == 1 && numFields < 2) {
6401:     PetscCall(DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values));
6402:     PetscFunctionReturn(PETSC_SUCCESS);
6403:   }
6404:   /* Get points */
6405:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, ornt, &numPoints, &points, &clSection, &clPoints, &clp));
6406:   /* Get sizes */
6407:   asize = 0;
6408:   for (PetscInt p = 0; p < numPoints * 2; p += 2) {
6409:     PetscInt dof;
6410:     PetscCall(PetscSectionGetDof(section, points[p], &dof));
6411:     asize += dof;
6412:   }
6413:   if (values) {
6414:     const PetscScalar *vArray;
6415:     PetscInt           size;

6417:     if (*values) {
6418:       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);
6419:     } else PetscCall(DMGetWorkArray(dm, asize, MPIU_SCALAR, values));
6420:     if (useClPerm) PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, asize, &perm));
6421:     PetscCall(VecGetArrayRead(v, &vArray));
6422:     /* Get values */
6423:     if (numFields > 0) PetscCall(DMPlexVecGetClosure_Fields_Static(dm, section, numPoints, points, numFields, perm, vArray, &size, *values));
6424:     else PetscCall(DMPlexVecGetClosure_Static(dm, section, numPoints, points, perm, vArray, &size, *values));
6425:     PetscCheck(asize == size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Section size %" PetscInt_FMT " does not match Vec closure size %" PetscInt_FMT, asize, size);
6426:     /* Cleanup array */
6427:     PetscCall(VecRestoreArrayRead(v, &vArray));
6428:   }
6429:   if (csize) *csize = asize;
6430:   /* Cleanup points */
6431:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
6432:   PetscFunctionReturn(PETSC_SUCCESS);
6433: }

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

6438:   Not collective

6440:   Input Parameters:
6441: + dm      - The `DM`
6442: . section - The section describing the layout in `v`, or `NULL` to use the default section
6443: . v       - The local vector
6444: - point   - The point in the `DM`

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

6451:   Level: intermediate

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

6458:   A typical use could be
6459: .vb
6460:    values = NULL;
6461:    PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values));
6462:    for (cl = 0; cl < clSize; ++cl) {
6463:      <Compute on closure>
6464:    }
6465:    PetscCall(DMPlexVecRestoreClosure(dm, NULL, v, p, &clSize, &values));
6466: .ve
6467:   or
6468: .vb
6469:    PetscMalloc1(clMaxSize, &values);
6470:    for (p = pStart; p < pEnd; ++p) {
6471:      clSize = clMaxSize;
6472:      PetscCall(DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values));
6473:      for (cl = 0; cl < clSize; ++cl) {
6474:        <Compute on closure>
6475:      }
6476:    }
6477:    PetscFree(values);
6478: .ve

6480:   Fortran Notes:
6481:   The `csize` argument is not present in the Fortran binding.

6483:   `values` must be declared with
6484: .vb
6485:   PetscScalar,dimension(:),pointer   :: values
6486: .ve
6487:   and it will be allocated internally by PETSc to hold the values returned

6489: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecRestoreClosure()`, `DMPlexVecSetClosure()`, `DMPlexMatSetClosure()`
6490: @*/
6491: PetscErrorCode DMPlexVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
6492: {
6493:   PetscFunctionBeginHot;
6494:   PetscCall(DMPlexVecGetOrientedClosure_Internal(dm, section, PETSC_TRUE, v, point, 0, csize, values));
6495:   PetscFunctionReturn(PETSC_SUCCESS);
6496: }

6498: PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt depth, PetscInt *csize, PetscScalar *values[])
6499: {
6500:   DMLabel            depthLabel;
6501:   PetscSection       clSection;
6502:   IS                 clPoints;
6503:   PetscScalar       *array;
6504:   const PetscScalar *vArray;
6505:   PetscInt          *points = NULL;
6506:   const PetscInt    *clp, *perm = NULL;
6507:   PetscInt           mdepth, numFields, numPoints, Np = 0, p, clsize, size;

6509:   PetscFunctionBeginHot;
6511:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
6514:   PetscCall(DMPlexGetDepth(dm, &mdepth));
6515:   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
6516:   PetscCall(PetscSectionGetNumFields(section, &numFields));
6517:   if (mdepth == 1 && numFields < 2) {
6518:     PetscCall(DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values));
6519:     PetscFunctionReturn(PETSC_SUCCESS);
6520:   }
6521:   /* Get points */
6522:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &numPoints, &points, &clSection, &clPoints, &clp));
6523:   for (clsize = 0, p = 0; p < Np; p++) {
6524:     PetscInt dof;
6525:     PetscCall(PetscSectionGetDof(section, points[2 * p], &dof));
6526:     clsize += dof;
6527:   }
6528:   PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, clsize, &perm));
6529:   /* Filter points */
6530:   for (p = 0; p < numPoints * 2; p += 2) {
6531:     PetscInt dep;

6533:     PetscCall(DMLabelGetValue(depthLabel, points[p], &dep));
6534:     if (dep != depth) continue;
6535:     points[Np * 2 + 0] = points[p];
6536:     points[Np * 2 + 1] = points[p + 1];
6537:     ++Np;
6538:   }
6539:   /* Get array */
6540:   if (!values || !*values) {
6541:     PetscInt asize = 0, dof;

6543:     for (p = 0; p < Np * 2; p += 2) {
6544:       PetscCall(PetscSectionGetDof(section, points[p], &dof));
6545:       asize += dof;
6546:     }
6547:     if (!values) {
6548:       PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
6549:       if (csize) *csize = asize;
6550:       PetscFunctionReturn(PETSC_SUCCESS);
6551:     }
6552:     PetscCall(DMGetWorkArray(dm, asize, MPIU_SCALAR, &array));
6553:   } else {
6554:     array = *values;
6555:   }
6556:   PetscCall(VecGetArrayRead(v, &vArray));
6557:   /* Get values */
6558:   if (numFields > 0) PetscCall(DMPlexVecGetClosure_Fields_Static(dm, section, Np, points, numFields, perm, vArray, &size, array));
6559:   else PetscCall(DMPlexVecGetClosure_Static(dm, section, Np, points, perm, vArray, &size, array));
6560:   /* Cleanup points */
6561:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
6562:   /* Cleanup array */
6563:   PetscCall(VecRestoreArrayRead(v, &vArray));
6564:   if (!*values) {
6565:     if (csize) *csize = size;
6566:     *values = array;
6567:   } else {
6568:     PetscCheck(size <= *csize, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %" PetscInt_FMT " < actual size %" PetscInt_FMT, *csize, size);
6569:     *csize = size;
6570:   }
6571:   PetscFunctionReturn(PETSC_SUCCESS);
6572: }

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

6577:   Not collective

6579:   Input Parameters:
6580: + dm      - The `DM`
6581: . section - The section describing the layout in `v`, or `NULL` to use the default section
6582: . v       - The local vector
6583: . point   - The point in the `DM`
6584: . csize   - The number of values in the closure, or `NULL`
6585: - values  - The array of values

6587:   Level: intermediate

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

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

6595: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecGetClosure()`, `DMPlexVecSetClosure()`, `DMPlexMatSetClosure()`
6596: @*/
6597: PetscErrorCode DMPlexVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
6598: {
6599:   PetscInt size = 0;

6601:   PetscFunctionBegin;
6602:   /* Should work without recalculating size */
6603:   PetscCall(DMRestoreWorkArray(dm, size, MPIU_SCALAR, (void *)values));
6604:   *values = NULL;
6605:   PetscFunctionReturn(PETSC_SUCCESS);
6606: }

6608: static inline void add(PetscScalar *x, PetscScalar y)
6609: {
6610:   *x += y;
6611: }
6612: static inline void insert(PetscScalar *x, PetscScalar y)
6613: {
6614:   *x = y;
6615: }

6617: 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[])
6618: {
6619:   PetscInt        cdof;  /* The number of constraints on this point */
6620:   const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6621:   PetscScalar    *a;
6622:   PetscInt        off, cind = 0, k;

6624:   PetscFunctionBegin;
6625:   PetscCall(PetscSectionGetConstraintDof(section, point, &cdof));
6626:   PetscCall(PetscSectionGetOffset(section, point, &off));
6627:   a = &array[off];
6628:   if (!cdof || setBC) {
6629:     if (clperm) {
6630:       if (perm) {
6631:         for (k = 0; k < dof; ++k) fuse(&a[k], values[clperm[offset + perm[k]]] * (flip ? flip[perm[k]] : 1.));
6632:       } else {
6633:         for (k = 0; k < dof; ++k) fuse(&a[k], values[clperm[offset + k]] * (flip ? flip[k] : 1.));
6634:       }
6635:     } else {
6636:       if (perm) {
6637:         for (k = 0; k < dof; ++k) fuse(&a[k], values[offset + perm[k]] * (flip ? flip[perm[k]] : 1.));
6638:       } else {
6639:         for (k = 0; k < dof; ++k) fuse(&a[k], values[offset + k] * (flip ? flip[k] : 1.));
6640:       }
6641:     }
6642:   } else {
6643:     PetscCall(PetscSectionGetConstraintIndices(section, point, &cdofs));
6644:     if (clperm) {
6645:       if (perm) {
6646:         for (k = 0; k < dof; ++k) {
6647:           if ((cind < cdof) && (k == cdofs[cind])) {
6648:             ++cind;
6649:             continue;
6650:           }
6651:           fuse(&a[k], values[clperm[offset + perm[k]]] * (flip ? flip[perm[k]] : 1.));
6652:         }
6653:       } else {
6654:         for (k = 0; k < dof; ++k) {
6655:           if ((cind < cdof) && (k == cdofs[cind])) {
6656:             ++cind;
6657:             continue;
6658:           }
6659:           fuse(&a[k], values[clperm[offset + k]] * (flip ? flip[k] : 1.));
6660:         }
6661:       }
6662:     } else {
6663:       if (perm) {
6664:         for (k = 0; k < dof; ++k) {
6665:           if ((cind < cdof) && (k == cdofs[cind])) {
6666:             ++cind;
6667:             continue;
6668:           }
6669:           fuse(&a[k], values[offset + perm[k]] * (flip ? flip[perm[k]] : 1.));
6670:         }
6671:       } else {
6672:         for (k = 0; k < dof; ++k) {
6673:           if ((cind < cdof) && (k == cdofs[cind])) {
6674:             ++cind;
6675:             continue;
6676:           }
6677:           fuse(&a[k], values[offset + k] * (flip ? flip[k] : 1.));
6678:         }
6679:       }
6680:     }
6681:   }
6682:   PetscFunctionReturn(PETSC_SUCCESS);
6683: }

6685: 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[])
6686: {
6687:   PetscInt        cdof;  /* The number of constraints on this point */
6688:   const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6689:   PetscScalar    *a;
6690:   PetscInt        off, cind = 0, k;

6692:   PetscFunctionBegin;
6693:   PetscCall(PetscSectionGetConstraintDof(section, point, &cdof));
6694:   PetscCall(PetscSectionGetOffset(section, point, &off));
6695:   a = &array[off];
6696:   if (cdof) {
6697:     PetscCall(PetscSectionGetConstraintIndices(section, point, &cdofs));
6698:     if (clperm) {
6699:       if (perm) {
6700:         for (k = 0; k < dof; ++k) {
6701:           if ((cind < cdof) && (k == cdofs[cind])) {
6702:             fuse(&a[k], values[clperm[offset + perm[k]]] * (flip ? flip[perm[k]] : 1.));
6703:             cind++;
6704:           }
6705:         }
6706:       } else {
6707:         for (k = 0; k < dof; ++k) {
6708:           if ((cind < cdof) && (k == cdofs[cind])) {
6709:             fuse(&a[k], values[clperm[offset + k]] * (flip ? flip[k] : 1.));
6710:             cind++;
6711:           }
6712:         }
6713:       }
6714:     } else {
6715:       if (perm) {
6716:         for (k = 0; k < dof; ++k) {
6717:           if ((cind < cdof) && (k == cdofs[cind])) {
6718:             fuse(&a[k], values[offset + perm[k]] * (flip ? flip[perm[k]] : 1.));
6719:             cind++;
6720:           }
6721:         }
6722:       } else {
6723:         for (k = 0; k < dof; ++k) {
6724:           if ((cind < cdof) && (k == cdofs[cind])) {
6725:             fuse(&a[k], values[offset + k] * (flip ? flip[k] : 1.));
6726:             cind++;
6727:           }
6728:         }
6729:       }
6730:     }
6731:   }
6732:   PetscFunctionReturn(PETSC_SUCCESS);
6733: }

6735: 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[])
6736: {
6737:   PetscScalar    *a;
6738:   PetscInt        fdof, foff, fcdof, foffset = *offset;
6739:   const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
6740:   PetscInt        cind = 0, b;

6742:   PetscFunctionBegin;
6743:   PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
6744:   PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &fcdof));
6745:   PetscCall(PetscSectionGetFieldOffset(section, point, f, &foff));
6746:   a = &array[foff];
6747:   if (!fcdof || setBC) {
6748:     if (clperm) {
6749:       if (perm) {
6750:         for (b = 0; b < fdof; b++) fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6751:       } else {
6752:         for (b = 0; b < fdof; b++) fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
6753:       }
6754:     } else {
6755:       if (perm) {
6756:         for (b = 0; b < fdof; b++) fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
6757:       } else {
6758:         for (b = 0; b < fdof; b++) fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
6759:       }
6760:     }
6761:   } else {
6762:     PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
6763:     if (clperm) {
6764:       if (perm) {
6765:         for (b = 0; b < fdof; b++) {
6766:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6767:             ++cind;
6768:             continue;
6769:           }
6770:           fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6771:         }
6772:       } else {
6773:         for (b = 0; b < fdof; b++) {
6774:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6775:             ++cind;
6776:             continue;
6777:           }
6778:           fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
6779:         }
6780:       }
6781:     } else {
6782:       if (perm) {
6783:         for (b = 0; b < fdof; b++) {
6784:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6785:             ++cind;
6786:             continue;
6787:           }
6788:           fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
6789:         }
6790:       } else {
6791:         for (b = 0; b < fdof; b++) {
6792:           if ((cind < fcdof) && (b == fcdofs[cind])) {
6793:             ++cind;
6794:             continue;
6795:           }
6796:           fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
6797:         }
6798:       }
6799:     }
6800:   }
6801:   *offset += fdof;
6802:   PetscFunctionReturn(PETSC_SUCCESS);
6803: }

6805: 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[])
6806: {
6807:   PetscScalar    *a;
6808:   PetscInt        fdof, foff, fcdof, foffset = *offset;
6809:   const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
6810:   PetscInt        Nc, cind = 0, ncind = 0, b;
6811:   PetscBool       ncSet, fcSet;

6813:   PetscFunctionBegin;
6814:   PetscCall(PetscSectionGetFieldComponents(section, f, &Nc));
6815:   PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
6816:   PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &fcdof));
6817:   PetscCall(PetscSectionGetFieldOffset(section, point, f, &foff));
6818:   a = &array[foff];
6819:   if (fcdof) {
6820:     /* We just override fcdof and fcdofs with Ncc and comps */
6821:     PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
6822:     if (clperm) {
6823:       if (perm) {
6824:         if (comps) {
6825:           for (b = 0; b < fdof; b++) {
6826:             ncSet = fcSet = PETSC_FALSE;
6827:             if (b % Nc == comps[ncind]) {
6828:               ncind = (ncind + 1) % Ncc;
6829:               ncSet = PETSC_TRUE;
6830:             }
6831:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6832:               ++cind;
6833:               fcSet = PETSC_TRUE;
6834:             }
6835:             if (ncSet && fcSet) fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6836:           }
6837:         } else {
6838:           for (b = 0; b < fdof; b++) {
6839:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6840:               fuse(&a[b], values[clperm[foffset + perm[b]]] * (flip ? flip[perm[b]] : 1.));
6841:               ++cind;
6842:             }
6843:           }
6844:         }
6845:       } else {
6846:         if (comps) {
6847:           for (b = 0; b < fdof; b++) {
6848:             ncSet = fcSet = PETSC_FALSE;
6849:             if (b % Nc == comps[ncind]) {
6850:               ncind = (ncind + 1) % Ncc;
6851:               ncSet = PETSC_TRUE;
6852:             }
6853:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6854:               ++cind;
6855:               fcSet = PETSC_TRUE;
6856:             }
6857:             if (ncSet && fcSet) fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
6858:           }
6859:         } else {
6860:           for (b = 0; b < fdof; b++) {
6861:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6862:               fuse(&a[b], values[clperm[foffset + b]] * (flip ? flip[b] : 1.));
6863:               ++cind;
6864:             }
6865:           }
6866:         }
6867:       }
6868:     } else {
6869:       if (perm) {
6870:         if (comps) {
6871:           for (b = 0; b < fdof; b++) {
6872:             ncSet = fcSet = PETSC_FALSE;
6873:             if (b % Nc == comps[ncind]) {
6874:               ncind = (ncind + 1) % Ncc;
6875:               ncSet = PETSC_TRUE;
6876:             }
6877:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6878:               ++cind;
6879:               fcSet = PETSC_TRUE;
6880:             }
6881:             if (ncSet && fcSet) fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
6882:           }
6883:         } else {
6884:           for (b = 0; b < fdof; b++) {
6885:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6886:               fuse(&a[b], values[foffset + perm[b]] * (flip ? flip[perm[b]] : 1.));
6887:               ++cind;
6888:             }
6889:           }
6890:         }
6891:       } else {
6892:         if (comps) {
6893:           for (b = 0; b < fdof; b++) {
6894:             ncSet = fcSet = PETSC_FALSE;
6895:             if (b % Nc == comps[ncind]) {
6896:               ncind = (ncind + 1) % Ncc;
6897:               ncSet = PETSC_TRUE;
6898:             }
6899:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6900:               ++cind;
6901:               fcSet = PETSC_TRUE;
6902:             }
6903:             if (ncSet && fcSet) fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
6904:           }
6905:         } else {
6906:           for (b = 0; b < fdof; b++) {
6907:             if ((cind < fcdof) && (b == fcdofs[cind])) {
6908:               fuse(&a[b], values[foffset + b] * (flip ? flip[b] : 1.));
6909:               ++cind;
6910:             }
6911:           }
6912:         }
6913:       }
6914:     }
6915:   }
6916:   *offset += fdof;
6917:   PetscFunctionReturn(PETSC_SUCCESS);
6918: }

6920: static inline PetscErrorCode DMPlexVecSetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
6921: {
6922:   PetscScalar    *array;
6923:   const PetscInt *cone, *coneO;
6924:   PetscInt        pStart, pEnd, p, numPoints, off, dof;

6926:   PetscFunctionBeginHot;
6927:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6928:   PetscCall(DMPlexGetConeSize(dm, point, &numPoints));
6929:   PetscCall(DMPlexGetCone(dm, point, &cone));
6930:   PetscCall(DMPlexGetConeOrientation(dm, point, &coneO));
6931:   PetscCall(VecGetArray(v, &array));
6932:   for (p = 0, off = 0; p <= numPoints; ++p, off += dof) {
6933:     const PetscInt cp = !p ? point : cone[p - 1];
6934:     const PetscInt o  = !p ? 0 : coneO[p - 1];

6936:     if ((cp < pStart) || (cp >= pEnd)) {
6937:       dof = 0;
6938:       continue;
6939:     }
6940:     PetscCall(PetscSectionGetDof(section, cp, &dof));
6941:     /* ADD_VALUES */
6942:     {
6943:       const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6944:       PetscScalar    *a;
6945:       PetscInt        cdof, coff, cind = 0, k;

6947:       PetscCall(PetscSectionGetConstraintDof(section, cp, &cdof));
6948:       PetscCall(PetscSectionGetOffset(section, cp, &coff));
6949:       a = &array[coff];
6950:       if (!cdof) {
6951:         if (o >= 0) {
6952:           for (k = 0; k < dof; ++k) a[k] += values[off + k];
6953:         } else {
6954:           for (k = 0; k < dof; ++k) a[k] += values[off + dof - k - 1];
6955:         }
6956:       } else {
6957:         PetscCall(PetscSectionGetConstraintIndices(section, cp, &cdofs));
6958:         if (o >= 0) {
6959:           for (k = 0; k < dof; ++k) {
6960:             if ((cind < cdof) && (k == cdofs[cind])) {
6961:               ++cind;
6962:               continue;
6963:             }
6964:             a[k] += values[off + k];
6965:           }
6966:         } else {
6967:           for (k = 0; k < dof; ++k) {
6968:             if ((cind < cdof) && (k == cdofs[cind])) {
6969:               ++cind;
6970:               continue;
6971:             }
6972:             a[k] += values[off + dof - k - 1];
6973:           }
6974:         }
6975:       }
6976:     }
6977:   }
6978:   PetscCall(VecRestoreArray(v, &array));
6979:   PetscFunctionReturn(PETSC_SUCCESS);
6980: }

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

6985:   Not collective

6987:   Input Parameters:
6988: + dm      - The `DM`
6989: . section - The section describing the layout in `v`, or `NULL` to use the default section
6990: . v       - The local vector
6991: . point   - The point in the `DM`
6992: . values  - The array of values
6993: - mode    - The insert mode. One of `INSERT_ALL_VALUES`, `ADD_ALL_VALUES`, `INSERT_VALUES`, `ADD_VALUES`, `INSERT_BC_VALUES`, and `ADD_BC_VALUES`,
6994:             where `INSERT_ALL_VALUES` and `ADD_ALL_VALUES` also overwrite boundary conditions.

6996:   Level: intermediate

6998:   Note:
6999:   Usually the input arrays were obtained with `DMPlexVecGetClosure()`

7001:   Fortran Note:
7002:   `values` must be declared with
7003: .vb
7004:   PetscScalar,dimension(:),pointer   :: values
7005: .ve

7007: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexVecGetClosure()`, `DMPlexMatSetClosure()`
7008: @*/
7009: PetscErrorCode DMPlexVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
7010: {
7011:   PetscSection    clSection;
7012:   IS              clPoints;
7013:   PetscScalar    *array;
7014:   PetscInt       *points = NULL;
7015:   const PetscInt *clp, *clperm = NULL;
7016:   PetscInt        depth, numFields, numPoints, p, clsize;

7018:   PetscFunctionBeginHot;
7020:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
7023:   PetscCall(DMPlexGetDepth(dm, &depth));
7024:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7025:   if (depth == 1 && numFields < 2 && mode == ADD_VALUES) {
7026:     PetscCall(DMPlexVecSetClosure_Depth1_Static(dm, section, v, point, values, mode));
7027:     PetscFunctionReturn(PETSC_SUCCESS);
7028:   }
7029:   /* Get points */
7030:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &numPoints, &points, &clSection, &clPoints, &clp));
7031:   for (clsize = 0, p = 0; p < numPoints; p++) {
7032:     PetscInt dof;
7033:     PetscCall(PetscSectionGetDof(section, points[2 * p], &dof));
7034:     clsize += dof;
7035:   }
7036:   PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, clsize, &clperm));
7037:   /* Get array */
7038:   PetscCall(VecGetArray(v, &array));
7039:   /* Get values */
7040:   if (numFields > 0) {
7041:     PetscInt offset = 0, f;
7042:     for (f = 0; f < numFields; ++f) {
7043:       const PetscInt    **perms = NULL;
7044:       const PetscScalar **flips = NULL;

7046:       PetscCall(PetscSectionGetFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7047:       switch (mode) {
7048:       case INSERT_VALUES:
7049:         for (p = 0; p < numPoints; p++) {
7050:           const PetscInt     point = points[2 * p];
7051:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7052:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7053:           PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, clperm, values, &offset, array));
7054:         }
7055:         break;
7056:       case INSERT_ALL_VALUES:
7057:         for (p = 0; p < numPoints; p++) {
7058:           const PetscInt     point = points[2 * p];
7059:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7060:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7061:           PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, clperm, values, &offset, array));
7062:         }
7063:         break;
7064:       case INSERT_BC_VALUES:
7065:         for (p = 0; p < numPoints; p++) {
7066:           const PetscInt     point = points[2 * p];
7067:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7068:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7069:           PetscCall(updatePointFieldsBC_private(section, point, perm, flip, f, -1, NULL, insert, clperm, values, &offset, array));
7070:         }
7071:         break;
7072:       case ADD_VALUES:
7073:         for (p = 0; p < numPoints; p++) {
7074:           const PetscInt     point = points[2 * p];
7075:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7076:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7077:           PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, clperm, values, &offset, array));
7078:         }
7079:         break;
7080:       case ADD_ALL_VALUES:
7081:         for (p = 0; p < numPoints; p++) {
7082:           const PetscInt     point = points[2 * p];
7083:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7084:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7085:           PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, clperm, values, &offset, array));
7086:         }
7087:         break;
7088:       case ADD_BC_VALUES:
7089:         for (p = 0; p < numPoints; p++) {
7090:           const PetscInt     point = points[2 * p];
7091:           const PetscInt    *perm  = perms ? perms[p] : NULL;
7092:           const PetscScalar *flip  = flips ? flips[p] : NULL;
7093:           PetscCall(updatePointFieldsBC_private(section, point, perm, flip, f, -1, NULL, add, clperm, values, &offset, array));
7094:         }
7095:         break;
7096:       default:
7097:         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
7098:       }
7099:       PetscCall(PetscSectionRestoreFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7100:     }
7101:   } else {
7102:     PetscInt            dof, off;
7103:     const PetscInt    **perms = NULL;
7104:     const PetscScalar **flips = NULL;

7106:     PetscCall(PetscSectionGetPointSyms(section, numPoints, points, &perms, &flips));
7107:     switch (mode) {
7108:     case INSERT_VALUES:
7109:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7110:         const PetscInt     point = points[2 * p];
7111:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7112:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7113:         PetscCall(PetscSectionGetDof(section, point, &dof));
7114:         PetscCall(updatePoint_private(section, point, dof, insert, PETSC_FALSE, perm, flip, clperm, values, off, array));
7115:       }
7116:       break;
7117:     case INSERT_ALL_VALUES:
7118:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7119:         const PetscInt     point = points[2 * p];
7120:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7121:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7122:         PetscCall(PetscSectionGetDof(section, point, &dof));
7123:         PetscCall(updatePoint_private(section, point, dof, insert, PETSC_TRUE, perm, flip, clperm, values, off, array));
7124:       }
7125:       break;
7126:     case INSERT_BC_VALUES:
7127:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7128:         const PetscInt     point = points[2 * p];
7129:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7130:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7131:         PetscCall(PetscSectionGetDof(section, point, &dof));
7132:         PetscCall(updatePointBC_private(section, point, dof, insert, perm, flip, clperm, values, off, array));
7133:       }
7134:       break;
7135:     case ADD_VALUES:
7136:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7137:         const PetscInt     point = points[2 * p];
7138:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7139:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7140:         PetscCall(PetscSectionGetDof(section, point, &dof));
7141:         PetscCall(updatePoint_private(section, point, dof, add, PETSC_FALSE, perm, flip, clperm, values, off, array));
7142:       }
7143:       break;
7144:     case ADD_ALL_VALUES:
7145:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7146:         const PetscInt     point = points[2 * p];
7147:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7148:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7149:         PetscCall(PetscSectionGetDof(section, point, &dof));
7150:         PetscCall(updatePoint_private(section, point, dof, add, PETSC_TRUE, perm, flip, clperm, values, off, array));
7151:       }
7152:       break;
7153:     case ADD_BC_VALUES:
7154:       for (p = 0, off = 0; p < numPoints; p++, off += dof) {
7155:         const PetscInt     point = points[2 * p];
7156:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7157:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7158:         PetscCall(PetscSectionGetDof(section, point, &dof));
7159:         PetscCall(updatePointBC_private(section, point, dof, add, perm, flip, clperm, values, off, array));
7160:       }
7161:       break;
7162:     default:
7163:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
7164:     }
7165:     PetscCall(PetscSectionRestorePointSyms(section, numPoints, points, &perms, &flips));
7166:   }
7167:   /* Cleanup points */
7168:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
7169:   /* Cleanup array */
7170:   PetscCall(VecRestoreArray(v, &array));
7171:   PetscFunctionReturn(PETSC_SUCCESS);
7172: }

7174: /* Check whether the given point is in the label. If not, update the offset to skip this point */
7175: static inline PetscErrorCode CheckPoint_Private(DMLabel label, PetscInt labelId, PetscSection section, PetscInt point, PetscInt f, PetscInt *offset, PetscBool *contains)
7176: {
7177:   PetscFunctionBegin;
7178:   *contains = PETSC_TRUE;
7179:   if (label) {
7180:     PetscInt fdof;

7182:     PetscCall(DMLabelStratumHasPoint(label, labelId, point, contains));
7183:     if (!*contains) {
7184:       PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
7185:       *offset += fdof;
7186:       PetscFunctionReturn(PETSC_SUCCESS);
7187:     }
7188:   }
7189:   PetscFunctionReturn(PETSC_SUCCESS);
7190: }

7192: /* Unlike DMPlexVecSetClosure(), this uses plex-native closure permutation, not a user-specified permutation such as DMPlexSetClosurePermutationTensor(). */
7193: 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)
7194: {
7195:   PetscSection    clSection;
7196:   IS              clPoints;
7197:   PetscScalar    *array;
7198:   PetscInt       *points = NULL;
7199:   const PetscInt *clp;
7200:   PetscInt        numFields, numPoints, p;
7201:   PetscInt        offset = 0, f;

7203:   PetscFunctionBeginHot;
7205:   if (!section) PetscCall(DMGetLocalSection(dm, &section));
7208:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7209:   /* Get points */
7210:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &numPoints, &points, &clSection, &clPoints, &clp));
7211:   /* Get array */
7212:   PetscCall(VecGetArray(v, &array));
7213:   /* Get values */
7214:   for (f = 0; f < numFields; ++f) {
7215:     const PetscInt    **perms = NULL;
7216:     const PetscScalar **flips = NULL;
7217:     PetscBool           contains;

7219:     if (!fieldActive[f]) {
7220:       for (p = 0; p < numPoints * 2; p += 2) {
7221:         PetscInt fdof;
7222:         PetscCall(PetscSectionGetFieldDof(section, points[p], f, &fdof));
7223:         offset += fdof;
7224:       }
7225:       continue;
7226:     }
7227:     PetscCall(PetscSectionGetFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7228:     switch (mode) {
7229:     case INSERT_VALUES:
7230:       for (p = 0; p < numPoints; p++) {
7231:         const PetscInt     point = points[2 * p];
7232:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7233:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7234:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7235:         if (!contains) continue;
7236:         PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, NULL, values, &offset, array));
7237:       }
7238:       break;
7239:     case INSERT_ALL_VALUES:
7240:       for (p = 0; p < numPoints; p++) {
7241:         const PetscInt     point = points[2 * p];
7242:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7243:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7244:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7245:         if (!contains) continue;
7246:         PetscCall(updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, NULL, values, &offset, array));
7247:       }
7248:       break;
7249:     case INSERT_BC_VALUES:
7250:       for (p = 0; p < numPoints; p++) {
7251:         const PetscInt     point = points[2 * p];
7252:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7253:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7254:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7255:         if (!contains) continue;
7256:         PetscCall(updatePointFieldsBC_private(section, point, perm, flip, f, Ncc, comps, insert, NULL, values, &offset, array));
7257:       }
7258:       break;
7259:     case ADD_VALUES:
7260:       for (p = 0; p < numPoints; p++) {
7261:         const PetscInt     point = points[2 * p];
7262:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7263:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7264:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7265:         if (!contains) continue;
7266:         PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, NULL, values, &offset, array));
7267:       }
7268:       break;
7269:     case ADD_ALL_VALUES:
7270:       for (p = 0; p < numPoints; p++) {
7271:         const PetscInt     point = points[2 * p];
7272:         const PetscInt    *perm  = perms ? perms[p] : NULL;
7273:         const PetscScalar *flip  = flips ? flips[p] : NULL;
7274:         PetscCall(CheckPoint_Private(label, labelId, section, point, f, &offset, &contains));
7275:         if (!contains) continue;
7276:         PetscCall(updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, NULL, values, &offset, array));
7277:       }
7278:       break;
7279:     default:
7280:       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
7281:     }
7282:     PetscCall(PetscSectionRestoreFieldPointSyms(section, f, numPoints, points, &perms, &flips));
7283:   }
7284:   /* Cleanup points */
7285:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &numPoints, &points, &clSection, &clPoints, &clp));
7286:   /* Cleanup array */
7287:   PetscCall(VecRestoreArray(v, &array));
7288:   PetscFunctionReturn(PETSC_SUCCESS);
7289: }

7291: static PetscErrorCode DMPlexPrintMatSetValues(PetscViewer viewer, Mat A, PetscInt point, PetscInt numRIndices, const PetscInt rindices[], PetscInt numCIndices, const PetscInt cindices[], const PetscScalar values[])
7292: {
7293:   PetscMPIInt rank;
7294:   PetscInt    i, j;

7296:   PetscFunctionBegin;
7297:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
7298:   PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]mat for point %" PetscInt_FMT "\n", rank, point));
7299:   for (i = 0; i < numRIndices; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]mat row indices[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, i, rindices[i]));
7300:   for (i = 0; i < numCIndices; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]mat col indices[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, i, cindices[i]));
7301:   numCIndices = numCIndices ? numCIndices : numRIndices;
7302:   if (!values) PetscFunctionReturn(PETSC_SUCCESS);
7303:   for (i = 0; i < numRIndices; i++) {
7304:     PetscCall(PetscViewerASCIIPrintf(viewer, "[%d]", rank));
7305:     for (j = 0; j < numCIndices; j++) {
7306: #if defined(PETSC_USE_COMPLEX)
7307:       PetscCall(PetscViewerASCIIPrintf(viewer, " (%g,%g)", (double)PetscRealPart(values[i * numCIndices + j]), (double)PetscImaginaryPart(values[i * numCIndices + j])));
7308: #else
7309:       PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)values[i * numCIndices + j]));
7310: #endif
7311:     }
7312:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
7313:   }
7314:   PetscFunctionReturn(PETSC_SUCCESS);
7315: }

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

7320:   Input Parameters:
7321: + section - The section for this data layout
7322: . islocal - Is the section (and thus indices being requested) local or global?
7323: . point   - The point contributing dofs with these indices
7324: . off     - The global offset of this point
7325: . loff    - The local offset of each field
7326: . setBC   - The flag determining whether to include indices of boundary values
7327: . perm    - A permutation of the dofs on this point, or NULL
7328: - indperm - A permutation of the entire indices array, or NULL

7330:   Output Parameter:
7331: . indices - Indices for dofs on this point

7333:   Level: developer

7335:   Note: The indices could be local or global, depending on the value of 'off'.
7336: */
7337: PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection section, PetscBool islocal, PetscInt point, PetscInt off, PetscInt *loff, PetscBool setBC, const PetscInt perm[], const PetscInt indperm[], PetscInt indices[])
7338: {
7339:   PetscInt        dof;   /* The number of unknowns on this point */
7340:   PetscInt        cdof;  /* The number of constraints on this point */
7341:   const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
7342:   PetscInt        cind = 0, k;

7344:   PetscFunctionBegin;
7345:   PetscCheck(islocal || !setBC, PetscObjectComm((PetscObject)section), PETSC_ERR_ARG_INCOMP, "setBC incompatible with global indices; use a local section or disable setBC");
7346:   PetscCall(PetscSectionGetDof(section, point, &dof));
7347:   PetscCall(PetscSectionGetConstraintDof(section, point, &cdof));
7348:   if (!cdof || setBC) {
7349:     for (k = 0; k < dof; ++k) {
7350:       const PetscInt preind = perm ? *loff + perm[k] : *loff + k;
7351:       const PetscInt ind    = indperm ? indperm[preind] : preind;

7353:       indices[ind] = off + k;
7354:     }
7355:   } else {
7356:     PetscCall(PetscSectionGetConstraintIndices(section, point, &cdofs));
7357:     for (k = 0; k < dof; ++k) {
7358:       const PetscInt preind = perm ? *loff + perm[k] : *loff + k;
7359:       const PetscInt ind    = indperm ? indperm[preind] : preind;

7361:       if ((cind < cdof) && (k == cdofs[cind])) {
7362:         /* Insert check for returning constrained indices */
7363:         indices[ind] = -(off + k + 1);
7364:         ++cind;
7365:       } else {
7366:         indices[ind] = off + k - (islocal ? 0 : cind);
7367:       }
7368:     }
7369:   }
7370:   *loff += dof;
7371:   PetscFunctionReturn(PETSC_SUCCESS);
7372: }

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

7377:  Input Parameters:
7378: + section - a section (global or local)
7379: - islocal - `PETSC_TRUE` if requesting local indices (i.e., section is local); `PETSC_FALSE` for global
7380: . point - point within section
7381: . off - The offset of this point in the (local or global) indexed space - should match islocal and (usually) the section
7382: . foffs - array of length numFields containing the offset in canonical point ordering (the location in indices) of each field
7383: . setBC - identify constrained (boundary condition) points via involution.
7384: . perms - perms[f][permsoff][:] is a permutation of dofs within each field
7385: . permsoff - offset
7386: - indperm - index permutation

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

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

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

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

7405:  Example:
7406:  Suppose a point contains one field with three components, and for which the unconstrained indices are {10, 11, 12}.
7407:  When the middle component is constrained, we get the array {10, -12, 12} for (islocal=TRUE, setBC=FALSE).
7408:  Note that -12 is the involution of 11, so the user can involute negative indices to recover local indices.
7409:  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.

7411:  Level: developer
7412: */
7413: 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[])
7414: {
7415:   PetscInt numFields, foff, f;

7417:   PetscFunctionBegin;
7418:   PetscCheck(islocal || !setBC, PetscObjectComm((PetscObject)section), PETSC_ERR_ARG_INCOMP, "setBC incompatible with global indices; use a local section or disable setBC");
7419:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7420:   for (f = 0, foff = 0; f < numFields; ++f) {
7421:     PetscInt        fdof, cfdof;
7422:     const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
7423:     PetscInt        cind = 0, b;
7424:     const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;

7426:     PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
7427:     PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &cfdof));
7428:     if (!cfdof || setBC) {
7429:       for (b = 0; b < fdof; ++b) {
7430:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7431:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7433:         indices[ind] = off + foff + b;
7434:       }
7435:     } else {
7436:       PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
7437:       for (b = 0; b < fdof; ++b) {
7438:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7439:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7441:         if ((cind < cfdof) && (b == fcdofs[cind])) {
7442:           indices[ind] = -(off + foff + b + 1);
7443:           ++cind;
7444:         } else {
7445:           indices[ind] = off + foff + b - (islocal ? 0 : cind);
7446:         }
7447:       }
7448:     }
7449:     foff += (setBC || islocal ? fdof : (fdof - cfdof));
7450:     foffs[f] += fdof;
7451:   }
7452:   PetscFunctionReturn(PETSC_SUCCESS);
7453: }

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

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

7460:  Notes:
7461:  The semantics of this function relate to that of setBC=FALSE in DMPlexGetIndicesPointFields_Internal.
7462:  Since this function uses global indices, setBC=TRUE would be invalid, so no such argument exists.
7463: */
7464: static PetscErrorCode DMPlexGetIndicesPointFieldsSplit_Internal(PetscSection section, PetscSection globalSection, PetscInt point, PetscInt foffs[], const PetscInt ***perms, PetscInt permsoff, const PetscInt indperm[], PetscInt indices[])
7465: {
7466:   PetscInt numFields, foff, f;

7468:   PetscFunctionBegin;
7469:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7470:   for (f = 0; f < numFields; ++f) {
7471:     PetscInt        fdof, cfdof;
7472:     const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
7473:     PetscInt        cind = 0, b;
7474:     const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;

7476:     PetscCall(PetscSectionGetFieldDof(section, point, f, &fdof));
7477:     PetscCall(PetscSectionGetFieldConstraintDof(section, point, f, &cfdof));
7478:     PetscCall(PetscSectionGetFieldOffset(globalSection, point, f, &foff));
7479:     if (!cfdof) {
7480:       for (b = 0; b < fdof; ++b) {
7481:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7482:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7484:         indices[ind] = foff + b;
7485:       }
7486:     } else {
7487:       PetscCall(PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs));
7488:       for (b = 0; b < fdof; ++b) {
7489:         const PetscInt preind = perm ? foffs[f] + perm[b] : foffs[f] + b;
7490:         const PetscInt ind    = indperm ? indperm[preind] : preind;

7492:         if ((cind < cfdof) && (b == fcdofs[cind])) {
7493:           indices[ind] = -(foff + b + 1);
7494:           ++cind;
7495:         } else {
7496:           indices[ind] = foff + b - cind;
7497:         }
7498:       }
7499:     }
7500:     foffs[f] += fdof;
7501:   }
7502:   PetscFunctionReturn(PETSC_SUCCESS);
7503: }

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

7509:   PetscFunctionBegin;
7510:   PetscCall(PetscSectionGetNumFields(section, &numFields));
7511:   PetscCall(PetscSectionGetChart(section, &sStart, &sEnd));
7512:   PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
7513:   for (PetscInt p = 0; p < nPoints; p++) {
7514:     PetscInt     b       = pnts[2 * p];
7515:     PetscInt     bSecDof = 0, bOff;
7516:     PetscInt     cSecDof = 0;
7517:     PetscSection indices_section;

7519:     if (b >= sStart && b < sEnd) PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7520:     if (!bSecDof) continue;
7521:     if (b >= cStart && b < cEnd) PetscCall(PetscSectionGetDof(cSec, b, &cSecDof));
7522:     indices_section = cSecDof > 0 ? cSec : section;
7523:     if (numFields) {
7524:       PetscInt fStart[32], fEnd[32];

7526:       fStart[0] = 0;
7527:       fEnd[0]   = 0;
7528:       for (PetscInt f = 0; f < numFields; f++) {
7529:         PetscInt fDof = 0;

7531:         PetscCall(PetscSectionGetFieldDof(indices_section, b, f, &fDof));
7532:         fStart[f + 1] = fStart[f] + fDof;
7533:         fEnd[f + 1]   = fStart[f + 1];
7534:       }
7535:       PetscCall(PetscSectionGetOffset(indices_section, b, &bOff));
7536:       // only apply permutations on one side
7537:       PetscCall(DMPlexGetIndicesPointFields_Internal(indices_section, PETSC_TRUE, b, bOff, fEnd, PETSC_TRUE, perms, perms ? p : -1, NULL, tmpIndices));
7538:       for (PetscInt f = 0; f < numFields; f++) {
7539:         for (PetscInt i = fStart[f]; i < fEnd[f]; i++) { indices[fieldOffsets[f]++] = (cSecDof > 0) ? tmpIndices[i] : -(tmpIndices[i] + 1); }
7540:       }
7541:     } else {
7542:       PetscInt bEnd = 0;

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

7547:       for (PetscInt i = 0; i < bEnd; i++) indices[fieldOffsets[0]++] = (cSecDof > 0) ? tmpIndices[i] : -(tmpIndices[i] + 1);
7548:     }
7549:   }
7550:   PetscFunctionReturn(PETSC_SUCCESS);
7551: }

7553: 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[])
7554: {
7555:   Mat             cMat;
7556:   PetscSection    aSec, cSec;
7557:   IS              aIS;
7558:   PetscInt        aStart = -1, aEnd = -1;
7559:   PetscInt        sStart = -1, sEnd = -1;
7560:   PetscInt        cStart = -1, cEnd = -1;
7561:   const PetscInt *anchors;
7562:   PetscInt        numFields, p;
7563:   PetscInt        newNumPoints = 0, newNumIndices = 0;
7564:   PetscInt       *newPoints, *indices, *newIndices, *tmpIndices, *tmpNewIndices;
7565:   PetscInt        oldOffsets[32];
7566:   PetscInt        newOffsets[32];
7567:   PetscInt        oldOffsetsCopy[32];
7568:   PetscInt        newOffsetsCopy[32];
7569:   PetscScalar    *modMat         = NULL;
7570:   PetscBool       anyConstrained = PETSC_FALSE;

7572:   PetscFunctionBegin;
7575:   PetscCall(PetscSectionGetNumFields(section, &numFields));

7577:   PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
7578:   /* if there are point-to-point constraints */
7579:   if (aSec) {
7580:     PetscCall(PetscArrayzero(newOffsets, 32));
7581:     PetscCall(PetscArrayzero(oldOffsets, 32));
7582:     PetscCall(ISGetIndices(aIS, &anchors));
7583:     PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
7584:     PetscCall(PetscSectionGetChart(section, &sStart, &sEnd));
7585:     /* figure out how many points are going to be in the new element matrix
7586:      * (we allow double counting, because it's all just going to be summed
7587:      * into the global matrix anyway) */
7588:     for (p = 0; p < 2 * numPoints; p += 2) {
7589:       PetscInt b    = points[p];
7590:       PetscInt bDof = 0, bSecDof = 0;

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

7595:       for (PetscInt f = 0; f < numFields; f++) {
7596:         PetscInt fDof = 0;

7598:         PetscCall(PetscSectionGetFieldDof(section, b, f, &fDof));
7599:         oldOffsets[f + 1] += fDof;
7600:       }
7601:       if (b >= aStart && b < aEnd) PetscCall(PetscSectionGetDof(aSec, b, &bDof));
7602:       if (bDof) {
7603:         /* this point is constrained */
7604:         /* it is going to be replaced by its anchors */
7605:         PetscInt bOff, q;

7607:         PetscCall(PetscSectionGetOffset(aSec, b, &bOff));
7608:         for (q = 0; q < bDof; q++) {
7609:           PetscInt a    = anchors[bOff + q];
7610:           PetscInt aDof = 0;

7612:           if (a >= sStart && a < sEnd) PetscCall(PetscSectionGetDof(section, a, &aDof));
7613:           if (aDof) {
7614:             anyConstrained = PETSC_TRUE;
7615:             newNumPoints += 1;
7616:           }
7617:           newNumIndices += aDof;
7618:           for (PetscInt f = 0; f < numFields; ++f) {
7619:             PetscInt fDof = 0;

7621:             if (a >= sStart && a < sEnd) PetscCall(PetscSectionGetFieldDof(section, a, f, &fDof));
7622:             newOffsets[f + 1] += fDof;
7623:           }
7624:         }
7625:       } else {
7626:         /* this point is not constrained */
7627:         newNumPoints++;
7628:         newNumIndices += bSecDof;
7629:         for (PetscInt f = 0; f < numFields; ++f) {
7630:           PetscInt fDof;

7632:           PetscCall(PetscSectionGetFieldDof(section, b, f, &fDof));
7633:           newOffsets[f + 1] += fDof;
7634:         }
7635:       }
7636:     }
7637:   }
7638:   if (!anyConstrained) {
7639:     if (outNumPoints) *outNumPoints = 0;
7640:     if (outNumIndices) *outNumIndices = 0;
7641:     if (outPoints) *outPoints = NULL;
7642:     if (outMat) *outMat = NULL;
7643:     if (aSec) PetscCall(ISRestoreIndices(aIS, &anchors));
7644:     PetscFunctionReturn(PETSC_SUCCESS);
7645:   }

7647:   if (outNumPoints) *outNumPoints = newNumPoints;
7648:   if (outNumIndices) *outNumIndices = newNumIndices;

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

7653:   if (!outPoints && !outMat) {
7654:     if (offsets) {
7655:       for (PetscInt f = 0; f <= numFields; f++) offsets[f] = newOffsets[f];
7656:     }
7657:     if (aSec) PetscCall(ISRestoreIndices(aIS, &anchors));
7658:     PetscFunctionReturn(PETSC_SUCCESS);
7659:   }

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

7664:   PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
7665:   PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));

7667:   /* output arrays */
7668:   PetscCall(DMGetWorkArray(dm, 2 * newNumPoints, MPIU_INT, &newPoints));
7669:   PetscCall(PetscArrayzero(newPoints, 2 * newNumPoints));

7671:   // get the new Points
7672:   for (PetscInt p = 0, newP = 0; p < numPoints; p++) {
7673:     PetscInt b    = points[2 * p];
7674:     PetscInt bDof = 0, bSecDof = 0, bOff;

7676:     if (b >= sStart && b < sEnd) PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7677:     if (!bSecDof) continue;
7678:     if (b >= aStart && b < aEnd) PetscCall(PetscSectionGetDof(aSec, b, &bDof));
7679:     if (bDof) {
7680:       PetscCall(PetscSectionGetOffset(aSec, b, &bOff));
7681:       for (PetscInt q = 0; q < bDof; q++) {
7682:         PetscInt a = anchors[bOff + q], aDof = 0;

7684:         if (a >= sStart && a < sEnd) PetscCall(PetscSectionGetDof(section, a, &aDof));
7685:         if (aDof) {
7686:           newPoints[2 * newP]     = a;
7687:           newPoints[2 * newP + 1] = 0; // orientations are accounted for in constructing the matrix, newly added points are in default orientation
7688:           newP++;
7689:         }
7690:       }
7691:     } else {
7692:       newPoints[2 * newP]     = b;
7693:       newPoints[2 * newP + 1] = points[2 * p + 1];
7694:       newP++;
7695:     }
7696:   }

7698:   if (outMat) {
7699:     PetscScalar *tmpMat;
7700:     PetscCall(PetscArraycpy(oldOffsetsCopy, oldOffsets, 32));
7701:     PetscCall(PetscArraycpy(newOffsetsCopy, newOffsets, 32));

7703:     PetscCall(DMGetWorkArray(dm, numIndices, MPIU_INT, &indices));
7704:     PetscCall(DMGetWorkArray(dm, numIndices, MPIU_INT, &tmpIndices));
7705:     PetscCall(DMGetWorkArray(dm, newNumIndices, MPIU_INT, &newIndices));
7706:     PetscCall(DMGetWorkArray(dm, newNumIndices, MPIU_INT, &tmpNewIndices));

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

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

7714:     PetscCall(DMGetWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &modMat));
7715:     PetscCall(DMGetWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &tmpMat));
7716:     PetscCall(PetscArrayzero(modMat, newNumIndices * numIndices));
7717:     // for each field, insert the anchor modification into modMat
7718:     for (PetscInt f = 0; f < PetscMax(1, numFields); f++) {
7719:       PetscInt fStart    = oldOffsets[f];
7720:       PetscInt fNewStart = newOffsets[f];
7721:       for (PetscInt p = 0, newP = 0, o = fStart, oNew = fNewStart; p < numPoints; p++) {
7722:         PetscInt b    = points[2 * p];
7723:         PetscInt bDof = 0, bSecDof = 0, bOff;

7725:         if (b >= sStart && b < sEnd) {
7726:           if (numFields) {
7727:             PetscCall(PetscSectionGetFieldDof(section, b, f, &bSecDof));
7728:           } else {
7729:             PetscCall(PetscSectionGetDof(section, b, &bSecDof));
7730:           }
7731:         }
7732:         if (!bSecDof) continue;
7733:         if (b >= aStart && b < aEnd) PetscCall(PetscSectionGetDof(aSec, b, &bDof));
7734:         if (bDof) {
7735:           PetscCall(PetscSectionGetOffset(aSec, b, &bOff));
7736:           for (PetscInt q = 0; q < bDof; q++, newP++) {
7737:             PetscInt a = anchors[bOff + q], aDof = 0;

7739:             if (a >= sStart && a < sEnd) {
7740:               if (numFields) {
7741:                 PetscCall(PetscSectionGetFieldDof(section, a, f, &aDof));
7742:               } else {
7743:                 PetscCall(PetscSectionGetDof(section, a, &aDof));
7744:               }
7745:             }
7746:             if (aDof) {
7747:               PetscCall(MatGetValues(cMat, bSecDof, &indices[o], aDof, &newIndices[oNew], tmpMat));
7748:               for (PetscInt d = 0; d < bSecDof; d++) {
7749:                 for (PetscInt e = 0; e < aDof; e++) modMat[(o + d) * newNumIndices + oNew + e] = tmpMat[d * aDof + e];
7750:               }
7751:             }
7752:             oNew += aDof;
7753:           }
7754:         } else {
7755:           // Insert the identity matrix in this block
7756:           for (PetscInt d = 0; d < bSecDof; d++) modMat[(o + d) * newNumIndices + oNew + d] = 1;
7757:           oNew += bSecDof;
7758:           newP++;
7759:         }
7760:         o += bSecDof;
7761:       }
7762:     }

7764:     *outMat = modMat;

7766:     PetscCall(DMRestoreWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &tmpMat));
7767:     PetscCall(DMRestoreWorkArray(dm, newNumIndices, MPIU_INT, &tmpNewIndices));
7768:     PetscCall(DMRestoreWorkArray(dm, newNumIndices, MPIU_INT, &newIndices));
7769:     PetscCall(DMRestoreWorkArray(dm, numIndices, MPIU_INT, &tmpIndices));
7770:     PetscCall(DMRestoreWorkArray(dm, numIndices, MPIU_INT, &indices));
7771:   }
7772:   PetscCall(ISRestoreIndices(aIS, &anchors));

7774:   /* output */
7775:   if (outPoints) {
7776:     *outPoints = newPoints;
7777:   } else {
7778:     PetscCall(DMRestoreWorkArray(dm, 2 * newNumPoints, MPIU_INT, &newPoints));
7779:   }
7780:   for (PetscInt f = 0; f <= numFields; f++) offsets[f] = newOffsets[f];
7781:   PetscFunctionReturn(PETSC_SUCCESS);
7782: }

7784: PETSC_INTERN PetscErrorCode DMPlexAnchorsModifyMat_Internal(DM dm, PetscSection section, PetscInt numPoints, PetscInt numIndices, const PetscInt points[], const PetscInt ***perms, PetscInt numRows, PetscInt numCols, const PetscScalar values[], PetscInt *outNumPoints, PetscInt *outNumIndices, PetscInt *outPoints[], PetscScalar *outValues[], PetscInt offsets[], PetscBool multiplyRight, PetscBool multiplyLeft)
7785: {
7786:   PetscScalar *modMat        = NULL;
7787:   PetscInt     newNumIndices = -1;

7789:   PetscFunctionBegin;
7790:   /* If M is the matrix represented by values, get the matrix C such that we will add M * C (or, if multiplyLeft, C^T * M * C) into the global matrix.
7791:      modMat is that matrix C */
7792:   PetscCall(DMPlexAnchorsGetSubMatModification(dm, section, numPoints, numIndices, points, perms, outNumPoints, &newNumIndices, outPoints, offsets, outValues ? &modMat : NULL));
7793:   if (outNumIndices) *outNumIndices = newNumIndices;
7794:   if (modMat) {
7795:     const PetscScalar *newValues = values;

7797:     if (multiplyRight) {
7798:       PetscScalar *newNewValues = NULL;
7799:       PetscBLASInt M, N, K;
7800:       PetscScalar  a = 1.0, b = 0.0;

7802:       PetscCheck(numCols == numIndices, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "values matrix has the wrong number of columns: %" PetscInt_FMT ", expected %" PetscInt_FMT, numCols, numIndices);

7804:       PetscCall(PetscBLASIntCast(newNumIndices, &M));
7805:       PetscCall(PetscBLASIntCast(numRows, &N));
7806:       PetscCall(PetscBLASIntCast(numIndices, &K));
7807:       PetscCall(DMGetWorkArray(dm, numRows * newNumIndices, MPIU_SCALAR, &newNewValues));
7808:       // row-major to column-major conversion, right multiplication becomes left multiplication
7809:       PetscCallBLAS("BLASgemm", BLASgemm_("N", "N", &M, &N, &K, &a, modMat, &M, newValues, &K, &b, newNewValues, &M));
7810:       numCols   = newNumIndices;
7811:       newValues = newNewValues;
7812:     }

7814:     if (multiplyLeft) {
7815:       PetscScalar *newNewValues = NULL;
7816:       PetscBLASInt M, N, K;
7817:       PetscScalar  a = 1.0, b = 0.0;

7819:       PetscCheck(numRows == numIndices, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "values matrix has the wrong number of rows: %" PetscInt_FMT ", expected %" PetscInt_FMT, numRows, numIndices);

7821:       PetscCall(PetscBLASIntCast(numCols, &M));
7822:       PetscCall(PetscBLASIntCast(newNumIndices, &N));
7823:       PetscCall(PetscBLASIntCast(numIndices, &K));
7824:       PetscCall(DMGetWorkArray(dm, newNumIndices * numCols, MPIU_SCALAR, &newNewValues));
7825:       // row-major to column-major conversion, left multiplication becomes right multiplication
7826:       PetscCallBLAS("BLASgemm", BLASgemm_("N", "T", &M, &N, &K, &a, newValues, &M, modMat, &N, &b, newNewValues, &M));
7827:       if (newValues != values) PetscCall(DMRestoreWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &newValues));
7828:       newValues = newNewValues;
7829:     }
7830:     *outValues = (PetscScalar *)newValues;
7831:     PetscCall(DMRestoreWorkArray(dm, numIndices * newNumIndices, MPIU_SCALAR, &modMat));
7832:   }
7833:   PetscFunctionReturn(PETSC_SUCCESS);
7834: }

7836: PETSC_INTERN PetscErrorCode DMPlexAnchorsModifyMat(DM dm, PetscSection section, PetscInt numPoints, PetscInt numIndices, const PetscInt points[], const PetscInt ***perms, const PetscScalar values[], PetscInt *outNumPoints, PetscInt *outNumIndices, PetscInt *outPoints[], PetscScalar *outValues[], PetscInt offsets[], PetscBool multiplyLeft)
7837: {
7838:   PetscFunctionBegin;
7839:   PetscCall(DMPlexAnchorsModifyMat_Internal(dm, section, numPoints, numIndices, points, perms, numIndices, numIndices, values, outNumPoints, outNumIndices, outPoints, outValues, offsets, PETSC_TRUE, multiplyLeft));
7840:   PetscFunctionReturn(PETSC_SUCCESS);
7841: }

7843: static PetscErrorCode DMPlexGetClosureIndicesSize_Internal(DM dm, PetscSection section, PetscInt point, PetscInt *closureSize)
7844: {
7845:   /* Closure ordering */
7846:   PetscSection    clSection;
7847:   IS              clPoints;
7848:   const PetscInt *clp;
7849:   PetscInt       *points;
7850:   PetscInt        Ncl, Ni = 0;

7852:   PetscFunctionBeginHot;
7853:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &Ncl, &points, &clSection, &clPoints, &clp));
7854:   for (PetscInt p = 0; p < Ncl * 2; p += 2) {
7855:     PetscInt dof;

7857:     PetscCall(PetscSectionGetDof(section, points[p], &dof));
7858:     Ni += dof;
7859:   }
7860:   PetscCall(DMPlexRestoreCompressedClosure(dm, section, point, &Ncl, &points, &clSection, &clPoints, &clp));
7861:   *closureSize = Ni;
7862:   PetscFunctionReturn(PETSC_SUCCESS);
7863: }

7865: static PetscErrorCode DMPlexGetClosureIndices_Internal(DM dm, PetscSection section, PetscSection idxSection, PetscInt point, PetscBool useClPerm, PetscInt *numRows, PetscInt *numCols, PetscInt *indices[], PetscInt outOffsets[], PetscScalar *values[], PetscBool multiplyRight, PetscBool multiplyLeft)
7866: {
7867:   /* Closure ordering */
7868:   PetscSection    clSection;
7869:   IS              clPoints;
7870:   const PetscInt *clp;
7871:   PetscInt       *points;
7872:   const PetscInt *clperm = NULL;
7873:   /* Dof permutation and sign flips */
7874:   const PetscInt    **perms[32] = {NULL};
7875:   const PetscScalar **flips[32] = {NULL};
7876:   PetscScalar        *valCopy   = NULL;
7877:   /* Hanging node constraints */
7878:   PetscInt    *pointsC = NULL;
7879:   PetscScalar *valuesC = NULL;
7880:   PetscInt     NclC, NiC;

7882:   PetscInt *idx;
7883:   PetscInt  Nf, Ncl, Ni = 0, offsets[32], p, f;
7884:   PetscBool isLocal = (section == idxSection) ? PETSC_TRUE : PETSC_FALSE;
7885:   PetscInt  idxStart, idxEnd;
7886:   PetscInt  nRows, nCols;

7888:   PetscFunctionBeginHot;
7892:   PetscAssertPointer(numRows, 6);
7893:   PetscAssertPointer(numCols, 7);
7894:   if (indices) PetscAssertPointer(indices, 8);
7895:   if (outOffsets) PetscAssertPointer(outOffsets, 9);
7896:   if (values) PetscAssertPointer(values, 10);
7897:   PetscCall(PetscSectionGetNumFields(section, &Nf));
7898:   PetscCheck(Nf <= 31, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %" PetscInt_FMT " limited to 31", Nf);
7899:   PetscCall(PetscArrayzero(offsets, 32));
7900:   /* 1) Get points in closure */
7901:   PetscCall(DMPlexGetCompressedClosure(dm, section, point, 0, &Ncl, &points, &clSection, &clPoints, &clp));
7902:   if (useClPerm) {
7903:     PetscInt depth, clsize;
7904:     PetscCall(DMPlexGetPointDepth(dm, point, &depth));
7905:     for (clsize = 0, p = 0; p < Ncl; p++) {
7906:       PetscInt dof;
7907:       PetscCall(PetscSectionGetDof(section, points[2 * p], &dof));
7908:       clsize += dof;
7909:     }
7910:     PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject)dm, depth, clsize, &clperm));
7911:   }
7912:   /* 2) Get number of indices on these points and field offsets from section */
7913: