xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 3924b6122239b9f35da0ebe3c7fa55fdb9db4750) !
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 = DMGetDefaultSection(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 = VecScatterCreate(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 = VecScatterCreate(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     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
282     ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr);
283     off  = off/sdim;
284     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
285   }
286   if (!dof) {
287     ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
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     ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
292   } else {
293     ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr);
294     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
295     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
296   }
297   *nv  = q;
298   PetscFunctionReturn(0);
299 }
300 
301 /*
302    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
303    Higher order meshes are also supported
304 */
305 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
306 {
307   DMLabel              label;
308   PetscSection         coordSection,parentSection;
309   Vec                  coordinates,hovec;
310   const PetscScalar    *array;
311   PetscInt             bf,p,sdim,dim,depth,novl;
312   PetscInt             cStart,cEnd,cEndInterior,vStart,vEnd,nvert;
313   PetscMPIInt          size;
314   PetscBool            localized,isascii;
315   PetscBool            enable_mfem,enable_boundary,enable_ncmesh;
316   PetscBT              pown,vown;
317   PetscErrorCode       ierr;
318   PetscContainer       glvis_container;
319   PetscBool            cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
320   const char           *fmt;
321   char                 emark[64] = "",bmark[64] = "";
322 
323   PetscFunctionBegin;
324   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
325   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
326   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
327   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
328   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
329   if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
330   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
331 
332   /* get container: determines if a process visualizes is portion of the data or not */
333   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
334   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
335   {
336     PetscViewerGLVisInfo glvis_info;
337     ierr    = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
338     enabled = glvis_info->enabled;
339     fmt     = glvis_info->fmt;
340   }
341   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
342   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
343   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
344   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
345   ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr);
346   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
347   if (periodic && !localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
348   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
349   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
350   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
351   if (!coordinates) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");
352 
353   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
354      DMPlex does not currently support HO meshes, so there's no API for this */
355   ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr);
356 
357   /*
358      a couple of sections of the mesh specification are disabled
359        - boundary: unless we want to visualize boundary attributes or we have a periodic mesh
360                    the boundary is not needed for proper mesh visualization
361        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
362                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
363   */
364   enable_boundary = periodic;
365   enable_ncmesh   = PETSC_FALSE;
366   enable_mfem     = PETSC_FALSE;
367   /* I'm tired of problems with negative values in the markers, disable them */
368   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
369   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr);
370   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);
371   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);
372   ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),NULL);CHKERRQ(ierr);
373   ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),NULL);CHKERRQ(ierr);
374   ierr = PetscOptionsEnd();CHKERRQ(ierr);
375   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
376   if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
377   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
378   if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
379                                                              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
380   if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
381                                                            "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
382   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
383     if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth);
384     cellvertex = PETSC_TRUE;
385   }
386 
387   /* Identify possible cells in the overlap */
388   novl = 0;
389   pown = NULL;
390   if (size > 1) {
391     IS             globalNum = NULL;
392     const PetscInt *gNum;
393     PetscBool      ovl  = PETSC_FALSE;
394 
395     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
396     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
397     for (p=cStart; p<cEnd; p++) {
398       if (gNum[p-cStart] < 0) {
399         ovl = PETSC_TRUE;
400         novl++;
401       }
402     }
403     if (ovl) {
404       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
405          TODO: garbage collector? attach pown to dm?  */
406       ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr);
407       for (p=cStart; p<cEnd; p++) {
408         if (gNum[p-cStart] < 0) continue;
409         else {
410           ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
411         }
412       }
413     }
414     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
415   }
416 
417   /* return if this process is disabled */
418   if (!enabled) {
419     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
420     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
421     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
422     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
423     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
424     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
425     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
426     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
427     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
428     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
429     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
430     PetscFunctionReturn(0);
431   }
432 
433   if (enable_mfem) {
434     if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
435       PetscInt    vpc = 0;
436       char        fec[64];
437       int         vids[8] = {0,1,2,3,4,5,6,7};
438       int         hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3};
439       int         quadv[8] = {0,1,3,2}, triv[3] = {0,1,2};
440       int         *dof = NULL;
441       PetscScalar *array,*ptr;
442 
443       ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr);
444       if (cEnd-cStart) {
445         PetscInt fpc;
446 
447         ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
448         switch(dim) {
449           case 1:
450             vpc = 2;
451             dof = hexv;
452             break;
453           case 2:
454             switch (fpc) {
455               case 3:
456                 vpc = 3;
457                 dof = triv;
458                 break;
459               case 4:
460                 vpc = 4;
461                 dof = quadv;
462                 break;
463               default:
464                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
465                 break;
466             }
467             break;
468           case 3:
469             switch (fpc) {
470               case 4: /* TODO: still need to understand L2 ordering for tets */
471                 vpc = 4;
472                 dof = tetv;
473                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
474                 break;
475               case 6:
476                 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc);
477                 vpc = 8;
478                 dof = hexv;
479                 break;
480               case 8:
481                 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
482                 vpc = 8;
483                 dof = hexv;
484                 break;
485               default:
486                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
487                 break;
488             }
489             break;
490           default:
491             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
492             break;
493         }
494         ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr);
495       }
496       if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
497       ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr);
498       ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr);
499       ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr);
500       ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr);
501       ierr = VecGetArray(hovec,&array);CHKERRQ(ierr);
502       ptr  = array;
503       for (p=cStart;p<cEnd;p++) {
504         PetscInt    csize,v,d;
505         PetscScalar *vals = NULL;
506 
507         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
508         ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
509         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);
510         for (v=0;v<vpc;v++) {
511           for (d=0;d<sdim;d++) {
512             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
513           }
514         }
515         ptr += vpc*sdim;
516         ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
517       }
518       ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr);
519     }
520   }
521 
522 
523   /* header */
524   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
525 
526   /* topological dimension */
527   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
528   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
529 
530   /* elements */
531   ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr);
532   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
533   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
534   for (p=cStart;p<cEnd;p++) {
535     int      vids[8];
536     PetscInt i,nv = 0,cid = -1,mid = 1;
537 
538     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
539     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
540     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
541     ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
542     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
543     for (i=0;i<nv;i++) {
544       ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr);
545     }
546     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
547   }
548 
549   /* boundary */
550   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
551   if (!enable_boundary) {
552     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
553   } else {
554     DMLabel  perLabel;
555     PetscBT  bfaces;
556     PetscInt fStart,fEnd,fEndInterior,*fcells;
557     PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0;
558     PetscInt fv1[]     = {0,1},
559              fv2tri[]  = {0,1,
560                           1,2,
561                           2,0},
562              fv2quad[] = {0,1,
563                           1,2,
564                           2,3,
565                           3,0},
566              fv3tet[]  = {0,1,2,
567                           0,3,1,
568                           0,2,3,
569                           2,1,3},
570              fv3hex[]  = {0,1,2,3,
571                        4,5,6,7,
572                        0,3,5,4,
573                        2,1,7,6,
574                        3,2,6,5,
575                        0,4,7,1};
576 
577     /* determine orientation of boundary mesh */
578     if (cEnd-cStart) {
579       ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
580       switch(dim) {
581         case 1:
582           if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
583           faces = fv1;
584           vpf = 1;
585           vpc = 2;
586           break;
587         case 2:
588           switch (fpc) {
589             case 3:
590               faces = fv2tri;
591               vpf   = 2;
592               vpc   = 3;
593               break;
594             case 4:
595               faces = fv2quad;
596               vpf   = 2;
597               vpc   = 4;
598               break;
599             default:
600               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
601               break;
602           }
603           break;
604         case 3:
605           switch (fpc) {
606             case 4:
607               faces = fv3tet;
608               vpf   = 3;
609               vpc   = 4;
610               break;
611             case 6:
612               faces = fv3hex;
613               vpf   = 4;
614               vpc   = 8;
615               break;
616             default:
617               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
618               break;
619           }
620           break;
621         default:
622           SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
623           break;
624       }
625     }
626     ierr = DMPlexGetHybridBounds(dm,NULL,&fEndInterior,NULL,NULL);CHKERRQ(ierr);
627     ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
628     fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
629     ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr);
630     ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr);
631     ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr);
632     ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
633     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
634       ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr);
635       ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
636       ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr);
637       for (p=cStart;p<cEnd;p++) {
638         PetscInt dof;
639         ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
640         if (dof) {
641           PetscInt    v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
642           PetscScalar *vals = NULL;
643 
644           if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
645           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);
646           ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
647           ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
648           for (v=0;v<cellClosureSize;v++)
649             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
650               vidxs = cellClosure + 2*v;
651               break;
652             }
653           if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
654           for (v=0;v<vpc;v++) {
655             PetscInt s;
656 
657             for (s=0;s<sdim;s++) {
658               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+vpc*sdim])>PETSC_MACHINE_EPSILON) {
659                 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr);
660               }
661             }
662           }
663           ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
664           ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
665         }
666       }
667       if (dim > 1) {
668         PetscInt eEnd,eStart,eEndInterior;
669 
670         ierr = DMPlexGetHybridBounds(dm,NULL,NULL,&eEndInterior,NULL);CHKERRQ(ierr);
671         ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr);
672         eEnd = eEndInterior < 0 ? eEnd : eEndInterior;
673         for (p=eStart;p<eEnd;p++) {
674           const PetscInt *cone;
675           PetscInt       coneSize,i;
676           PetscBool      ispe = PETSC_TRUE;
677 
678           ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
679           ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
680           for (i=0;i<coneSize;i++) {
681             PetscInt v;
682 
683             ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
684             ispe = (PetscBool)(ispe && (v==2));
685           }
686           if (ispe && coneSize) {
687             ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
688           }
689         }
690         if (dim > 2) {
691           for (p=fStart;p<fEnd;p++) {
692             const PetscInt *cone;
693             PetscInt       coneSize,i;
694             PetscBool      ispe = PETSC_TRUE;
695 
696             ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
697             ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
698             for (i=0;i<coneSize;i++) {
699               PetscInt v;
700 
701               ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
702               ispe = (PetscBool)(ispe && (v==2));
703             }
704             if (ispe && coneSize) {
705               ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
706             }
707           }
708         }
709       }
710     }
711     for (p=fStart;p<fEnd;p++) {
712       const PetscInt *support;
713       PetscInt       supportSize;
714       PetscBool      isbf = PETSC_FALSE;
715 
716       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
717       if (pown) {
718         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
719         PetscInt  i;
720 
721         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
722         for (i=0;i<supportSize;i++) {
723           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
724           else has_ghost = PETSC_TRUE;
725         }
726         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
727       } else {
728         isbf = (PetscBool)(supportSize == 1);
729       }
730       if (!isbf && perLabel) {
731         const PetscInt *cone;
732         PetscInt       coneSize,i;
733 
734         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
735         ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
736         isbf = PETSC_TRUE;
737         for (i=0;i<coneSize;i++) {
738           PetscInt v,d;
739 
740           ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
741           ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr);
742           isbf = (PetscBool)(isbf && v != d);
743         }
744       }
745       if (isbf) {
746         ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr);
747       }
748     }
749     /* count boundary faces */
750     for (p=fStart,bf=0;p<fEnd;p++) {
751       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
752         const PetscInt *support;
753         PetscInt       supportSize,c;
754 
755         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
756         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
757         if (pown) {
758           for (c=0;c<supportSize;c++) {
759             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
760               bf++;
761             }
762           }
763         } else bf += supportSize;
764       }
765     }
766     ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr);
767     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
768     for (p=fStart;p<fEnd;p++) {
769       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
770         const PetscInt *support;
771         PetscInt       supportSize,c,nc = 0;
772 
773         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
774         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
775         if (pown) {
776           for (c=0;c<supportSize;c++) {
777             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
778               fcells[nc++] = support[c];
779             }
780           }
781         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
782         for (c=0;c<nc;c++) {
783           const    PetscInt *cone;
784           int      vids[8];
785           PetscInt i,cell,cl,nv,cid = -1,mid = -1;
786 
787           cell = fcells[c];
788           ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr);
789           for (cl=0;cl<fpc;cl++)
790             if (cone[cl] == p)
791               break;
792 
793           /* face material id and type */
794           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
795           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
796           /* vertex ids */
797           ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
798           for (i=0;i<vpf;i++) {
799             ierr = PetscViewerASCIIPrintf(viewer," %D",vids[faces[cl*vpf+i]]);CHKERRQ(ierr);
800           }
801           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
802         }
803         bf = bf-nc;
804       }
805     }
806     if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
807     ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr);
808     ierr = PetscFree(fcells);CHKERRQ(ierr);
809   }
810 
811   /* mark owned vertices */
812   vown = NULL;
813   if (pown) {
814     ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
815     for (p=cStart;p<cEnd;p++) {
816       PetscInt i,closureSize,*closure = NULL;
817 
818       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
819       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
820       for (i=0;i<closureSize;i++) {
821         const PetscInt pp = closure[2*i];
822 
823         if (pp >= vStart && pp < vEnd) {
824           ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr);
825         }
826       }
827       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
828     }
829   }
830 
831   /* vertex_parents (Non-conforming meshes) */
832   parentSection  = NULL;
833   if (enable_ncmesh) {
834     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
835   }
836   if (parentSection) {
837     PetscInt vp,gvp;
838 
839     for (vp=0,p=vStart;p<vEnd;p++) {
840       DMLabel  dlabel;
841       PetscInt parent,depth;
842 
843       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
844       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
845       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
846       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
847       if (parent != p) vp++;
848     }
849     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
850     if (gvp) {
851       PetscInt  maxsupp;
852       PetscBool *skip = NULL;
853 
854       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
855       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
856       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
857       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
858       for (p=vStart;p<vEnd;p++) {
859         DMLabel  dlabel;
860         PetscInt parent;
861 
862         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
863         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
864         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
865         if (parent != p) {
866           int            vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
867           PetscInt       i,nv,ssize,n,numChildren,depth = -1;
868           const PetscInt *children;
869 
870           ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr);
871           switch (ssize) {
872             case 2: /* edge */
873               nv   = 0;
874               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
875               ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
876               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
877               for (i=0;i<nv;i++) {
878                 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr);
879               }
880               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
881               vp--;
882               break;
883             case 4: /* face */
884               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
885               for (n=0;n<numChildren;n++) {
886                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
887                 if (!depth) {
888                   const PetscInt *hvsupp,*hesupp,*cone;
889                   PetscInt       hvsuppSize,hesuppSize,coneSize;
890                   PetscInt       hv = children[n],he = -1,f;
891 
892                   ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
893                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
894                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
895                   for (i=0;i<hvsuppSize;i++) {
896                     PetscInt ep;
897                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
898                     if (ep != hvsupp[i]) {
899                       he = hvsupp[i];
900                     } else {
901                       skip[i] = PETSC_TRUE;
902                     }
903                   }
904                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
905                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
906                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
907                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
908                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
909                   for (f=0;f<hesuppSize;f++) {
910                     PetscInt j;
911 
912                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
913                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
914                     for (j=0;j<coneSize;j++) {
915                       PetscInt k;
916                       for (k=0;k<hvsuppSize;k++) {
917                         if (hvsupp[k] == cone[j]) {
918                           skip[k] = PETSC_TRUE;
919                           break;
920                         }
921                       }
922                     }
923                   }
924                   for (i=0;i<hvsuppSize;i++) {
925                     if (!skip[i]) {
926                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
927                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
928                     }
929                   }
930                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
931                   for (i=0;i<2;i++) {
932                     ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)(vids[i]-vStart));CHKERRQ(ierr);
933                   }
934                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
935                   vp--;
936                 }
937               }
938               break;
939             default:
940               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize);
941           }
942         }
943       }
944       ierr = PetscFree(skip);CHKERRQ(ierr);
945     }
946     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
947   }
948   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
949   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
950 
951   /* vertices */
952   if (hovec) { /* higher-order meshes */
953     const char *fec;
954     PetscInt   i,n;
955 
956     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
957     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
958     ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
959     ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr);
960     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
961     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr);
962     ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
963     ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
964     ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr);
965     ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr);
966     if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
967     for (i=0;i<n/sdim;i++) {
968       PetscInt s;
969       for (s=0;s<sdim;s++) {
970         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr);
971       }
972       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
973     }
974     ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr);
975   } else {
976     ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr);
977     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
978     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr);
979     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
980     ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
981     for (p=0;p<nvert/sdim;p++) {
982       PetscInt s;
983       for (s=0;s<sdim;s++) {
984         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr);
985       }
986       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
987     }
988     ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
989   }
990   PetscFunctionReturn(0);
991 }
992 
993 /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
994 PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
995 {
996   PetscErrorCode ierr;
997   PetscBool      isglvis,isascii;
998 
999   PetscFunctionBegin;
1000   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1001   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1002   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
1003   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1004   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
1005   if (isglvis) {
1006     PetscViewer          view;
1007     PetscViewerGLVisType type;
1008 
1009     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
1010     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
1011     if (view) { /* in the socket case, it may happen that the connection failed */
1012       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
1013         PetscMPIInt size,rank;
1014 
1015         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
1016         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1017         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
1018       }
1019       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
1020       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
1021       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
1022         PetscInt    dim;
1023         const char* name;
1024 
1025         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
1026         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1027         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
1028         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
1029       }
1030     }
1031   } else {
1032     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
1033   }
1034   PetscFunctionReturn(0);
1035 }
1036