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