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