xref: /petsc/src/dm/impls/plex/plexglvis.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 #include <petsc/private/glvisviewerimpl.h>
2 #include <petsc/private/petscimpl.h>
3 #include <petsc/private/dmpleximpl.h>
4 #include <petscbt.h>
5 #include <petscdmplex.h>
6 #include <petscsf.h>
7 #include <petscds.h>
8 
9 typedef struct {
10   PetscInt    nf;
11   VecScatter *scctx;
12 } GLVisViewerCtx;
13 
14 static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) {
15   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
16   PetscInt        i;
17 
18   PetscFunctionBegin;
19   for (i = 0; i < ctx->nf; i++) PetscCall(VecScatterDestroy(&ctx->scctx[i]));
20   PetscCall(PetscFree(ctx->scctx));
21   PetscCall(PetscFree(vctx));
22   PetscFunctionReturn(0);
23 }
24 
25 static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) {
26   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
27   PetscInt        f;
28 
29   PetscFunctionBegin;
30   for (f = 0; f < nf; f++) {
31     PetscCall(VecScatterBegin(ctx->scctx[f], (Vec)oX, (Vec)oXfield[f], INSERT_VALUES, SCATTER_FORWARD));
32     PetscCall(VecScatterEnd(ctx->scctx[f], (Vec)oX, (Vec)oXfield[f], INSERT_VALUES, SCATTER_FORWARD));
33   }
34   PetscFunctionReturn(0);
35 }
36 
37 /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */
38 PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer) {
39   DM              dm = (DM)odm;
40   Vec             xlocal, xfield, *Ufield;
41   PetscDS         ds;
42   IS              globalNum, isfield;
43   PetscBT         vown;
44   char          **fieldname = NULL, **fec_type = NULL;
45   const PetscInt *gNum;
46   PetscInt       *nlocal, *bs, *idxs, *dims;
47   PetscInt        f, maxfields, nfields, c, totc, totdofs, Nv, cum, i;
48   PetscInt        dim, cStart, cEnd, vStart, vEnd;
49   GLVisViewerCtx *ctx;
50   PetscSection    s;
51 
52   PetscFunctionBegin;
53   PetscCall(DMGetDimension(dm, &dim));
54   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
55   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
56   PetscCall(PetscObjectQuery((PetscObject)dm, "_glvis_plex_gnum", (PetscObject *)&globalNum));
57   if (!globalNum) {
58     PetscCall(DMPlexCreateCellNumbering_Internal(dm, PETSC_TRUE, &globalNum));
59     PetscCall(PetscObjectCompose((PetscObject)dm, "_glvis_plex_gnum", (PetscObject)globalNum));
60     PetscCall(PetscObjectDereference((PetscObject)globalNum));
61   }
62   PetscCall(ISGetIndices(globalNum, &gNum));
63   PetscCall(PetscBTCreate(vEnd - vStart, &vown));
64   for (c = cStart, totc = 0; c < cEnd; c++) {
65     if (gNum[c - cStart] >= 0) {
66       PetscInt i, numPoints, *points = NULL;
67 
68       totc++;
69       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &numPoints, &points));
70       for (i = 0; i < numPoints * 2; i += 2) {
71         if ((points[i] >= vStart) && (points[i] < vEnd)) PetscCall(PetscBTSet(vown, points[i] - vStart));
72       }
73       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &numPoints, &points));
74     }
75   }
76   for (f = 0, Nv = 0; f < vEnd - vStart; f++)
77     if (PetscLikely(PetscBTLookup(vown, f))) Nv++;
78 
79   PetscCall(DMCreateLocalVector(dm, &xlocal));
80   PetscCall(VecGetLocalSize(xlocal, &totdofs));
81   PetscCall(DMGetLocalSection(dm, &s));
82   PetscCall(PetscSectionGetNumFields(s, &nfields));
83   for (f = 0, maxfields = 0; f < nfields; f++) {
84     PetscInt bs;
85 
86     PetscCall(PetscSectionGetFieldComponents(s, f, &bs));
87     maxfields += bs;
88   }
89   PetscCall(PetscCalloc7(maxfields, &fieldname, maxfields, &nlocal, maxfields, &bs, maxfields, &dims, maxfields, &fec_type, totdofs, &idxs, maxfields, &Ufield));
90   PetscCall(PetscNew(&ctx));
91   PetscCall(PetscCalloc1(maxfields, &ctx->scctx));
92   PetscCall(DMGetDS(dm, &ds));
93   if (ds) {
94     for (f = 0; f < nfields; f++) {
95       const char *fname;
96       char        name[256];
97       PetscObject disc;
98       size_t      len;
99 
100       PetscCall(PetscSectionGetFieldName(s, f, &fname));
101       PetscCall(PetscStrlen(fname, &len));
102       if (len) {
103         PetscCall(PetscStrcpy(name, fname));
104       } else {
105         PetscCall(PetscSNPrintf(name, 256, "Field%" PetscInt_FMT, f));
106       }
107       PetscCall(PetscDSGetDiscretization(ds, f, &disc));
108       if (disc) {
109         PetscClassId id;
110         PetscInt     Nc;
111         char         fec[64];
112 
113         PetscCall(PetscObjectGetClassId(disc, &id));
114         if (id == PETSCFE_CLASSID) {
115           PetscFE            fem = (PetscFE)disc;
116           PetscDualSpace     sp;
117           PetscDualSpaceType spname;
118           PetscInt           order;
119           PetscBool          islag, continuous, H1 = PETSC_TRUE;
120 
121           PetscCall(PetscFEGetNumComponents(fem, &Nc));
122           PetscCall(PetscFEGetDualSpace(fem, &sp));
123           PetscCall(PetscDualSpaceGetType(sp, &spname));
124           PetscCall(PetscStrcmp(spname, PETSCDUALSPACELAGRANGE, &islag));
125           PetscCheck(islag, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dual space");
126           PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &continuous));
127           PetscCall(PetscDualSpaceGetOrder(sp, &order));
128           if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */
129             PetscCall(PetscSNPrintf(fec, 64, "FiniteElementCollection: H1_%" PetscInt_FMT "D_P1", dim));
130           } else {
131             PetscCheck(continuous || !order, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Discontinuous space visualization currently unsupported for order %" PetscInt_FMT, order);
132             H1 = PETSC_FALSE;
133             PetscCall(PetscSNPrintf(fec, 64, "FiniteElementCollection: L2_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, order));
134           }
135           PetscCall(PetscStrallocpy(name, &fieldname[ctx->nf]));
136           bs[ctx->nf]   = Nc;
137           dims[ctx->nf] = dim;
138           if (H1) {
139             nlocal[ctx->nf] = Nc * Nv;
140             PetscCall(PetscStrallocpy(fec, &fec_type[ctx->nf]));
141             PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nv * Nc, &xfield));
142             for (i = 0, cum = 0; i < vEnd - vStart; i++) {
143               PetscInt j, off;
144 
145               if (PetscUnlikely(!PetscBTLookup(vown, i))) continue;
146               PetscCall(PetscSectionGetFieldOffset(s, i + vStart, f, &off));
147               for (j = 0; j < Nc; j++) idxs[cum++] = off + j;
148             }
149             PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal), Nv * Nc, idxs, PETSC_USE_POINTER, &isfield));
150           } else {
151             nlocal[ctx->nf] = Nc * totc;
152             PetscCall(PetscStrallocpy(fec, &fec_type[ctx->nf]));
153             PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nc * totc, &xfield));
154             for (i = 0, cum = 0; i < cEnd - cStart; i++) {
155               PetscInt j, off;
156 
157               if (PetscUnlikely(gNum[i] < 0)) continue;
158               PetscCall(PetscSectionGetFieldOffset(s, i + cStart, f, &off));
159               for (j = 0; j < Nc; j++) idxs[cum++] = off + j;
160             }
161             PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal), totc * Nc, idxs, PETSC_USE_POINTER, &isfield));
162           }
163           PetscCall(VecScatterCreate(xlocal, isfield, xfield, NULL, &ctx->scctx[ctx->nf]));
164           PetscCall(VecDestroy(&xfield));
165           PetscCall(ISDestroy(&isfield));
166           ctx->nf++;
167         } else if (id == PETSCFV_CLASSID) {
168           PetscInt c;
169 
170           PetscCall(PetscFVGetNumComponents((PetscFV)disc, &Nc));
171           PetscCall(PetscSNPrintf(fec, 64, "FiniteElementCollection: L2_%" PetscInt_FMT "D_P0", dim));
172           for (c = 0; c < Nc; c++) {
173             char comp[256];
174             PetscCall(PetscSNPrintf(comp, 256, "%s-Comp%" PetscInt_FMT, name, c));
175             PetscCall(PetscStrallocpy(comp, &fieldname[ctx->nf]));
176             bs[ctx->nf]     = 1; /* Does PetscFV support components with different block size? */
177             nlocal[ctx->nf] = totc;
178             dims[ctx->nf]   = dim;
179             PetscCall(PetscStrallocpy(fec, &fec_type[ctx->nf]));
180             PetscCall(VecCreateSeq(PETSC_COMM_SELF, totc, &xfield));
181             for (i = 0, cum = 0; i < cEnd - cStart; i++) {
182               PetscInt off;
183 
184               if (PetscUnlikely(gNum[i]) < 0) continue;
185               PetscCall(PetscSectionGetFieldOffset(s, i + cStart, f, &off));
186               idxs[cum++] = off + c;
187             }
188             PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal), totc, idxs, PETSC_USE_POINTER, &isfield));
189             PetscCall(VecScatterCreate(xlocal, isfield, xfield, NULL, &ctx->scctx[ctx->nf]));
190             PetscCall(VecDestroy(&xfield));
191             PetscCall(ISDestroy(&isfield));
192             ctx->nf++;
193           }
194         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
195       } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Missing discretization for field %" PetscInt_FMT, f);
196     }
197   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Needs a DS attached to the DM");
198   PetscCall(PetscBTDestroy(&vown));
199   PetscCall(VecDestroy(&xlocal));
200   PetscCall(ISRestoreIndices(globalNum, &gNum));
201 
202   /* create work vectors */
203   for (f = 0; f < ctx->nf; f++) {
204     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm), nlocal[f], PETSC_DECIDE, &Ufield[f]));
205     PetscCall(PetscObjectSetName((PetscObject)Ufield[f], fieldname[f]));
206     PetscCall(VecSetBlockSize(Ufield[f], bs[f]));
207     PetscCall(VecSetDM(Ufield[f], dm));
208   }
209 
210   /* customize the viewer */
211   PetscCall(PetscViewerGLVisSetFields(viewer, ctx->nf, (const char **)fec_type, dims, DMPlexSampleGLVisFields_Private, (PetscObject *)Ufield, ctx, DestroyGLVisViewerCtx_Private));
212   for (f = 0; f < ctx->nf; f++) {
213     PetscCall(PetscFree(fieldname[f]));
214     PetscCall(PetscFree(fec_type[f]));
215     PetscCall(VecDestroy(&Ufield[f]));
216   }
217   PetscCall(PetscFree7(fieldname, nlocal, bs, dims, fec_type, idxs, Ufield));
218   PetscFunctionReturn(0);
219 }
220 
221 typedef enum {
222   MFEM_POINT = 0,
223   MFEM_SEGMENT,
224   MFEM_TRIANGLE,
225   MFEM_SQUARE,
226   MFEM_TETRAHEDRON,
227   MFEM_CUBE,
228   MFEM_PRISM,
229   MFEM_UNDEF
230 } MFEM_cid;
231 
232 MFEM_cid mfem_table_cid[4][7] = {
233   {MFEM_POINT, MFEM_UNDEF, MFEM_UNDEF,   MFEM_UNDEF,    MFEM_UNDEF,       MFEM_UNDEF, MFEM_UNDEF},
234   {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF,    MFEM_UNDEF,       MFEM_UNDEF, MFEM_UNDEF},
235   {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_TRIANGLE, MFEM_SQUARE,      MFEM_UNDEF, MFEM_UNDEF},
236   {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF,    MFEM_TETRAHEDRON, MFEM_PRISM, MFEM_CUBE }
237 };
238 
239 MFEM_cid mfem_table_cid_unint[4][9] = {
240   {MFEM_POINT, MFEM_UNDEF, MFEM_UNDEF,   MFEM_UNDEF,    MFEM_UNDEF,       MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_UNDEF},
241   {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF,    MFEM_UNDEF,       MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_UNDEF},
242   {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_TRIANGLE, MFEM_SQUARE,      MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_UNDEF},
243   {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF,    MFEM_TETRAHEDRON, MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_CUBE }
244 };
245 
246 static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid) {
247   DMLabel  dlabel;
248   PetscInt depth, csize, pdepth, dim;
249 
250   PetscFunctionBegin;
251   PetscCall(DMPlexGetDepthLabel(dm, &dlabel));
252   PetscCall(DMLabelGetValue(dlabel, p, &pdepth));
253   PetscCall(DMPlexGetConeSize(dm, p, &csize));
254   PetscCall(DMPlexGetDepth(dm, &depth));
255   PetscCall(DMGetDimension(dm, &dim));
256   if (label) {
257     PetscCall(DMLabelGetValue(label, p, mid));
258     *mid = *mid - minl + 1; /* MFEM does not like negative markers */
259   } else *mid = 1;
260   if (depth >= 0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
261     PetscCheck(dim >= 0 && dim <= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT, dim);
262     PetscCheck(csize <= 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "Found cone size %" PetscInt_FMT " for point %" PetscInt_FMT, csize, p);
263     PetscCheck(depth == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Found depth %" PetscInt_FMT " for point %" PetscInt_FMT ". You should interpolate the mesh first", depth, p);
264     *cid = mfem_table_cid_unint[dim][csize];
265   } else {
266     PetscCheck(csize <= 6, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cone size %" PetscInt_FMT " for point %" PetscInt_FMT, csize, p);
267     PetscCheck(pdepth >= 0 && pdepth <= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Depth %" PetscInt_FMT " for point %" PetscInt_FMT, csize, p);
268     *cid = mfem_table_cid[pdepth][csize];
269   }
270   PetscFunctionReturn(0);
271 }
272 
273 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[]) {
274   PetscInt dim, sdim, dof = 0, off = 0, i, q, vStart, vEnd, numPoints, *points = NULL;
275 
276   PetscFunctionBegin;
277   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
278   PetscCall(DMGetDimension(dm, &dim));
279   sdim = dim;
280   if (csec) {
281     PetscInt sStart, sEnd;
282 
283     PetscCall(DMGetCoordinateDim(dm, &sdim));
284     PetscCall(PetscSectionGetChart(csec, &sStart, &sEnd));
285     PetscCall(PetscSectionGetOffset(csec, vStart, &off));
286     off = off / sdim;
287     if (p >= sStart && p < sEnd) PetscCall(PetscSectionGetDof(csec, p, &dof));
288   }
289   if (!dof) {
290     PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &numPoints, &points));
291     for (i = 0, q = 0; i < numPoints * 2; i += 2)
292       if ((points[i] >= vStart) && (points[i] < vEnd)) vids[q++] = points[i] - vStart + off;
293     PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &numPoints, &points));
294   } else {
295     PetscCall(PetscSectionGetOffset(csec, p, &off));
296     PetscCall(PetscSectionGetDof(csec, p, &dof));
297     for (q = 0; q < dof / sdim; q++) vids[q] = off / sdim + q;
298   }
299   *nv = q;
300   PetscFunctionReturn(0);
301 }
302 
303 static PetscErrorCode GLVisCreateFE(PetscFE femIn, char name[32], PetscFE *fem, IS *perm) {
304   DM              K;
305   PetscSpace      P;
306   PetscDualSpace  Q;
307   PetscQuadrature q, fq;
308   PetscInt        dim, deg, dof;
309   DMPolytopeType  ptype;
310   PetscBool       isSimplex, isTensor;
311   PetscBool       continuity = PETSC_FALSE;
312   PetscDTNodeType nodeType   = PETSCDTNODES_GAUSSJACOBI;
313   PetscBool       endpoint   = PETSC_TRUE;
314   MPI_Comm        comm;
315 
316   PetscFunctionBegin;
317   PetscCall(PetscObjectGetComm((PetscObject)femIn, &comm));
318   PetscCall(PetscFEGetBasisSpace(femIn, &P));
319   PetscCall(PetscFEGetDualSpace(femIn, &Q));
320   PetscCall(PetscDualSpaceGetDM(Q, &K));
321   PetscCall(DMGetDimension(K, &dim));
322   PetscCall(PetscSpaceGetDegree(P, &deg, NULL));
323   PetscCall(PetscSpaceGetNumComponents(P, &dof));
324   PetscCall(DMPlexGetCellType(K, 0, &ptype));
325   switch (ptype) {
326   case DM_POLYTOPE_QUADRILATERAL:
327   case DM_POLYTOPE_HEXAHEDRON: isSimplex = PETSC_FALSE; break;
328   default: isSimplex = PETSC_TRUE; break;
329   }
330   isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
331   if (isSimplex) deg = PetscMin(deg, 3); /* Permutation not coded for degree higher than 3 */
332   /* Create space */
333   PetscCall(PetscSpaceCreate(comm, &P));
334   PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
335   PetscCall(PetscSpacePolynomialSetTensor(P, isTensor));
336   PetscCall(PetscSpaceSetNumComponents(P, dof));
337   PetscCall(PetscSpaceSetNumVariables(P, dim));
338   PetscCall(PetscSpaceSetDegree(P, deg, PETSC_DETERMINE));
339   PetscCall(PetscSpaceSetUp(P));
340   /* Create dual space */
341   PetscCall(PetscDualSpaceCreate(comm, &Q));
342   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
343   PetscCall(PetscDualSpaceLagrangeSetTensor(Q, isTensor));
344   PetscCall(PetscDualSpaceLagrangeSetContinuity(Q, continuity));
345   PetscCall(PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0));
346   PetscCall(PetscDualSpaceSetNumComponents(Q, dof));
347   PetscCall(PetscDualSpaceSetOrder(Q, deg));
348   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
349   PetscCall(PetscDualSpaceSetDM(Q, K));
350   PetscCall(DMDestroy(&K));
351   PetscCall(PetscDualSpaceSetUp(Q));
352   /* Create quadrature */
353   if (isSimplex) {
354     PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1, +1, &q));
355     PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, deg + 1, -1, +1, &fq));
356   } else {
357     PetscCall(PetscDTGaussTensorQuadrature(dim, 1, deg + 1, -1, +1, &q));
358     PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, deg + 1, -1, +1, &fq));
359   }
360   /* Create finite element */
361   PetscCall(PetscFECreate(comm, fem));
362   PetscCall(PetscSNPrintf(name, 32, "L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, deg));
363   PetscCall(PetscObjectSetName((PetscObject)*fem, name));
364   PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
365   PetscCall(PetscFESetNumComponents(*fem, dof));
366   PetscCall(PetscFESetBasisSpace(*fem, P));
367   PetscCall(PetscFESetDualSpace(*fem, Q));
368   PetscCall(PetscFESetQuadrature(*fem, q));
369   PetscCall(PetscFESetFaceQuadrature(*fem, fq));
370   PetscCall(PetscFESetUp(*fem));
371 
372   /* Both MFEM and PETSc are lexicographic, but PLEX stores the swapped cone */
373   *perm = NULL;
374   if (isSimplex && dim == 3) {
375     PetscInt celldofs, *pidx;
376 
377     PetscCall(PetscDualSpaceGetDimension(Q, &celldofs));
378     celldofs /= dof;
379     PetscCall(PetscMalloc1(celldofs, &pidx));
380     switch (celldofs) {
381     case 4:
382       pidx[0] = 2;
383       pidx[1] = 0;
384       pidx[2] = 1;
385       pidx[3] = 3;
386       break;
387     case 10:
388       pidx[0] = 5;
389       pidx[1] = 3;
390       pidx[2] = 0;
391       pidx[3] = 4;
392       pidx[4] = 1;
393       pidx[5] = 2;
394       pidx[6] = 8;
395       pidx[7] = 6;
396       pidx[8] = 7;
397       pidx[9] = 9;
398       break;
399     case 20:
400       pidx[0]  = 9;
401       pidx[1]  = 7;
402       pidx[2]  = 4;
403       pidx[3]  = 0;
404       pidx[4]  = 8;
405       pidx[5]  = 5;
406       pidx[6]  = 1;
407       pidx[7]  = 6;
408       pidx[8]  = 2;
409       pidx[9]  = 3;
410       pidx[10] = 15;
411       pidx[11] = 13;
412       pidx[12] = 10;
413       pidx[13] = 14;
414       pidx[14] = 11;
415       pidx[15] = 12;
416       pidx[16] = 18;
417       pidx[17] = 16;
418       pidx[18] = 17;
419       pidx[19] = 19;
420       break;
421     default: SETERRQ(comm, PETSC_ERR_SUP, "Unhandled degree,dof pair %" PetscInt_FMT ",%" PetscInt_FMT, deg, celldofs); break;
422     }
423     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dof, celldofs, pidx, PETSC_OWN_POINTER, perm));
424   }
425 
426   /* Cleanup */
427   PetscCall(PetscSpaceDestroy(&P));
428   PetscCall(PetscDualSpaceDestroy(&Q));
429   PetscCall(PetscQuadratureDestroy(&q));
430   PetscCall(PetscQuadratureDestroy(&fq));
431   PetscFunctionReturn(0);
432 }
433 
434 /*
435    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
436    Higher order meshes are also supported
437 */
438 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) {
439   DMLabel            label;
440   PetscSection       coordSection, coordSectionCell, parentSection, hoSection = NULL;
441   Vec                coordinates, coordinatesCell, hovec;
442   const PetscScalar *array;
443   PetscInt           bf, p, sdim, dim, depth, novl, minl;
444   PetscInt           cStart, cEnd, vStart, vEnd, nvert;
445   PetscMPIInt        size;
446   PetscBool          localized, isascii;
447   PetscBool          enable_mfem, enable_boundary, enable_ncmesh, view_ovl = PETSC_FALSE;
448   PetscBT            pown, vown;
449   PetscContainer     glvis_container;
450   PetscBool          cellvertex = PETSC_FALSE, enabled = PETSC_TRUE;
451   PetscBool          enable_emark, enable_bmark;
452   const char        *fmt;
453   char               emark[64] = "", bmark[64] = "";
454 
455   PetscFunctionBegin;
456   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
457   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
458   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
459   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer must be of type VIEWERASCII");
460   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
461   PetscCheck(size <= 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Use single sequential viewers for parallel visualization");
462   PetscCall(DMGetDimension(dm, &dim));
463   PetscCall(DMPlexGetDepth(dm, &depth));
464 
465   /* get container: determines if a process visualizes is portion of the data or not */
466   PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
467   PetscCheck(glvis_container, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis container");
468   {
469     PetscViewerGLVisInfo glvis_info;
470     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info));
471     enabled = glvis_info->enabled;
472     fmt     = glvis_info->fmt;
473   }
474 
475   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh */
476   PetscCall(PetscObjectQuery((PetscObject)dm, "_glvis_mesh_coords", (PetscObject *)&hovec));
477   PetscCall(PetscObjectReference((PetscObject)hovec));
478   if (!hovec) {
479     DM           cdm;
480     PetscFE      disc;
481     PetscClassId classid;
482 
483     PetscCall(DMGetCoordinateDM(dm, &cdm));
484     PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&disc));
485     PetscCall(PetscObjectGetClassId((PetscObject)disc, &classid));
486     if (classid == PETSCFE_CLASSID) {
487       DM      hocdm;
488       PetscFE hodisc;
489       Vec     vec;
490       Mat     mat;
491       char    name[32], fec_type[64];
492       IS      perm = NULL;
493 
494       PetscCall(GLVisCreateFE(disc, name, &hodisc, &perm));
495       PetscCall(DMClone(cdm, &hocdm));
496       PetscCall(DMSetField(hocdm, 0, NULL, (PetscObject)hodisc));
497       PetscCall(PetscFEDestroy(&hodisc));
498       PetscCall(DMCreateDS(hocdm));
499 
500       PetscCall(DMGetCoordinates(dm, &vec));
501       PetscCall(DMCreateGlobalVector(hocdm, &hovec));
502       PetscCall(DMCreateInterpolation(cdm, hocdm, &mat, NULL));
503       PetscCall(MatInterpolate(mat, vec, hovec));
504       PetscCall(MatDestroy(&mat));
505       PetscCall(DMGetLocalSection(hocdm, &hoSection));
506       PetscCall(PetscSectionSetClosurePermutation(hoSection, (PetscObject)hocdm, depth, perm));
507       PetscCall(ISDestroy(&perm));
508       PetscCall(DMDestroy(&hocdm));
509       PetscCall(PetscSNPrintf(fec_type, sizeof(fec_type), "FiniteElementCollection: %s", name));
510       PetscCall(PetscObjectSetName((PetscObject)hovec, fec_type));
511     }
512   }
513 
514   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
515   PetscCall(DMPlexGetGhostCellStratum(dm, &p, NULL));
516   if (p >= 0) cEnd = p;
517   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
518   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
519   PetscCall(DMGetCoordinateSection(dm, &coordSection));
520   PetscCall(DMGetCoordinateDim(dm, &sdim));
521   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
522   PetscCheck(coordinates || hovec, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Missing local coordinates vector");
523 
524   /*
525      a couple of sections of the mesh specification are disabled
526        - boundary: the boundary is not needed for proper mesh visualization unless we want to visualize boundary attributes or we have high-order coordinates in 3D (topologically)
527        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
528                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
529   */
530   enable_boundary = PETSC_FALSE;
531   enable_ncmesh   = PETSC_FALSE;
532   enable_mfem     = PETSC_FALSE;
533   enable_emark    = PETSC_FALSE;
534   enable_bmark    = PETSC_FALSE;
535   /* I'm tired of problems with negative values in the markers, disable them */
536   PetscOptionsBegin(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->prefix, "GLVis PetscViewer DMPlex Options", "PetscViewer");
537   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary", "Enable boundary section in mesh representation", NULL, enable_boundary, &enable_boundary, NULL));
538   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh", "Enable vertex_parents section in mesh representation (allows derefinement)", NULL, enable_ncmesh, &enable_ncmesh, NULL));
539   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem", "Dump a mesh that can be used with MFEM's FiniteElementSpaces", NULL, enable_mfem, &enable_mfem, NULL));
540   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_overlap", "Include overlap region in local meshes", NULL, view_ovl, &view_ovl, NULL));
541   PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_emarker", "String for the material id label", NULL, emark, emark, sizeof(emark), &enable_emark));
542   PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_bmarker", "String for the boundary id label", NULL, bmark, bmark, sizeof(bmark), &enable_bmark));
543   PetscOptionsEnd();
544   if (enable_bmark) enable_boundary = PETSC_TRUE;
545 
546   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
547   PetscCheck(!enable_ncmesh || size == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not supported in parallel");
548   PetscCheck(!enable_boundary || depth < 0 || dim == depth, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG,
549              "Mesh must be interpolated. "
550              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
551   PetscCheck(!enable_ncmesh || depth < 0 || dim == depth, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG,
552              "Mesh must be interpolated. "
553              "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
554   if (depth >= 0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
555     PetscCheck(depth == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported depth %" PetscInt_FMT ". You should interpolate the mesh first", depth);
556     cellvertex = PETSC_TRUE;
557   }
558 
559   /* Identify possible cells in the overlap */
560   novl = 0;
561   pown = NULL;
562   if (size > 1) {
563     IS              globalNum = NULL;
564     const PetscInt *gNum;
565     PetscBool       ovl = PETSC_FALSE;
566 
567     PetscCall(PetscObjectQuery((PetscObject)dm, "_glvis_plex_gnum", (PetscObject *)&globalNum));
568     if (!globalNum) {
569       if (view_ovl) {
570         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), cEnd - cStart, 0, 1, &globalNum));
571       } else {
572         PetscCall(DMPlexCreateCellNumbering_Internal(dm, PETSC_TRUE, &globalNum));
573       }
574       PetscCall(PetscObjectCompose((PetscObject)dm, "_glvis_plex_gnum", (PetscObject)globalNum));
575       PetscCall(PetscObjectDereference((PetscObject)globalNum));
576     }
577     PetscCall(ISGetIndices(globalNum, &gNum));
578     for (p = cStart; p < cEnd; p++) {
579       if (gNum[p - cStart] < 0) {
580         ovl = PETSC_TRUE;
581         novl++;
582       }
583     }
584     if (ovl) {
585       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
586          TODO: garbage collector? attach pown to dm?  */
587       PetscCall(PetscBTCreate(cEnd - cStart, &pown));
588       for (p = cStart; p < cEnd; p++) {
589         if (gNum[p - cStart] < 0) continue;
590         else PetscCall(PetscBTSet(pown, p - cStart));
591       }
592     }
593     PetscCall(ISRestoreIndices(globalNum, &gNum));
594   }
595 
596   /* vertex_parents (Non-conforming meshes) */
597   parentSection = NULL;
598   if (enable_ncmesh) {
599     PetscCall(DMPlexGetTree(dm, &parentSection, NULL, NULL, NULL, NULL));
600     enable_ncmesh = (PetscBool)(enable_ncmesh && parentSection);
601   }
602   /* return if this process is disabled */
603   if (!enabled) {
604     PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh %s\n", enable_ncmesh ? "v1.1" : "v1.0"));
605     PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n"));
606     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim));
607     PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
608     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
609     PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
610     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
611     PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
612     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
613     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
614     PetscCall(PetscBTDestroy(&pown));
615     PetscCall(VecDestroy(&hovec));
616     PetscFunctionReturn(0);
617   }
618 
619   if (enable_mfem) {
620     if (localized && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
621       PetscInt     vpc = 0;
622       char         fec[64];
623       PetscInt     vids[8] = {0, 1, 2, 3, 4, 5, 6, 7};
624       PetscInt     hexv[8] = {0, 1, 3, 2, 4, 5, 7, 6}, tetv[4] = {0, 1, 2, 3};
625       PetscInt     quadv[8] = {0, 1, 3, 2}, triv[3] = {0, 1, 2};
626       PetscInt    *dof = NULL;
627       PetscScalar *array, *ptr;
628 
629       PetscCall(PetscSNPrintf(fec, sizeof(fec), "FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P1", dim));
630       if (cEnd - cStart) {
631         PetscInt fpc;
632 
633         PetscCall(DMPlexGetConeSize(dm, cStart, &fpc));
634         switch (dim) {
635         case 1:
636           vpc = 2;
637           dof = hexv;
638           break;
639         case 2:
640           switch (fpc) {
641           case 3:
642             vpc = 3;
643             dof = triv;
644             break;
645           case 4:
646             vpc = 4;
647             dof = quadv;
648             break;
649           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: faces per cell %" PetscInt_FMT, fpc);
650           }
651           break;
652         case 3:
653           switch (fpc) {
654           case 4: /* TODO: still need to understand L2 ordering for tets */
655             vpc = 4;
656             dof = tetv;
657             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled tethraedral case");
658           case 6:
659             PetscCheck(!cellvertex, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: vertices per cell %" PetscInt_FMT, fpc);
660             vpc = 8;
661             dof = hexv;
662             break;
663           case 8:
664             PetscCheck(cellvertex, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: faces per cell %" PetscInt_FMT, fpc);
665             vpc = 8;
666             dof = hexv;
667             break;
668           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: faces per cell %" PetscInt_FMT, fpc);
669           }
670           break;
671         default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unhandled dim");
672         }
673         PetscCall(DMPlexReorderCell(dm, cStart, vids));
674       }
675       PetscCheck(dof, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing dofs");
676       PetscCall(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart - novl) * vpc * sdim, &hovec));
677       PetscCall(PetscObjectSetName((PetscObject)hovec, fec));
678       PetscCall(VecGetArray(hovec, &array));
679       ptr = array;
680       for (p = cStart; p < cEnd; p++) {
681         PetscInt     csize, v, d;
682         PetscScalar *vals = NULL;
683 
684         if (PetscUnlikely(pown && !PetscBTLookup(pown, p - cStart))) continue;
685         PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, p, &csize, &vals));
686         PetscCheck(csize == vpc * sdim || csize == vpc * sdim * 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported closure size %" PetscInt_FMT " (vpc %" PetscInt_FMT ", sdim %" PetscInt_FMT ")", csize, vpc, sdim);
687         for (v = 0; v < vpc; v++) {
688           for (d = 0; d < sdim; d++) ptr[sdim * dof[v] + d] = vals[sdim * vids[v] + d];
689         }
690         ptr += vpc * sdim;
691         PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, p, &csize, &vals));
692       }
693       PetscCall(VecRestoreArray(hovec, &array));
694     }
695   }
696   /* if we have high-order coordinates in 3D, we need to specify the boundary */
697   if (hovec && dim == 3) enable_boundary = PETSC_TRUE;
698 
699   /* header */
700   PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh %s\n", enable_ncmesh ? "v1.1" : "v1.0"));
701 
702   /* topological dimension */
703   PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n"));
704   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim));
705 
706   /* elements */
707   minl  = 1;
708   label = NULL;
709   if (enable_emark) {
710     PetscInt lminl = PETSC_MAX_INT;
711 
712     PetscCall(DMGetLabel(dm, emark, &label));
713     if (label) {
714       IS       vals;
715       PetscInt ldef;
716 
717       PetscCall(DMLabelGetDefaultValue(label, &ldef));
718       PetscCall(DMLabelGetValueIS(label, &vals));
719       PetscCall(ISGetMinMax(vals, &lminl, NULL));
720       PetscCall(ISDestroy(&vals));
721       lminl = PetscMin(ldef, lminl);
722     }
723     PetscCall(MPIU_Allreduce(&lminl, &minl, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
724     if (minl == PETSC_MAX_INT) minl = 1;
725   }
726   PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
727   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", cEnd - cStart - novl));
728   for (p = cStart; p < cEnd; p++) {
729     PetscInt vids[8];
730     PetscInt i, nv = 0, cid = -1, mid = 1;
731 
732     if (PetscUnlikely(pown && !PetscBTLookup(pown, p - cStart))) continue;
733     PetscCall(DMPlexGetPointMFEMCellID_Internal(dm, label, minl, p, &mid, &cid));
734     PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm, p, (localized && !hovec) ? coordSection : NULL, &nv, vids));
735     PetscCall(DMPlexReorderCell(dm, p, vids));
736     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT, mid, cid));
737     for (i = 0; i < nv; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, vids[i]));
738     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
739   }
740 
741   /* boundary */
742   PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
743   if (!enable_boundary) {
744     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
745   } else {
746     DMLabel  perLabel;
747     PetscBT  bfaces;
748     PetscInt fStart, fEnd, *fcells;
749 
750     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
751     PetscCall(PetscBTCreate(fEnd - fStart, &bfaces));
752     PetscCall(DMPlexGetMaxSizes(dm, NULL, &p));
753     PetscCall(PetscMalloc1(p, &fcells));
754     PetscCall(DMGetLabel(dm, "glvis_periodic_cut", &perLabel));
755     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
756       PetscCall(DMCreateLabel(dm, "glvis_periodic_cut"));
757       PetscCall(DMGetLabel(dm, "glvis_periodic_cut", &perLabel));
758       PetscCall(DMLabelSetDefaultValue(perLabel, 1));
759       PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell));
760       PetscCall(DMGetCellCoordinatesLocal(dm, &coordinatesCell));
761       for (p = cStart; p < cEnd; p++) {
762         DMPolytopeType cellType;
763         PetscInt       dof;
764 
765         PetscCall(DMPlexGetCellType(dm, p, &cellType));
766         PetscCall(PetscSectionGetDof(coordSectionCell, p, &dof));
767         if (dof) {
768           PetscInt     uvpc, v, csize, csizeCell, cellClosureSize, *cellClosure = NULL, *vidxs = NULL;
769           PetscScalar *vals = NULL, *valsCell = NULL;
770 
771           uvpc = DMPolytopeTypeGetNumVertices(cellType);
772           PetscCheck(dof % sdim == 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Incompatible number of cell dofs %" PetscInt_FMT " and space dimension %" PetscInt_FMT, dof, sdim);
773           PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, p, &csize, &vals));
774           PetscCall(DMPlexVecGetClosure(dm, coordSectionCell, coordinatesCell, p, &csizeCell, &valsCell));
775           PetscCheck(csize == csizeCell, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cell %" PetscInt_FMT " has invalid localized coordinates", p);
776           PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cellClosureSize, &cellClosure));
777           for (v = 0; v < cellClosureSize; v++)
778             if (cellClosure[2 * v] >= vStart && cellClosure[2 * v] < vEnd) {
779               vidxs = cellClosure + 2 * v;
780               break;
781             }
782           PetscCheck(vidxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing vertices");
783           for (v = 0; v < uvpc; v++) {
784             PetscInt s;
785 
786             for (s = 0; s < sdim; s++) {
787               if (PetscAbsScalar(vals[v * sdim + s] - valsCell[v * sdim + s]) > PETSC_MACHINE_EPSILON) PetscCall(DMLabelSetValue(perLabel, vidxs[2 * v], 2));
788             }
789           }
790           PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cellClosureSize, &cellClosure));
791           PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, p, &csize, &vals));
792           PetscCall(DMPlexVecRestoreClosure(dm, coordSectionCell, coordinatesCell, p, &csizeCell, &valsCell));
793         }
794       }
795       if (dim > 1) {
796         PetscInt eEnd, eStart;
797 
798         PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
799         for (p = eStart; p < eEnd; p++) {
800           const PetscInt *cone;
801           PetscInt        coneSize, i;
802           PetscBool       ispe = PETSC_TRUE;
803 
804           PetscCall(DMPlexGetCone(dm, p, &cone));
805           PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
806           for (i = 0; i < coneSize; i++) {
807             PetscInt v;
808 
809             PetscCall(DMLabelGetValue(perLabel, cone[i], &v));
810             ispe = (PetscBool)(ispe && (v == 2));
811           }
812           if (ispe && coneSize) {
813             PetscInt        ch, numChildren;
814             const PetscInt *children;
815 
816             PetscCall(DMLabelSetValue(perLabel, p, 2));
817             PetscCall(DMPlexGetTreeChildren(dm, p, &numChildren, &children));
818             for (ch = 0; ch < numChildren; ch++) PetscCall(DMLabelSetValue(perLabel, children[ch], 2));
819           }
820         }
821         if (dim > 2) {
822           for (p = fStart; p < fEnd; p++) {
823             const PetscInt *cone;
824             PetscInt        coneSize, i;
825             PetscBool       ispe = PETSC_TRUE;
826 
827             PetscCall(DMPlexGetCone(dm, p, &cone));
828             PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
829             for (i = 0; i < coneSize; i++) {
830               PetscInt v;
831 
832               PetscCall(DMLabelGetValue(perLabel, cone[i], &v));
833               ispe = (PetscBool)(ispe && (v == 2));
834             }
835             if (ispe && coneSize) {
836               PetscInt        ch, numChildren;
837               const PetscInt *children;
838 
839               PetscCall(DMLabelSetValue(perLabel, p, 2));
840               PetscCall(DMPlexGetTreeChildren(dm, p, &numChildren, &children));
841               for (ch = 0; ch < numChildren; ch++) PetscCall(DMLabelSetValue(perLabel, children[ch], 2));
842             }
843           }
844         }
845       }
846     }
847     for (p = fStart; p < fEnd; p++) {
848       const PetscInt *support;
849       PetscInt        supportSize;
850       PetscBool       isbf = PETSC_FALSE;
851 
852       PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
853       if (pown) {
854         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
855         PetscInt  i;
856 
857         PetscCall(DMPlexGetSupport(dm, p, &support));
858         for (i = 0; i < supportSize; i++) {
859           if (PetscLikely(PetscBTLookup(pown, support[i] - cStart))) has_owned = PETSC_TRUE;
860           else has_ghost = PETSC_TRUE;
861         }
862         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
863       } else {
864         isbf = (PetscBool)(supportSize == 1);
865       }
866       if (!isbf && perLabel) {
867         const PetscInt *cone;
868         PetscInt        coneSize, i;
869 
870         PetscCall(DMPlexGetCone(dm, p, &cone));
871         PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
872         isbf = PETSC_TRUE;
873         for (i = 0; i < coneSize; i++) {
874           PetscInt v, d;
875 
876           PetscCall(DMLabelGetValue(perLabel, cone[i], &v));
877           PetscCall(DMLabelGetDefaultValue(perLabel, &d));
878           isbf = (PetscBool)(isbf && v != d);
879         }
880       }
881       if (isbf) PetscCall(PetscBTSet(bfaces, p - fStart));
882     }
883     /* count boundary faces */
884     for (p = fStart, bf = 0; p < fEnd; p++) {
885       if (PetscUnlikely(PetscBTLookup(bfaces, p - fStart))) {
886         const PetscInt *support;
887         PetscInt        supportSize, c;
888 
889         PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
890         PetscCall(DMPlexGetSupport(dm, p, &support));
891         for (c = 0; c < supportSize; c++) {
892           const PetscInt *cone;
893           PetscInt        cell, cl, coneSize;
894 
895           cell = support[c];
896           if (pown && PetscUnlikely(!PetscBTLookup(pown, cell - cStart))) continue;
897           PetscCall(DMPlexGetCone(dm, cell, &cone));
898           PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
899           for (cl = 0; cl < coneSize; cl++) {
900             if (cone[cl] == p) {
901               bf += 1;
902               break;
903             }
904           }
905         }
906       }
907     }
908     minl  = 1;
909     label = NULL;
910     if (enable_bmark) {
911       PetscInt lminl = PETSC_MAX_INT;
912 
913       PetscCall(DMGetLabel(dm, bmark, &label));
914       if (label) {
915         IS       vals;
916         PetscInt ldef;
917 
918         PetscCall(DMLabelGetDefaultValue(label, &ldef));
919         PetscCall(DMLabelGetValueIS(label, &vals));
920         PetscCall(ISGetMinMax(vals, &lminl, NULL));
921         PetscCall(ISDestroy(&vals));
922         lminl = PetscMin(ldef, lminl);
923       }
924       PetscCall(MPIU_Allreduce(&lminl, &minl, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
925       if (minl == PETSC_MAX_INT) minl = 1;
926     }
927     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", bf));
928     for (p = fStart; p < fEnd; p++) {
929       if (PetscUnlikely(PetscBTLookup(bfaces, p - fStart))) {
930         const PetscInt *support;
931         PetscInt        supportSize, c, nc = 0;
932 
933         PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
934         PetscCall(DMPlexGetSupport(dm, p, &support));
935         if (pown) {
936           for (c = 0; c < supportSize; c++) {
937             if (PetscLikely(PetscBTLookup(pown, support[c] - cStart))) fcells[nc++] = support[c];
938           }
939         } else
940           for (c = 0; c < supportSize; c++) fcells[nc++] = support[c];
941         for (c = 0; c < nc; c++) {
942           const DMPolytopeType *faceTypes;
943           DMPolytopeType        cellType;
944           const PetscInt       *faceSizes, *cone;
945           PetscInt              vids[8], *faces, st, i, coneSize, cell, cl, nv, cid = -1, mid = -1;
946 
947           cell = fcells[c];
948           PetscCall(DMPlexGetCone(dm, cell, &cone));
949           PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
950           for (cl = 0; cl < coneSize; cl++)
951             if (cone[cl] == p) break;
952           if (cl == coneSize) continue;
953 
954           /* face material id and type */
955           PetscCall(DMPlexGetPointMFEMCellID_Internal(dm, label, minl, p, &mid, &cid));
956           PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT, mid, cid));
957           /* vertex ids */
958           PetscCall(DMPlexGetCellType(dm, cell, &cellType));
959           PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm, cell, (localized && !hovec) ? coordSection : NULL, &nv, vids));
960           PetscCall(DMPlexGetRawFaces_Internal(dm, cellType, vids, NULL, &faceTypes, &faceSizes, (const PetscInt **)&faces));
961           st = 0;
962           for (i = 0; i < cl; i++) st += faceSizes[i];
963           PetscCall(DMPlexInvertCell(faceTypes[cl], faces + st));
964           for (i = 0; i < faceSizes[cl]; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, faces[st + i]));
965           PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
966           PetscCall(DMPlexRestoreRawFaces_Internal(dm, cellType, vids, NULL, &faceTypes, &faceSizes, (const PetscInt **)&faces));
967           bf -= 1;
968         }
969       }
970     }
971     PetscCheck(!bf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Remaining boundary faces %" PetscInt_FMT, bf);
972     PetscCall(PetscBTDestroy(&bfaces));
973     PetscCall(PetscFree(fcells));
974   }
975 
976   /* mark owned vertices */
977   vown = NULL;
978   if (pown) {
979     PetscCall(PetscBTCreate(vEnd - vStart, &vown));
980     for (p = cStart; p < cEnd; p++) {
981       PetscInt i, closureSize, *closure = NULL;
982 
983       if (PetscUnlikely(!PetscBTLookup(pown, p - cStart))) continue;
984       PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure));
985       for (i = 0; i < closureSize; i++) {
986         const PetscInt pp = closure[2 * i];
987 
988         if (pp >= vStart && pp < vEnd) PetscCall(PetscBTSet(vown, pp - vStart));
989       }
990       PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure));
991     }
992   }
993 
994   if (parentSection) {
995     PetscInt vp, gvp;
996 
997     for (vp = 0, p = vStart; p < vEnd; p++) {
998       DMLabel  dlabel;
999       PetscInt parent, depth;
1000 
1001       if (PetscUnlikely(vown && !PetscBTLookup(vown, p - vStart))) continue;
1002       PetscCall(DMPlexGetDepthLabel(dm, &dlabel));
1003       PetscCall(DMLabelGetValue(dlabel, p, &depth));
1004       PetscCall(DMPlexGetTreeParent(dm, p, &parent, NULL));
1005       if (parent != p) vp++;
1006     }
1007     PetscCall(MPIU_Allreduce(&vp, &gvp, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
1008     if (gvp) {
1009       PetscInt   maxsupp;
1010       PetscBool *skip = NULL;
1011 
1012       PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertex_parents\n"));
1013       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", vp));
1014       PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxsupp));
1015       PetscCall(PetscMalloc1(maxsupp, &skip));
1016       for (p = vStart; p < vEnd; p++) {
1017         DMLabel  dlabel;
1018         PetscInt parent;
1019 
1020         if (PetscUnlikely(vown && !PetscBTLookup(vown, p - vStart))) continue;
1021         PetscCall(DMPlexGetDepthLabel(dm, &dlabel));
1022         PetscCall(DMPlexGetTreeParent(dm, p, &parent, NULL));
1023         if (parent != p) {
1024           PetscInt        vids[8] = {-1, -1, -1, -1, -1, -1, -1, -1}; /* silent overzealous clang static analyzer */
1025           PetscInt        i, nv, ssize, n, numChildren, depth = -1;
1026           const PetscInt *children;
1027 
1028           PetscCall(DMPlexGetConeSize(dm, parent, &ssize));
1029           switch (ssize) {
1030           case 2: /* edge */
1031             nv = 0;
1032             PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm, parent, localized ? coordSection : NULL, &nv, vids));
1033             PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, p - vStart));
1034             for (i = 0; i < nv; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, vids[i]));
1035             PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1036             vp--;
1037             break;
1038           case 4: /* face */
1039             PetscCall(DMPlexGetTreeChildren(dm, parent, &numChildren, &children));
1040             for (n = 0; n < numChildren; n++) {
1041               PetscCall(DMLabelGetValue(dlabel, children[n], &depth));
1042               if (!depth) {
1043                 const PetscInt *hvsupp, *hesupp, *cone;
1044                 PetscInt        hvsuppSize, hesuppSize, coneSize;
1045                 PetscInt        hv = children[n], he = -1, f;
1046 
1047                 PetscCall(PetscArrayzero(skip, maxsupp));
1048                 PetscCall(DMPlexGetSupportSize(dm, hv, &hvsuppSize));
1049                 PetscCall(DMPlexGetSupport(dm, hv, &hvsupp));
1050                 for (i = 0; i < hvsuppSize; i++) {
1051                   PetscInt ep;
1052                   PetscCall(DMPlexGetTreeParent(dm, hvsupp[i], &ep, NULL));
1053                   if (ep != hvsupp[i]) {
1054                     he = hvsupp[i];
1055                   } else {
1056                     skip[i] = PETSC_TRUE;
1057                   }
1058                 }
1059                 PetscCheck(he != -1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Vertex %" PetscInt_FMT " support size %" PetscInt_FMT ": hanging edge not found", hv, hvsuppSize);
1060                 PetscCall(DMPlexGetCone(dm, he, &cone));
1061                 vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
1062                 PetscCall(DMPlexGetSupportSize(dm, he, &hesuppSize));
1063                 PetscCall(DMPlexGetSupport(dm, he, &hesupp));
1064                 for (f = 0; f < hesuppSize; f++) {
1065                   PetscInt j;
1066 
1067                   PetscCall(DMPlexGetCone(dm, hesupp[f], &cone));
1068                   PetscCall(DMPlexGetConeSize(dm, hesupp[f], &coneSize));
1069                   for (j = 0; j < coneSize; j++) {
1070                     PetscInt k;
1071                     for (k = 0; k < hvsuppSize; k++) {
1072                       if (hvsupp[k] == cone[j]) {
1073                         skip[k] = PETSC_TRUE;
1074                         break;
1075                       }
1076                     }
1077                   }
1078                 }
1079                 for (i = 0; i < hvsuppSize; i++) {
1080                   if (!skip[i]) {
1081                     PetscCall(DMPlexGetCone(dm, hvsupp[i], &cone));
1082                     vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
1083                   }
1084                 }
1085                 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, hv - vStart));
1086                 for (i = 0; i < 2; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, vids[i] - vStart));
1087                 PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1088                 vp--;
1089               }
1090             }
1091             break;
1092           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Don't know how to deal with support size %" PetscInt_FMT, ssize);
1093           }
1094         }
1095       }
1096       PetscCall(PetscFree(skip));
1097     }
1098     PetscCheck(!vp, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected %" PetscInt_FMT " hanging vertices", vp);
1099   }
1100   PetscCall(PetscBTDestroy(&vown));
1101 
1102   /* vertices */
1103   if (hovec) { /* higher-order meshes */
1104     const char *fec;
1105     PetscInt    i, n, s;
1106     PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
1107     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", vEnd - vStart));
1108     PetscCall(PetscViewerASCIIPrintf(viewer, "nodes\n"));
1109     PetscCall(PetscObjectGetName((PetscObject)hovec, &fec));
1110     PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
1111     PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", fec));
1112     PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", sdim));
1113     PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: 1\n\n")); /*Ordering::byVDIM*/
1114     if (hoSection) {
1115       DM cdm;
1116 
1117       PetscCall(VecGetDM(hovec, &cdm));
1118       for (p = cStart; p < cEnd; p++) {
1119         PetscScalar *vals = NULL;
1120         PetscInt     csize;
1121 
1122         if (PetscUnlikely(pown && !PetscBTLookup(pown, p - cStart))) continue;
1123         PetscCall(DMPlexVecGetClosure(cdm, hoSection, hovec, p, &csize, &vals));
1124         PetscCheck(csize % sdim == 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Size of closure %" PetscInt_FMT " incompatible with space dimension %" PetscInt_FMT, csize, sdim);
1125         for (i = 0; i < csize / sdim; i++) {
1126           for (s = 0; s < sdim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, (double)PetscRealPart(vals[i * sdim + s])));
1127           PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1128         }
1129         PetscCall(DMPlexVecRestoreClosure(cdm, hoSection, hovec, p, &csize, &vals));
1130       }
1131     } else {
1132       PetscCall(VecGetArrayRead(hovec, &array));
1133       PetscCall(VecGetLocalSize(hovec, &n));
1134       PetscCheck(n % sdim == 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Size of local coordinate vector %" PetscInt_FMT " incompatible with space dimension %" PetscInt_FMT, n, sdim);
1135       for (i = 0; i < n / sdim; i++) {
1136         for (s = 0; s < sdim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, (double)PetscRealPart(array[i * sdim + s])));
1137         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1138       }
1139       PetscCall(VecRestoreArrayRead(hovec, &array));
1140     }
1141   } else {
1142     PetscCall(VecGetLocalSize(coordinates, &nvert));
1143     PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
1144     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", nvert / sdim));
1145     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
1146     PetscCall(VecGetArrayRead(coordinates, &array));
1147     for (p = 0; p < nvert / sdim; p++) {
1148       PetscInt s;
1149       for (s = 0; s < sdim; s++) {
1150         PetscReal v = PetscRealPart(array[p * sdim + s]);
1151 
1152         PetscCall(PetscViewerASCIIPrintf(viewer, fmt, PetscIsInfOrNanReal(v) ? 0.0 : (double)v));
1153       }
1154       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1155     }
1156     PetscCall(VecRestoreArrayRead(coordinates, &array));
1157   }
1158   PetscCall(PetscBTDestroy(&pown));
1159   PetscCall(VecDestroy(&hovec));
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) {
1164   PetscFunctionBegin;
1165   PetscCall(DMView_GLVis(dm, viewer, DMPlexView_GLVis_ASCII));
1166   PetscFunctionReturn(0);
1167 }
1168