xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 11a4995d88c6b984f3abb7b1ba77d48c471e7283)
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);
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;
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;
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 = PetscCalloc5(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&fec_type,totdofs,&idxs);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           if (H1) {
144             nlocal[ctx->nf] = Nc * Nv;
145             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
146             ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr);
147             for (i=0,cum=0;i<vEnd-vStart;i++) {
148               PetscInt j,off;
149 
150               if (PetscUnlikely(!PetscBTLookup(vown,i))) continue;
151               ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr);
152               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
153             }
154             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
155           } else {
156             nlocal[ctx->nf] = Nc * totc;
157             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
158             ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr);
159             for (i=0,cum=0;i<cEnd-cStart;i++) {
160               PetscInt j,off;
161 
162               if (PetscUnlikely(gNum[i] < 0)) continue;
163               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
164               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
165             }
166             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
167           }
168           ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
169           ierr = VecDestroy(&xfield);CHKERRQ(ierr);
170           ierr = ISDestroy(&isfield);CHKERRQ(ierr);
171           ctx->nf++;
172         } else if (id == PETSCFV_CLASSID) {
173           PetscInt c;
174 
175           ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr);
176           ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P0",dim);CHKERRQ(ierr);
177           for (c = 0; c < Nc; c++) {
178             char comp[256];
179             ierr = PetscSNPrintf(comp,256,"%s-Comp%d",name,c);CHKERRQ(ierr);
180             ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr);
181             bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */
182             nlocal[ctx->nf] = totc;
183             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
184             ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr);
185             for (i=0,cum=0;i<cEnd-cStart;i++) {
186               PetscInt off;
187 
188               if (PetscUnlikely(gNum[i])<0) continue;
189               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
190               idxs[cum++] = off + c;
191             }
192             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
193             ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
194             ierr = VecDestroy(&xfield);CHKERRQ(ierr);
195             ierr = ISDestroy(&isfield);CHKERRQ(ierr);
196             ctx->nf++;
197           }
198         } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f);
199       } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f);
200     }
201   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
202   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
203   ierr = VecDestroy(&xlocal);CHKERRQ(ierr);
204   ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
205 
206   /* customize the viewer */
207   ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fieldname,(const char**)fec_type,nlocal,bs,DMPlexSampleGLVisFields_Private,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
208   for (f=0;f<ctx->nf;f++) {
209     ierr = PetscFree(fieldname[f]);CHKERRQ(ierr);
210     ierr = PetscFree(fec_type[f]);CHKERRQ(ierr);
211   }
212   ierr = PetscFree5(fieldname,nlocal,bs,fec_type,idxs);CHKERRQ(ierr);
213   PetscFunctionReturn(0);
214 }
215 
216 typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_UNDEF} MFEM_cid;
217 
218 MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
219                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
220                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF, MFEM_UNDEF},
221                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } };
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "DMPlexGetPointMFEMCellID_Internal"
225 PETSC_STATIC_INLINE PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid)
226 {
227   DMLabel        dlabel;
228   PetscInt       depth,csize;
229   PetscErrorCode ierr;
230 
231   PetscFunctionBegin;
232   ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
233   ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
234   ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr);
235   if (label) {
236     ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr);
237   } else {
238     *mid = 1;
239   }
240   *cid = mfem_table_cid[depth][csize];
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "DMPlexGetPointMFEMVertexIDs_Internal"
246 PETSC_STATIC_INLINE PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscInt *nv, int vids[])
247 {
248   PetscInt       dim,i,q,vStart,vEnd,numPoints,*points = NULL;
249   PetscErrorCode ierr;
250 
251   PetscFunctionBegin;
252   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
253   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
254   ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
255   for (i=0,q=0;i<numPoints*2;i+= 2)
256     if ((points[i] >= vStart) && (points[i] < vEnd))
257       vids[q++] = points[i]-vStart;
258   ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
259   ierr = DMPlexInvertCell(dim,q,vids);CHKERRQ(ierr);
260   *nv  = q;
261   PetscFunctionReturn(0);
262 }
263 
264 /* ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR */
265 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
266 {
267   DM                   cdm;
268   DMLabel              label;
269   PetscSection         coordSection,parentSection;
270   Vec                  coordinates;
271   IS                   globalNum = NULL;
272   const PetscScalar    *array;
273   const PetscInt       *gNum;
274   PetscInt             bf,p,sdim,dim,depth,novl;
275   PetscInt             cStart,cEnd,cEndInterior,fStart,fEnd,fEndInterior,vStart,vEnd;
276   PetscMPIInt          commsize;
277   PetscBool            ovl,isascii,enable_boundary,enable_ncmesh;
278   PetscBT              pown;
279   PetscErrorCode       ierr;
280   PetscInt             *faces = NULL,fpc = 0,vpf = 0;
281   PetscInt             fv1[]     = {0,1},
282                        fv2tri[]  = {0,1,
283                                     1,2,
284                                     2,0},
285                        fv2quad[] = {0,1,
286                                     1,2,
287                                     2,3,
288                                     3,0},
289                        fv3tet[]  = {0,1,2,
290                                     0,3,1,
291                                     0,2,3,
292                                     2,1,3},
293                        fv3hex[]  = {0,1,2,3,
294                                  4,5,6,7,
295                                  0,3,5,4,
296                                  2,1,7,6,
297                                  3,2,6,5,
298                                  0,4,7,1};
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, &fEndInterior, 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 
326   /*
327      a couple of sections of the mesh specification are disabled
328        - boundary: unless we want to visualize boundary attributes, the boundary is not needed for proper mesh visualization
329        - 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
330   */
331   enable_boundary = PETSC_FALSE;
332   enable_ncmesh   = PETSC_FALSE;
333   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
334   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);
335   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);
336   ierr = PetscOptionsEnd();CHKERRQ(ierr);
337 
338   /* Identify possible cells in the overlap */
339   gNum = NULL;
340   novl = 0;
341   ovl  = PETSC_FALSE;
342   pown = NULL;
343   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
344   if (commsize > 1) {
345     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
346     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
347     for (p=cStart; p<cEnd; p++) {
348       if (gNum[p-cStart] < 0) {
349         ovl = PETSC_TRUE;
350         novl++;
351       }
352     }
353     if (ovl) {
354       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
355          TODO: garbage collector? attach pown to dm?  */
356       ierr = PetscBTCreate(PetscMax(cEnd-cStart,vEnd-vStart),&pown);CHKERRQ(ierr);
357     }
358   }
359 
360   if (!enabled) {
361     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
362     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
363     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
364     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
365     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
366     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
367     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
368     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
369     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
370     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
371     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
372     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
373     if (globalNum) {
374       ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
375     }
376     PetscFunctionReturn(0);
377   }
378 
379   /* header */
380   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
381 
382   /* topological dimension */
383   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
384   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
385 
386   /* elements */
387   /* TODO: is this the label we want to use for marking material IDs?
388      We should decide to have a single marker for these stuff
389      Proposal: DMSetMaterialLabel?
390   */
391   ierr = DMGetLabel(dm,"Cell sets",&label);CHKERRQ(ierr);
392   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
393   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
394   for (p=cStart;p<cEnd;p++) {
395     int      vids[8];
396     PetscInt i,nv = 0,cid = -1,mid = 1;
397 
398     if (ovl) {
399       if (gNum[p-cStart] < 0) continue;
400       else {
401         ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
402       }
403     }
404     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
405     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,&nv,vids);CHKERRQ(ierr);
406     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
407     for (i=0;i<nv;i++) {
408       ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
409     }
410     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
411   }
412 
413   /* determine orientation of boundary mesh */
414   if (cEnd-cStart) {
415     ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
416     switch(dim) {
417       case 1:
418         if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
419         faces = fv1;
420         vpf = 1;
421         break;
422       case 2:
423         switch (fpc) {
424           case 3:
425             faces = fv2tri;
426             vpf   = 2;
427             break;
428           case 4:
429             faces = fv2quad;
430             vpf   = 2;
431             break;
432           default:
433             SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
434             break;
435         }
436         break;
437       case 3:
438         switch (fpc) {
439           case 4:
440             faces = fv3tet;
441             vpf   = 3;
442             break;
443           case 6:
444             faces = fv3hex;
445             vpf   = 4;
446             break;
447           default:
448             SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
449             break;
450         }
451         break;
452       default:
453         SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
454         break;
455     }
456   }
457 
458   /* boundary */
459   ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
460   fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
461   if (!ovl) {
462     for (p=fStart,bf=0;p<fEnd;p++) {
463       PetscInt supportSize;
464 
465       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
466       if (supportSize == 1) bf++;
467     }
468   } else {
469     for (p=fStart,bf=0;p<fEnd;p++) {
470       const PetscInt *support;
471       PetscInt       i,supportSize;
472       PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
473 
474       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
475       ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
476       for (i=0;i<supportSize;i++) {
477         if (PetscBTLookup(pown,support[i]-cStart)) has_owned = PETSC_TRUE;
478         else has_ghost = PETSC_TRUE;
479       }
480       if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) bf++;
481     }
482   }
483   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
484   if (!enable_boundary) {
485     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
486   } else {
487     /* TODO: is this the label we want to use for marking boundary facets?
488        We should decide to have a single marker for these stuff
489        Proposal: DMSetBoundaryLabel?
490     */
491     ierr = DMGetLabel(dm,"marker",&label);CHKERRQ(ierr);
492     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
493     if (!ovl) {
494       for (p=fStart;p<fEnd;p++) {
495         PetscInt supportSize;
496 
497         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
498         if (supportSize == 1) {
499           const    PetscInt *support, *cone;
500           PetscInt i,c,v,cid = -1,mid = -1;
501           PetscInt cellClosureSize,*cellClosure = NULL;
502 
503           ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
504           ierr = DMPlexGetCone(dm,support[0],&cone);CHKERRQ(ierr);
505           for (c=0;c<fpc;c++)
506             if (cone[c] == p)
507               break;
508 
509           ierr = DMPlexGetTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
510           for (v=0;v<cellClosureSize;v++)
511             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
512               break;
513           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
514           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
515           for (i=0;i<vpf;i++) {
516             ierr = PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);CHKERRQ(ierr);
517           }
518           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
519           ierr = DMPlexRestoreTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
520         }
521       }
522     } else {
523       for (p=fStart;p<fEnd;p++) {
524         const PetscInt *support;
525         PetscInt       i,supportSize,validcell = -1;
526         PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
527 
528         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
529         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
530         for (i=0;i<supportSize;i++) {
531           if (PetscBTLookup(pown,support[i]-cStart)) {
532             has_owned = PETSC_TRUE;
533             validcell = support[i];
534           } else has_ghost = PETSC_TRUE;
535         }
536         if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) {
537           const    PetscInt *support, *cone;
538           PetscInt i,c,v,cid = -1,mid = -1;
539           PetscInt cellClosureSize,*cellClosure = NULL;
540 
541           ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
542           ierr = DMPlexGetCone(dm,validcell,&cone);CHKERRQ(ierr);
543           for (c=0;c<fpc;c++)
544             if (cone[c] == p)
545               break;
546 
547           ierr = DMPlexGetTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
548           for (v=0;v<cellClosureSize;v++)
549             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
550               break;
551           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
552           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
553           for (i=0;i<vpf;i++) {
554             ierr = PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);CHKERRQ(ierr);
555           }
556           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
557           ierr = DMPlexRestoreTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
558         }
559       }
560     }
561   }
562 
563   /* mark owned vertices */
564   if (ovl) {
565     ierr = PetscBTMemzero(vEnd-vStart,pown);CHKERRQ(ierr);
566     for (p=cStart;p<cEnd;p++) {
567       PetscInt i,closureSize,*closure = NULL;
568 
569       if (gNum[p-cStart] < 0) continue;
570       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
571       for (i=0;i<closureSize;i++) {
572         const PetscInt pp = closure[2*i];
573 
574         if (pp >= vStart && pp < vEnd) {
575           ierr = PetscBTSet(pown,pp-vStart);CHKERRQ(ierr);
576         }
577       }
578       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
579     }
580   }
581   if (globalNum) {
582     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
583   }
584 
585   /* vertex_parents (Non-conforming meshes) */
586   parentSection  = NULL;
587   if (enable_ncmesh) {
588     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
589   }
590   if (parentSection) {
591     PetscInt vp,gvp;
592 
593     for (vp=0,p=vStart;p<vEnd;p++) {
594       DMLabel  dlabel;
595       PetscInt parent,depth;
596 
597       if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
598       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
599       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
600       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
601       if (parent != p) vp++;
602     }
603     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
604     if (gvp) {
605       PetscInt  maxsupp;
606       PetscBool *skip = NULL;
607 
608       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
609       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
610       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
611       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
612       for (p=vStart;p<vEnd;p++) {
613         DMLabel  dlabel;
614         PetscInt parent;
615 
616         if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
617         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
618         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
619         if (parent != p) {
620           int            vids[8];
621           PetscInt       i,nv,size,n,numChildren,depth = -1;
622           const PetscInt *children;
623           ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr);
624           switch (size) {
625             case 2: /* edge */
626               nv   = 0;
627               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,&nv,vids);CHKERRQ(ierr);
628               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
629               for (i=0;i<nv;i++) {
630                 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
631               }
632               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
633               vp--;
634               break;
635             case 4: /* face */
636               ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
637               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
638               for (n=0;n<numChildren;n++) {
639                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
640                 if (!depth) {
641                   const PetscInt *hvsupp,*hesupp,*cone;
642                   PetscInt       hvsuppSize,hesuppSize,coneSize;
643                   PetscInt       hv = children[n],he,f;
644                   PetscBool      found = PETSC_FALSE;
645 
646                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
647                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
648                   for (i=0;i<hvsuppSize;i++) {
649                     PetscInt ep;
650                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
651                     if (ep != hvsupp[i]) {
652                       he = hvsupp[i];
653                       found = PETSC_TRUE;
654                     } else {
655                       skip[i] = PETSC_TRUE;
656                     }
657                   }
658                   if (!found) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
659                   he      = hvsupp[i];
660                   skip[i] = PETSC_TRUE;
661                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
662                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
663                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
664                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
665                   for (f=0;f<hesuppSize;f++) {
666                     PetscInt j;
667 
668                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
669                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
670                     for (j=0;j<coneSize;j++) {
671                       PetscInt k;
672                       for (k=0;k<hvsuppSize;k++) {
673                         if (hvsupp[k] == cone[j]) {
674                           skip[k] = PETSC_TRUE;
675                           break;
676                         }
677                       }
678                     }
679                   }
680                   for (i=0;i<hvsuppSize;i++) {
681                     if (!skip[i]) {
682                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
683                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
684                     }
685                   }
686                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
687                   for (i=0;i<2;i++) {
688                     ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);CHKERRQ(ierr);
689                   }
690                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
691                   vp--;
692                 }
693               }
694               break;
695             default:
696               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size);
697           }
698         }
699       }
700       ierr = PetscFree(skip);CHKERRQ(ierr);
701     }
702     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
703   }
704   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
705 
706   /* vertices */
707   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
708   ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr);
709   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
710   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
711   ierr = DMPlexGetDepthStratum(cdm,0,&vStart,&vEnd);CHKERRQ(ierr);
712   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
713   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
714   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
715   ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
716   for (p=vStart;p<vEnd;p++) {
717     PetscInt i,st,dof;
718 
719     ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
720     ierr = PetscSectionGetOffset(coordSection,p,&st);CHKERRQ(ierr);
721     for (i=st;i<st+dof-1;i++) {
722       ierr = PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[i]));CHKERRQ(ierr);
723     }
724     ierr = PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[i]));CHKERRQ(ierr);
725   }
726   ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
727   PetscFunctionReturn(0);
728 }
729 
730 /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
731 PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
732 {
733   PetscErrorCode ierr;
734   PetscBool      isglvis,isascii;
735 
736   PetscFunctionBegin;
737   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
738   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
739   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
740   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
741   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
742   if (isglvis) {
743     PetscViewer          view;
744     PetscViewerGLVisType type;
745 
746     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
747     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
748     if (view) { /* in the socket case, it may happen that the connection failed */
749       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
750         PetscMPIInt size,rank;
751         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
752         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
753         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
754       }
755       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
756       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
757       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
758         PetscInt    dim;
759         const char* name;
760 
761         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
762         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
763         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
764         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
765       }
766     }
767   } else {
768     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
769   }
770   PetscFunctionReturn(0);
771 }
772