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