xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 96d2aba5c40a07bc26ef1954a3f1cfdaabadb7a8)
1 #include <petsc/private/glvisviewerimpl.h>
2 #include <petsc/private/petscimpl.h>
3 #include <petsc/private/dmpleximpl.h>
4 #include <petscbt.h>
5 #include <petscdmplex.h>
6 #include <petscsf.h>
7 #include <petscds.h>
8 
9 typedef struct {
10   PetscInt   nf;
11   VecScatter *scctx;
12 } GLVisViewerCtx;
13 
14 static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
15 {
16   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
17   PetscInt       i;
18   PetscErrorCode ierr;
19 
20   PetscFunctionBegin;
21   for (i=0;i<ctx->nf;i++) {
22     ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr);
23   }
24   ierr = PetscFree(ctx->scctx);CHKERRQ(ierr);
25   ierr = PetscFree(vctx);CHKERRQ(ierr);
26   PetscFunctionReturn(0);
27 }
28 
29 static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
30 {
31   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
32   PetscInt       f;
33   PetscErrorCode ierr;
34 
35   PetscFunctionBegin;
36   for (f=0;f<nf;f++) {
37     ierr = VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
38     ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39   }
40   PetscFunctionReturn(0);
41 }
42 
43 /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */
44 PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer)
45 {
46   DM             dm = (DM)odm;
47   Vec            xlocal,xfield,*Ufield;
48   PetscDS        ds;
49   IS             globalNum,isfield;
50   PetscBT        vown;
51   char           **fieldname = NULL,**fec_type = NULL;
52   const PetscInt *gNum;
53   PetscInt       *nlocal,*bs,*idxs,*dims;
54   PetscInt       f,maxfields,nfields,c,totc,totdofs,Nv,cum,i;
55   PetscInt       dim,cStart,cEnd,cEndInterior,vStart,vEnd;
56   GLVisViewerCtx *ctx;
57   PetscSection   s;
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
62   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
63   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
64   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
65   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
66   ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
67   ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
68   ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
69   for (c = cStart, totc = 0; c < cEnd; c++) {
70     if (gNum[c-cStart] >= 0) {
71       PetscInt i,numPoints,*points = NULL;
72 
73       totc++;
74       ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
75       for (i=0;i<numPoints*2;i+= 2) {
76         if ((points[i] >= vStart) && (points[i] < vEnd)) {
77           ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr);
78         }
79       }
80       ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
81     }
82   }
83   for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++;
84 
85   ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr);
86   ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr);
87   ierr = DMGetDefaultSection(dm,&s);CHKERRQ(ierr);
88   ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr);
89   for (f=0,maxfields=0;f<nfields;f++) {
90     PetscInt bs;
91 
92     ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr);
93     maxfields += bs;
94   }
95   ierr = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);CHKERRQ(ierr);
96   ierr = PetscNew(&ctx);CHKERRQ(ierr);
97   ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr);
98   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
99   if (ds) {
100     for (f=0;f<nfields;f++) {
101       const char* fname;
102       char        name[256];
103       PetscObject disc;
104       size_t      len;
105 
106       ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr);
107       ierr = PetscStrlen(fname,&len);CHKERRQ(ierr);
108       if (len) {
109         ierr = PetscStrcpy(name,fname);CHKERRQ(ierr);
110       } else {
111         ierr = PetscSNPrintf(name,256,"Field%d",f);CHKERRQ(ierr);
112       }
113       ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
114       if (disc) {
115         PetscClassId id;
116         PetscInt     Nc;
117         char         fec[64];
118 
119         ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
120         if (id == PETSCFE_CLASSID) {
121           PetscFE            fem = (PetscFE)disc;
122           PetscDualSpace     sp;
123           PetscDualSpaceType spname;
124           PetscInt           order;
125           PetscBool          islag,continuous,H1 = PETSC_TRUE;
126 
127           ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr);
128           ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr);
129           ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr);
130           ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr);
131           if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
132           ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr);
133           if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported");
134           ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
135           if (continuous && order > 0) {
136             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%dD_P1",dim);CHKERRQ(ierr);
137           } else {
138             H1   = PETSC_FALSE;
139             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P%d",dim,order);CHKERRQ(ierr);
140           }
141           ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr);
142           bs[ctx->nf]   = Nc;
143           dims[ctx->nf] = dim;
144           if (H1) {
145             nlocal[ctx->nf] = Nc * Nv;
146             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
147             ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr);
148             for (i=0,cum=0;i<vEnd-vStart;i++) {
149               PetscInt j,off;
150 
151               if (PetscUnlikely(!PetscBTLookup(vown,i))) continue;
152               ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr);
153               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
154             }
155             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
156           } else {
157             nlocal[ctx->nf] = Nc * totc;
158             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
159             ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr);
160             for (i=0,cum=0;i<cEnd-cStart;i++) {
161               PetscInt j,off;
162 
163               if (PetscUnlikely(gNum[i] < 0)) continue;
164               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
165               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
166             }
167             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
168           }
169           ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
170           ierr = VecDestroy(&xfield);CHKERRQ(ierr);
171           ierr = ISDestroy(&isfield);CHKERRQ(ierr);
172           ctx->nf++;
173         } else if (id == PETSCFV_CLASSID) {
174           PetscInt c;
175 
176           ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr);
177           ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P0",dim);CHKERRQ(ierr);
178           for (c = 0; c < Nc; c++) {
179             char comp[256];
180             ierr = PetscSNPrintf(comp,256,"%s-Comp%d",name,c);CHKERRQ(ierr);
181             ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr);
182             bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */
183             nlocal[ctx->nf] = totc;
184             dims[ctx->nf] = dim;
185             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
186             ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr);
187             for (i=0,cum=0;i<cEnd-cStart;i++) {
188               PetscInt off;
189 
190               if (PetscUnlikely(gNum[i])<0) continue;
191               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
192               idxs[cum++] = off + c;
193             }
194             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
195             ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
196             ierr = VecDestroy(&xfield);CHKERRQ(ierr);
197             ierr = ISDestroy(&isfield);CHKERRQ(ierr);
198             ctx->nf++;
199           }
200         } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f);
201       } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f);
202     }
203   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
204   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
205   ierr = VecDestroy(&xlocal);CHKERRQ(ierr);
206   ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
207 
208   /* create work vectors */
209   for (f=0;f<ctx->nf;f++) {
210     ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr);
211     ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr);
212     ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr);
213     ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr);
214   }
215 
216   /* customize the viewer */
217   ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
218   for (f=0;f<ctx->nf;f++) {
219     ierr = PetscFree(fieldname[f]);CHKERRQ(ierr);
220     ierr = PetscFree(fec_type[f]);CHKERRQ(ierr);
221     ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr);
222   }
223   ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr);
224   PetscFunctionReturn(0);
225 }
226 
227 typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_UNDEF} MFEM_cid;
228 
229 MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
230                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
231                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF, MFEM_UNDEF},
232                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } };
233 
234 static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid)
235 {
236   DMLabel        dlabel;
237   PetscInt       depth,csize;
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
242   ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
243   ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr);
244   if (label) {
245     ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr);
246   } else *mid = 1;
247   *cid = mfem_table_cid[depth][csize];
248   PetscFunctionReturn(0);
249 }
250 
251 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[])
252 {
253   PetscInt       dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;
254   PetscErrorCode ierr;
255 
256   PetscFunctionBegin;
257   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
258   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
259   sdim = dim;
260   if (csec) {
261     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
262     ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr);
263     off  = off/sdim;
264     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
265   }
266   if (!dof) {
267     ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
268     for (i=0,q=0;i<numPoints*2;i+= 2)
269       if ((points[i] >= vStart) && (points[i] < vEnd))
270         vids[q++] = points[i]-vStart+off;
271     ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
272   } else {
273     ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr);
274     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
275     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
276   }
277   *nv  = q;
278   PetscFunctionReturn(0);
279 }
280 
281 /*
282    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
283    Higher order meshes are also supported
284 */
285 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
286 {
287   DMLabel              label;
288   PetscSection         coordSection,parentSection;
289   Vec                  coordinates,hovec;
290   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            periodic, enabled = PETSC_TRUE;
301   const char           *fmt;
302 
303   PetscFunctionBegin;
304   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
305   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
306   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
307   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
308   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
309   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
310   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
311   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
312   if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
313 
314   /* get container: determines if a process visualizes is portion of the data or not */
315   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
316   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
317   {
318     PetscViewerGLVisInfo glvis_info;
319     ierr    = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
320     enabled = glvis_info->enabled;
321     fmt     = glvis_info->fmt;
322   }
323   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
324   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
325   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
326   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
327   ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr);
328   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
329   if (periodic && !localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
330   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
331   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
332   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
333 
334   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
335      DMPlex does not currently support HO meshes, so there's no API for this */
336   ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr);
337 
338   /*
339      a couple of sections of the mesh specification are disabled
340        - boundary: unless we want to visualize boundary attributes or we have a periodic mesh
341                    the boundary is not needed for proper mesh visualization
342        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
343                          and be able to derefine the mesh
344   */
345   enable_boundary = periodic;
346   enable_ncmesh = PETSC_FALSE;
347   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
348   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);
349   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);
350   ierr = PetscOptionsEnd();CHKERRQ(ierr);
351 
352   /* Identify possible cells in the overlap */
353   gNum = NULL;
354   novl = 0;
355   ovl  = PETSC_FALSE;
356   pown = NULL;
357   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
358   if (commsize > 1) {
359     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
360     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
361     for (p=cStart; p<cEnd; p++) {
362       if (gNum[p-cStart] < 0) {
363         ovl = PETSC_TRUE;
364         novl++;
365       }
366     }
367     if (ovl) {
368       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
369          TODO: garbage collector? attach pown to dm?  */
370       ierr = PetscBTCreate(PetscMax(cEnd-cStart,vEnd-vStart),&pown);CHKERRQ(ierr);
371     }
372   }
373 
374   if (!enabled) {
375     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
376     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
377     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
378     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
379     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
380     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
381     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
382     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
383     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
384     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
385     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
386     if (globalNum) {
387       ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
388     }
389     PetscFunctionReturn(0);
390   }
391 
392   /* header */
393   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
394 
395   /* topological dimension */
396   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
397   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
398 
399   /* elements */
400   /* TODO: is this the label we want to use for marking material IDs?
401      We should decide to have a single marker for these stuff
402      Proposal: DMSetMaterialLabel?
403   */
404   ierr = DMGetLabel(dm,"Cell Sets",&label);CHKERRQ(ierr);
405   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
406   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
407   for (p=cStart;p<cEnd;p++) {
408     int      vids[8];
409     PetscInt i,nv = 0,cid = -1,mid = 1;
410 
411     if (ovl) {
412       if (gNum[p-cStart] < 0) continue;
413       else {
414         ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
415       }
416     }
417     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
418     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
419     ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
420     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
421     for (i=0;i<nv;i++) {
422       ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
423     }
424     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
425   }
426 
427   /* boundary */
428   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
429   if (!enable_boundary) {
430     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
431   } else {
432     DMLabel  perLabel;
433     PetscBT  bfaces;
434     PetscInt fStart,fEnd,fEndInterior,*fcells;
435     PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0;
436     PetscInt fv1[]     = {0,1},
437              fv2tri[]  = {0,1,
438                           1,2,
439                           2,0},
440              fv2quad[] = {0,1,
441                           1,2,
442                           2,3,
443                           3,0},
444              fv3tet[]  = {0,1,2,
445                           0,3,1,
446                           0,2,3,
447                           2,1,3},
448              fv3hex[]  = {0,1,2,3,
449                        4,5,6,7,
450                        0,3,5,4,
451                        2,1,7,6,
452                        3,2,6,5,
453                        0,4,7,1};
454 
455     /* determine orientation of boundary mesh */
456     if (cEnd-cStart) {
457       ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
458       switch(dim) {
459         case 1:
460           if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
461           faces = fv1;
462           vpf = 1;
463           vpc = 2;
464           break;
465         case 2:
466           switch (fpc) {
467             case 3:
468               faces = fv2tri;
469               vpf   = 2;
470               vpc   = 3;
471               break;
472             case 4:
473               faces = fv2quad;
474               vpf   = 2;
475               vpc   = 4;
476               break;
477             default:
478               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
479               break;
480           }
481           break;
482         case 3:
483           switch (fpc) {
484             case 4:
485               faces = fv3tet;
486               vpf   = 3;
487               vpc   = 4;
488               break;
489             case 6:
490               faces = fv3hex;
491               vpf   = 4;
492               vpc   = 8;
493               break;
494             default:
495               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
496               break;
497           }
498           break;
499         default:
500           SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
501           break;
502       }
503     }
504     ierr = DMPlexGetHybridBounds(dm,NULL,&fEndInterior,NULL,NULL);CHKERRQ(ierr);
505     ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
506     fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
507     ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr);
508     ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr);
509     ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr);
510     ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
511     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
512       ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr);
513       ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
514       ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr);
515       for (p=cStart;p<cEnd;p++) {
516         PetscInt dof;
517         ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
518         if (dof) {
519           PetscInt    v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
520           PetscScalar *vals = NULL;
521 
522           if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
523           if (dof/sdim != vpc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible number of cell dofs %D and vertices %D",dof/sdim,vpc);
524           ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
525           ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
526           for (v=0;v<cellClosureSize;v++)
527             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
528               vidxs = cellClosure + 2*v;
529               break;
530             }
531           if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
532           for (v=0;v<vpc;v++) {
533             PetscInt s;
534             for (s=0;s<sdim;s++) {
535               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+vpc*sdim])>PETSC_MACHINE_EPSILON) {
536                 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr);
537               }
538             }
539           }
540           ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
541           ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
542         }
543       }
544       if (dim > 1) {
545         PetscInt eEnd,eStart,eEndInterior;
546         ierr = DMPlexGetHybridBounds(dm,NULL,NULL,&eEndInterior,NULL);CHKERRQ(ierr);
547         ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr);
548         eEnd = eEndInterior < 0 ? eEnd : eEndInterior;
549         for (p=eStart;p<eEnd;p++) {
550           const PetscInt *cone;
551           PetscInt       coneSize,i;
552           PetscBool      ispe = PETSC_TRUE;
553 
554           ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
555           ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
556           for (i=0;i<coneSize;i++) {
557             PetscInt v;
558 
559             ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
560             ispe = (PetscBool)(ispe && (v==2));
561           }
562           if (ispe && coneSize) {
563             ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
564           }
565         }
566         if (dim > 2) {
567           for (p=fStart;p<fEnd;p++) {
568             const PetscInt *cone;
569             PetscInt       coneSize,i;
570             PetscBool      ispe = PETSC_TRUE;
571 
572             ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
573             ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
574             for (i=0;i<coneSize;i++) {
575               PetscInt v;
576 
577               ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
578               ispe = (PetscBool)(ispe && (v==2));
579             }
580             if (ispe && coneSize) {
581               ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
582             }
583           }
584         }
585       }
586     }
587     for (p=fStart;p<fEnd;p++) {
588       const PetscInt *support;
589       PetscInt       supportSize;
590       PetscBool      isbf = PETSC_FALSE;
591 
592       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
593       if (ovl) {
594         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
595         PetscInt  i;
596 
597         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
598         for (i=0;i<supportSize;i++) {
599           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
600           else has_ghost = PETSC_TRUE;
601         }
602         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
603       } else {
604         isbf = (PetscBool)(supportSize == 1);
605       }
606       if (!isbf && perLabel) {
607         const PetscInt *cone;
608         PetscInt       coneSize,i;
609 
610         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
611         ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
612         isbf = PETSC_TRUE;
613         for (i=0;i<coneSize;i++) {
614           PetscInt v,d;
615 
616           ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
617           ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr);
618           isbf = (PetscBool)(isbf && v != d);
619         }
620       }
621       if (isbf) {
622         ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr);
623       }
624     }
625     /* count boundary faces */
626     for (p=fStart,bf=0;p<fEnd;p++) {
627       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
628         const PetscInt *support;
629         PetscInt       supportSize,c;
630 
631         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
632         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
633         if (ovl) {
634           for (c=0;c<supportSize;c++) {
635             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
636               bf++;
637             }
638           }
639         } else bf += supportSize;
640       }
641     }
642     /* TODO: is this the label we want to use for marking boundary facets?
643        We should decide to have a single marker for these stuff
644        Proposal: DMSetBoundaryLabel?
645     */
646     ierr = DMGetLabel(dm,"marker",&label);CHKERRQ(ierr);
647     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
648     for (p=fStart;p<fEnd;p++) {
649       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
650         const PetscInt *support;
651         PetscInt       supportSize,c,nc = 0;
652 
653         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
654         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
655         if (ovl) {
656           for (c=0;c<supportSize;c++) {
657             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
658               fcells[nc++] = support[c];
659             }
660           }
661         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
662         for (c=0;c<nc;c++) {
663           const    PetscInt *cone;
664           int      vids[8];
665           PetscInt i,cell,cl,nv,cid = -1,mid = -1;
666 
667           cell = fcells[c];
668           ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr);
669           for (cl=0;cl<fpc;cl++)
670             if (cone[cl] == p)
671               break;
672 
673           /* face material id and type */
674           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
675           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
676           /* vertex ids */
677           ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
678           for (i=0;i<vpf;i++) {
679             ierr = PetscViewerASCIIPrintf(viewer," %D",vids[faces[cl*vpf+i]]);CHKERRQ(ierr);
680           }
681           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
682         }
683         bf = bf-nc;
684       }
685     }
686     if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
687     ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr);
688     ierr = PetscFree(fcells);CHKERRQ(ierr);
689   }
690 
691   /* mark owned vertices */
692   if (ovl) {
693     ierr = PetscBTMemzero(vEnd-vStart,pown);CHKERRQ(ierr);
694     for (p=cStart;p<cEnd;p++) {
695       PetscInt i,closureSize,*closure = NULL;
696 
697       if (gNum[p-cStart] < 0) continue;
698       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
699       for (i=0;i<closureSize;i++) {
700         const PetscInt pp = closure[2*i];
701 
702         if (pp >= vStart && pp < vEnd) {
703           ierr = PetscBTSet(pown,pp-vStart);CHKERRQ(ierr);
704         }
705       }
706       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
707     }
708   }
709   if (globalNum) {
710     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
711   }
712 
713   /* vertex_parents (Non-conforming meshes) */
714   parentSection  = NULL;
715   if (enable_ncmesh) {
716     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
717   }
718   if (parentSection) {
719     PetscInt vp,gvp;
720 
721     for (vp=0,p=vStart;p<vEnd;p++) {
722       DMLabel  dlabel;
723       PetscInt parent,depth;
724 
725       if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
726       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
727       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
728       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
729       if (parent != p) vp++;
730     }
731     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
732     if (gvp) {
733       PetscInt  maxsupp;
734       PetscBool *skip = NULL;
735 
736       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
737       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
738       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
739       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
740       for (p=vStart;p<vEnd;p++) {
741         DMLabel  dlabel;
742         PetscInt parent;
743 
744         if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
745         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
746         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
747         if (parent != p) {
748           int            vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
749           PetscInt       i,nv,size,n,numChildren,depth = -1;
750           const PetscInt *children;
751           ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr);
752           switch (size) {
753             case 2: /* edge */
754               nv   = 0;
755               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
756               ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
757               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
758               for (i=0;i<nv;i++) {
759                 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
760               }
761               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
762               vp--;
763               break;
764             case 4: /* face */
765               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
766               for (n=0;n<numChildren;n++) {
767                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
768                 if (!depth) {
769                   const PetscInt *hvsupp,*hesupp,*cone;
770                   PetscInt       hvsuppSize,hesuppSize,coneSize;
771                   PetscInt       hv = children[n],he = -1,f;
772 
773                   ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
774                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
775                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
776                   for (i=0;i<hvsuppSize;i++) {
777                     PetscInt ep;
778                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
779                     if (ep != hvsupp[i]) {
780                       he = hvsupp[i];
781                     } else {
782                       skip[i] = PETSC_TRUE;
783                     }
784                   }
785                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
786                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
787                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
788                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
789                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
790                   for (f=0;f<hesuppSize;f++) {
791                     PetscInt j;
792 
793                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
794                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
795                     for (j=0;j<coneSize;j++) {
796                       PetscInt k;
797                       for (k=0;k<hvsuppSize;k++) {
798                         if (hvsupp[k] == cone[j]) {
799                           skip[k] = PETSC_TRUE;
800                           break;
801                         }
802                       }
803                     }
804                   }
805                   for (i=0;i<hvsuppSize;i++) {
806                     if (!skip[i]) {
807                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
808                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
809                     }
810                   }
811                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
812                   for (i=0;i<2;i++) {
813                     ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);CHKERRQ(ierr);
814                   }
815                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
816                   vp--;
817                 }
818               }
819               break;
820             default:
821               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size);
822           }
823         }
824       }
825       ierr = PetscFree(skip);CHKERRQ(ierr);
826     }
827     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
828   }
829   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
830 
831   /* vertices */
832   if (hovec) { /* higher-order meshes */
833     const char *fec;
834     PetscInt   i,n;
835 
836     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
837     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
838     ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
839     ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr);
840     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
841     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr);
842     ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
843     ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
844     ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr);
845     ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr);
846     if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
847     for (i=0;i<n/sdim;i++) {
848       PetscInt s;
849       for (s=0;s<sdim;s++) {
850         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr);
851       }
852       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
853     }
854     ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr);
855   } else {
856     ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr);
857     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
858     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr);
859     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
860     ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
861     for (p=0;p<nvert/sdim;p++) {
862       PetscInt s;
863       for (s=0;s<sdim;s++) {
864         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr);
865       }
866       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
867     }
868     ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
869   }
870   PetscFunctionReturn(0);
871 }
872 
873 /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
874 PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
875 {
876   PetscErrorCode ierr;
877   PetscBool      isglvis,isascii;
878 
879   PetscFunctionBegin;
880   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
881   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
882   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
883   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
884   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
885   if (isglvis) {
886     PetscViewer          view;
887     PetscViewerGLVisType type;
888 
889     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
890     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
891     if (view) { /* in the socket case, it may happen that the connection failed */
892       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
893         PetscMPIInt size,rank;
894         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
895         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
896         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
897       }
898       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
899       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
900       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
901         PetscInt    dim;
902         const char* name;
903 
904         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
905         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
906         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
907         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
908       }
909     }
910   } else {
911     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
912   }
913   PetscFunctionReturn(0);
914 }
915