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