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