xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 19e8fbdabb07dcd9104f2d10b4dec626a2e1a063)
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 (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 {
247     *mid = 1;
248   }
249   *cid = mfem_table_cid[depth][csize];
250   PetscFunctionReturn(0);
251 }
252 
253 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[])
254 {
255   PetscInt       dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
260   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
261   sdim = dim;
262   if (csec) {
263     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
264     ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr);
265     off  = off/sdim;
266     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
267   }
268   if (!dof) {
269     ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
270     for (i=0,q=0;i<numPoints*2;i+= 2)
271       if ((points[i] >= vStart) && (points[i] < vEnd))
272         vids[q++] = points[i]-vStart+off;
273     ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
274   } else {
275     ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr);
276     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
277     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
278   }
279   ierr = DMPlexInvertCell(dim,q,vids);CHKERRQ(ierr);
280   *nv  = q;
281   PetscFunctionReturn(0);
282 }
283 
284 /* ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR */
285 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
286 {
287   DMLabel              label;
288   PetscSection         coordSection,parentSection;
289   Vec                  coordinates;
290   IS                   globalNum = NULL;
291   const PetscScalar    *array;
292   const PetscInt       *gNum;
293   PetscInt             bf,p,sdim,dim,depth,novl;
294   PetscInt             cStart,cEnd,cEndInterior,vStart,vEnd,nvert;
295   PetscMPIInt          commsize;
296   PetscBool            localized,ovl,isascii,enable_boundary,enable_ncmesh;
297   PetscBT              pown;
298   PetscErrorCode       ierr;
299   PetscContainer       glvis_container;
300   PetscBool            enabled = PETSC_TRUE;
301 
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
304   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
305   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
306   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
307   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
308   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
309   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
310   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
311   if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
312 
313   /* get container: determines if a process visualizes is portion of the data or not */
314   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
315   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
316   {
317     PetscViewerGLVisInfo glvis_info;
318     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
319     enabled = glvis_info->enabled;
320   }
321   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
322   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
323   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
324   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
325   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
326   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
327 
328   /*
329      a couple of sections of the mesh specification are disabled
330        - boundary: unless we want to visualize boundary attributes, the boundary is not needed for proper mesh visualization
331        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package and be able to derefine the mesh
332   */
333   enable_boundary = PETSC_FALSE;
334   enable_ncmesh   = PETSC_FALSE;
335   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
336   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation; useful for debugging purposes",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr);
337   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation; useful for debugging purposes",NULL,enable_ncmesh,&enable_ncmesh,NULL);CHKERRQ(ierr);
338   ierr = PetscOptionsEnd();CHKERRQ(ierr);
339 
340   /* Identify possible cells in the overlap */
341   gNum = NULL;
342   novl = 0;
343   ovl  = PETSC_FALSE;
344   pown = NULL;
345   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
346   if (commsize > 1) {
347     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
348     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
349     for (p=cStart; p<cEnd; p++) {
350       if (gNum[p-cStart] < 0) {
351         ovl = PETSC_TRUE;
352         novl++;
353       }
354     }
355     if (ovl) {
356       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
357          TODO: garbage collector? attach pown to dm?  */
358       ierr = PetscBTCreate(PetscMax(cEnd-cStart,vEnd-vStart),&pown);CHKERRQ(ierr);
359     }
360   }
361 
362   if (!enabled) {
363     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
364     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
365     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
366     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
367     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
368     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
369     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
370     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
371     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
372     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
373     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
374     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
375     if (globalNum) {
376       ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
377     }
378     PetscFunctionReturn(0);
379   }
380 
381   /* header */
382   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
383 
384   /* topological dimension */
385   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
386   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
387 
388   /* elements */
389   /* TODO: is this the label we want to use for marking material IDs?
390      We should decide to have a single marker for these stuff
391      Proposal: DMSetMaterialLabel?
392   */
393   ierr = DMGetLabel(dm,"Cell Sets",&label);CHKERRQ(ierr);
394   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
395   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
396   for (p=cStart;p<cEnd;p++) {
397     int      vids[8];
398     PetscInt i,nv = 0,cid = -1,mid = 1;
399 
400     if (ovl) {
401       if (gNum[p-cStart] < 0) continue;
402       else {
403         ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
404       }
405     }
406     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
407     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
408     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
409     for (i=0;i<nv;i++) {
410       ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
411     }
412     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
413   }
414 
415   /* boundary */
416   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
417   if (!enable_boundary) {
418     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
419   } else {
420     PetscInt fStart,fEnd,fEndInterior;
421     PetscInt *faces = NULL,fpc = 0,vpf = 0;
422     PetscInt fv1[]     = {0,1},
423              fv2tri[]  = {0,1,
424                           1,2,
425                           2,0},
426              fv2quad[] = {0,1,
427                           1,2,
428                           2,3,
429                           3,0},
430              fv3tet[]  = {0,1,2,
431                           0,3,1,
432                           0,2,3,
433                           2,1,3},
434              fv3hex[]  = {0,1,2,3,
435                        4,5,6,7,
436                        0,3,5,4,
437                        2,1,7,6,
438                        3,2,6,5,
439                        0,4,7,1};
440 
441     if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Boundary specification with periodic mesh");
442     /* determine orientation of boundary mesh */
443     if (cEnd-cStart) {
444       ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
445       switch(dim) {
446         case 1:
447           if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
448           faces = fv1;
449           vpf = 1;
450           break;
451         case 2:
452           switch (fpc) {
453             case 3:
454               faces = fv2tri;
455               vpf   = 2;
456               break;
457             case 4:
458               faces = fv2quad;
459               vpf   = 2;
460               break;
461             default:
462               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
463               break;
464           }
465           break;
466         case 3:
467           switch (fpc) {
468             case 4:
469               faces = fv3tet;
470               vpf   = 3;
471               break;
472             case 6:
473               faces = fv3hex;
474               vpf   = 4;
475               break;
476             default:
477               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
478               break;
479           }
480           break;
481         default:
482           SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
483           break;
484       }
485     }
486 
487     ierr = DMPlexGetHybridBounds(dm, NULL, &fEndInterior, NULL, NULL);CHKERRQ(ierr);
488     ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
489     fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
490     if (!ovl) {
491       for (p=fStart,bf=0;p<fEnd;p++) {
492         PetscInt supportSize;
493 
494         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
495         if (supportSize == 1) bf++;
496       }
497     } else {
498       for (p=fStart,bf=0;p<fEnd;p++) {
499         const PetscInt *support;
500         PetscInt       i,supportSize;
501         PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
502 
503         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
504         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
505         for (i=0;i<supportSize;i++) {
506           if (PetscBTLookup(pown,support[i]-cStart)) has_owned = PETSC_TRUE;
507           else has_ghost = PETSC_TRUE;
508         }
509         if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) bf++;
510       }
511     }
512     /* TODO: is this the label we want to use for marking boundary facets?
513        We should decide to have a single marker for these stuff
514        Proposal: DMSetBoundaryLabel?
515     */
516     ierr = DMGetLabel(dm,"marker",&label);CHKERRQ(ierr);
517     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
518     if (!ovl) {
519       for (p=fStart;p<fEnd;p++) {
520         PetscInt supportSize;
521 
522         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
523         if (supportSize == 1) {
524           const    PetscInt *support, *cone;
525           PetscInt i,c,v,cid = -1,mid = -1;
526           PetscInt cellClosureSize,*cellClosure = NULL;
527 
528           ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
529           ierr = DMPlexGetCone(dm,support[0],&cone);CHKERRQ(ierr);
530           for (c=0;c<fpc;c++)
531             if (cone[c] == p)
532               break;
533 
534           ierr = DMPlexGetTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
535           for (v=0;v<cellClosureSize;v++)
536             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
537               break;
538           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
539           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
540           for (i=0;i<vpf;i++) {
541             ierr = PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);CHKERRQ(ierr);
542           }
543           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
544           ierr = DMPlexRestoreTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
545         }
546       }
547     } else {
548       for (p=fStart;p<fEnd;p++) {
549         const PetscInt *support;
550         PetscInt       i,supportSize,validcell = -1;
551         PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
552 
553         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
554         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
555         for (i=0;i<supportSize;i++) {
556           if (PetscBTLookup(pown,support[i]-cStart)) {
557             has_owned = PETSC_TRUE;
558             validcell = support[i];
559           } else has_ghost = PETSC_TRUE;
560         }
561         if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) {
562           const    PetscInt *support, *cone;
563           PetscInt i,c,v,cid = -1,mid = -1;
564           PetscInt cellClosureSize,*cellClosure = NULL;
565 
566           ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
567           ierr = DMPlexGetCone(dm,validcell,&cone);CHKERRQ(ierr);
568           for (c=0;c<fpc;c++)
569             if (cone[c] == p)
570               break;
571 
572           ierr = DMPlexGetTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
573           for (v=0;v<cellClosureSize;v++)
574             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
575               break;
576           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
577           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
578           for (i=0;i<vpf;i++) {
579             ierr = PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);CHKERRQ(ierr);
580           }
581           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
582           ierr = DMPlexRestoreTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
583         }
584       }
585     }
586   }
587 
588   /* mark owned vertices */
589   if (ovl) {
590     ierr = PetscBTMemzero(vEnd-vStart,pown);CHKERRQ(ierr);
591     for (p=cStart;p<cEnd;p++) {
592       PetscInt i,closureSize,*closure = NULL;
593 
594       if (gNum[p-cStart] < 0) continue;
595       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
596       for (i=0;i<closureSize;i++) {
597         const PetscInt pp = closure[2*i];
598 
599         if (pp >= vStart && pp < vEnd) {
600           ierr = PetscBTSet(pown,pp-vStart);CHKERRQ(ierr);
601         }
602       }
603       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
604     }
605   }
606   if (globalNum) {
607     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
608   }
609 
610   /* vertex_parents (Non-conforming meshes) */
611   parentSection  = NULL;
612   if (enable_ncmesh) {
613     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
614   }
615   if (parentSection) {
616     PetscInt vp,gvp;
617 
618     for (vp=0,p=vStart;p<vEnd;p++) {
619       DMLabel  dlabel;
620       PetscInt parent,depth;
621 
622       if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
623       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
624       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
625       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
626       if (parent != p) vp++;
627     }
628     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
629     if (gvp) {
630       PetscInt  maxsupp;
631       PetscBool *skip = NULL;
632 
633       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
634       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
635       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
636       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
637       for (p=vStart;p<vEnd;p++) {
638         DMLabel  dlabel;
639         PetscInt parent;
640 
641         if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
642         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
643         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
644         if (parent != p) {
645           int            vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
646           PetscInt       i,nv,size,n,numChildren,depth = -1;
647           const PetscInt *children;
648           ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr);
649           switch (size) {
650             case 2: /* edge */
651               nv   = 0;
652               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
653               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
654               for (i=0;i<nv;i++) {
655                 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
656               }
657               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
658               vp--;
659               break;
660             case 4: /* face */
661               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
662               for (n=0;n<numChildren;n++) {
663                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
664                 if (!depth) {
665                   const PetscInt *hvsupp,*hesupp,*cone;
666                   PetscInt       hvsuppSize,hesuppSize,coneSize;
667                   PetscInt       hv = children[n],he = -1,f;
668 
669                   ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
670                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
671                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
672                   for (i=0;i<hvsuppSize;i++) {
673                     PetscInt ep;
674                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
675                     if (ep != hvsupp[i]) {
676                       he = hvsupp[i];
677                     } else {
678                       skip[i] = PETSC_TRUE;
679                     }
680                   }
681                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
682                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
683                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
684                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
685                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
686                   for (f=0;f<hesuppSize;f++) {
687                     PetscInt j;
688 
689                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
690                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
691                     for (j=0;j<coneSize;j++) {
692                       PetscInt k;
693                       for (k=0;k<hvsuppSize;k++) {
694                         if (hvsupp[k] == cone[j]) {
695                           skip[k] = PETSC_TRUE;
696                           break;
697                         }
698                       }
699                     }
700                   }
701                   for (i=0;i<hvsuppSize;i++) {
702                     if (!skip[i]) {
703                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
704                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
705                     }
706                   }
707                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
708                   for (i=0;i<2;i++) {
709                     ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);CHKERRQ(ierr);
710                   }
711                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
712                   vp--;
713                 }
714               }
715               break;
716             default:
717               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size);
718           }
719         }
720       }
721       ierr = PetscFree(skip);CHKERRQ(ierr);
722     }
723     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
724   }
725   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
726 
727   /* vertices */
728   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
729   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
730   ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr);
731   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
732   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr);
733   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
734   ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
735   for (p=0;p<nvert/sdim;p++) {
736     PetscInt s;
737     for (s=0;s<sdim;s++) {
738       ierr = PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr);
739     }
740     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
741   }
742   ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
743   PetscFunctionReturn(0);
744 }
745 
746 /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
747 PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
748 {
749   PetscErrorCode ierr;
750   PetscBool      isglvis,isascii;
751 
752   PetscFunctionBegin;
753   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
754   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
755   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
756   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
757   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
758   if (isglvis) {
759     PetscViewer          view;
760     PetscViewerGLVisType type;
761 
762     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
763     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
764     if (view) { /* in the socket case, it may happen that the connection failed */
765       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
766         PetscMPIInt size,rank;
767         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
768         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
769         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
770       }
771       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
772       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
773       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
774         PetscInt    dim;
775         const char* name;
776 
777         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
778         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
779         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
780         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
781       }
782     }
783   } else {
784     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
785   }
786   PetscFunctionReturn(0);
787 }
788