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