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