Actual source code: plexglvis.c

petsc-3.13.0 2020-03-29
Report Typos and Errors
  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>

  9: typedef struct {
 10:   PetscInt   nf;
 11:   VecScatter *scctx;
 12: } GLVisViewerCtx;

 14: static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
 15: {
 16:   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
 17:   PetscInt       i;

 21:   for (i=0;i<ctx->nf;i++) {
 22:     VecScatterDestroy(&ctx->scctx[i]);
 23:   }
 24:   PetscFree(ctx->scctx);
 25:   PetscFree(vctx);
 26:   return(0);
 27: }

 29: static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
 30: {
 31:   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
 32:   PetscInt       f;

 36:   for (f=0;f<nf;f++) {
 37:     VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);
 38:     VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);
 39:   }
 40:   return(0);
 41: }

 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;

 61:   DMGetDimension(dm,&dim);
 62:   DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
 63:   DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
 64:   PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);
 65:   if (!globalNum) {
 66:     DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);
 67:     PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);
 68:     PetscObjectDereference((PetscObject)globalNum);
 69:   }
 70:   ISGetIndices(globalNum,&gNum);
 71:   PetscBTCreate(vEnd-vStart,&vown);
 72:   for (c = cStart, totc = 0; c < cEnd; c++) {
 73:     if (gNum[c-cStart] >= 0) {
 74:       PetscInt i,numPoints,*points = NULL;

 76:       totc++;
 77:       DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);
 78:       for (i=0;i<numPoints*2;i+= 2) {
 79:         if ((points[i] >= vStart) && (points[i] < vEnd)) {
 80:           PetscBTSet(vown,points[i]-vStart);
 81:         }
 82:       }
 83:       DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);
 84:     }
 85:   }
 86:   for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++;

 88:   DMCreateLocalVector(dm,&xlocal);
 89:   VecGetLocalSize(xlocal,&totdofs);
 90:   DMGetLocalSection(dm,&s);
 91:   PetscSectionGetNumFields(s,&nfields);
 92:   for (f=0,maxfields=0;f<nfields;f++) {
 93:     PetscInt bs;

 95:     PetscSectionGetFieldComponents(s,f,&bs);
 96:     maxfields += bs;
 97:   }
 98:   PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);
 99:   PetscNew(&ctx);
100:   PetscCalloc1(maxfields,&ctx->scctx);
101:   DMGetDS(dm,&ds);
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;

109:       PetscSectionGetFieldName(s,f,&fname);
110:       PetscStrlen(fname,&len);
111:       if (len) {
112:         PetscStrcpy(name,fname);
113:       } else {
114:         PetscSNPrintf(name,256,"Field%D",f);
115:       }
116:       PetscDSGetDiscretization(ds,f,&disc);
117:       if (disc) {
118:         PetscClassId id;
119:         PetscInt     Nc;
120:         char         fec[64];

122:         PetscObjectGetClassId(disc, &id);
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;

130:           PetscFEGetNumComponents(fem,&Nc);
131:           PetscFEGetDualSpace(fem,&sp);
132:           PetscDualSpaceGetType(sp,&spname);
133:           PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);
134:           if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
135:           PetscDualSpaceLagrangeGetContinuity(sp,&continuous);
136:           PetscDualSpaceGetOrder(sp,&order);
137:           if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */
138:             PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);
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:             PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order);
143:           }
144:           PetscStrallocpy(name,&fieldname[ctx->nf]);
145:           bs[ctx->nf]   = Nc;
146:           dims[ctx->nf] = dim;
147:           if (H1) {
148:             nlocal[ctx->nf] = Nc * Nv;
149:             PetscStrallocpy(fec,&fec_type[ctx->nf]);
150:             VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);
151:             for (i=0,cum=0;i<vEnd-vStart;i++) {
152:               PetscInt j,off;

154:               if (PetscUnlikely(!PetscBTLookup(vown,i))) continue;
155:               PetscSectionGetFieldOffset(s,i+vStart,f,&off);
156:               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
157:             }
158:             ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);
159:           } else {
160:             nlocal[ctx->nf] = Nc * totc;
161:             PetscStrallocpy(fec,&fec_type[ctx->nf]);
162:             VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);
163:             for (i=0,cum=0;i<cEnd-cStart;i++) {
164:               PetscInt j,off;

166:               if (PetscUnlikely(gNum[i] < 0)) continue;
167:               PetscSectionGetFieldOffset(s,i+cStart,f,&off);
168:               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
169:             }
170:             ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);
171:           }
172:           VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);
173:           VecDestroy(&xfield);
174:           ISDestroy(&isfield);
175:           ctx->nf++;
176:         } else if (id == PETSCFV_CLASSID) {
177:           PetscInt c;

179:           PetscFVGetNumComponents((PetscFV)disc,&Nc);
180:           PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim);
181:           for (c = 0; c < Nc; c++) {
182:             char comp[256];
183:             PetscSNPrintf(comp,256,"%s-Comp%D",name,c);
184:             PetscStrallocpy(comp,&fieldname[ctx->nf]);
185:             bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */
186:             nlocal[ctx->nf] = totc;
187:             dims[ctx->nf] = dim;
188:             PetscStrallocpy(fec,&fec_type[ctx->nf]);
189:             VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);
190:             for (i=0,cum=0;i<cEnd-cStart;i++) {
191:               PetscInt off;

193:               if (PetscUnlikely(gNum[i])<0) continue;
194:               PetscSectionGetFieldOffset(s,i+cStart,f,&off);
195:               idxs[cum++] = off + c;
196:             }
197:             ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);
198:             VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);
199:             VecDestroy(&xfield);
200:             ISDestroy(&isfield);
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:   PetscBTDestroy(&vown);
208:   VecDestroy(&xlocal);
209:   ISRestoreIndices(globalNum,&gNum);

211:   /* create work vectors */
212:   for (f=0;f<ctx->nf;f++) {
213:     VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);
214:     PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);
215:     VecSetBlockSize(Ufield[f],bs[f]);
216:     VecSetDM(Ufield[f],dm);
217:   }

219:   /* customize the viewer */
220:   PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);
221:   for (f=0;f<ctx->nf;f++) {
222:     PetscFree(fieldname[f]);
223:     PetscFree(fec_type[f]);
224:     VecDestroy(&Ufield[f]);
225:   }
226:   PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);
227:   return(0);
228: }

230: typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid;

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 } };

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 } };

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;

249:   DMPlexGetDepthLabel(dm,&dlabel);
250:   DMLabelGetValue(dlabel,p,&pdepth);
251:   DMPlexGetConeSize(dm,p,&csize);
252:   DMPlexGetDepth(dm,&depth);
253:   DMGetDimension(dm,&dim);
254:   if (label) {
255:     DMLabelGetValue(label,p,mid);
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 defined PETSC_USE_DEBUG
260:     if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim);
261:     if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p);
262:     if (depth != 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p);
263: #endif
264:     *cid = mfem_table_cid_unint[dim][csize];
265:   } else {
266: #if defined PETSC_USE_DEBUG
267:     if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p);
268:     if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p);
269: #endif
270:     *cid = mfem_table_cid[pdepth][csize];
271:   }
272:   return(0);
273: }

275: static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[])
276: {
277:   PetscInt       dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;

281:   DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
282:   DMGetDimension(dm,&dim);
283:   sdim = dim;
284:   if (csec) {
285:     PetscInt sStart,sEnd;

287:     DMGetCoordinateDim(dm,&sdim);
288:     PetscSectionGetChart(csec,&sStart,&sEnd);
289:     PetscSectionGetOffset(csec,vStart,&off);
290:     off  = off/sdim;
291:     if (p >= sStart && p < sEnd) {
292:       PetscSectionGetDof(csec,p,&dof);
293:     }
294:   }
295:   if (!dof) {
296:     DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);
297:     for (i=0,q=0;i<numPoints*2;i+= 2)
298:       if ((points[i] >= vStart) && (points[i] < vEnd))
299:         vids[q++] = points[i]-vStart+off;
300:     DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);
301:   } else {
302:     PetscSectionGetOffset(csec,p,&off);
303:     PetscSectionGetDof(csec,p,&dof);
304:     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
305:   }
306:   *nv = q;
307:   return(0);
308: }

310: /*
311:    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
312:    Higher order meshes are also supported
313: */
314: static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
315: {
316:   DMLabel              label;
317:   PetscSection         coordSection,parentSection;
318:   Vec                  coordinates,hovec;
319:   const PetscScalar    *array;
320:   PetscInt             bf,p,sdim,dim,depth,novl,minl;
321:   PetscInt             cStart,cEnd,vStart,vEnd,nvert;
322:   PetscMPIInt          size;
323:   PetscBool            localized,isascii;
324:   PetscBool            enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE;
325:   PetscBT              pown,vown;
326:   PetscErrorCode       ierr;
327:   PetscContainer       glvis_container;
328:   PetscBool            cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
329:   PetscBool            enable_emark,enable_bmark;
330:   const char           *fmt;
331:   char                 emark[64] = "",bmark[64] = "";

336:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
337:   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
338:   MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
339:   if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
340:   DMGetDimension(dm,&dim);

342:   /* get container: determines if a process visualizes is portion of the data or not */
343:   PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
344:   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
345:   {
346:     PetscViewerGLVisInfo glvis_info;
347:     PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
348:     enabled = glvis_info->enabled;
349:     fmt     = glvis_info->fmt;
350:   }

352:   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
353:      DMPlex does not currently support HO meshes, so there's no API for this */
354:   PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);

356:   DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
357:   DMPlexGetGhostCellStratum(dm,&p,NULL);
358:   if (p >= 0) cEnd = p;
359:   DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
360:   DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);
361:   DMGetCoordinatesLocalized(dm,&localized);
362:   if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
363:   DMGetCoordinateSection(dm,&coordSection);
364:   DMGetCoordinateDim(dm,&sdim);
365:   DMGetCoordinatesLocal(dm,&coordinates);
366:   if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");

368:   /*
369:      a couple of sections of the mesh specification are disabled
370:        - 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)
371:        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
372:                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
373:   */
374:   enable_boundary = PETSC_FALSE;
375:   enable_ncmesh   = PETSC_FALSE;
376:   enable_mfem     = PETSC_FALSE;
377:   enable_emark    = PETSC_FALSE;
378:   enable_bmark    = PETSC_FALSE;
379:   /* I'm tired of problems with negative values in the markers, disable them */
380:   PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");
381:   PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);
382:   PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL);
383:   PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL);
384:   PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL);
385:   PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);
386:   PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);
387:   PetscOptionsEnd();
388:   if (enable_bmark) enable_boundary = PETSC_TRUE;

390:   MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
391:   if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
392:   DMPlexGetDepth(dm,&depth);
393:   if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
394:                                                              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
395:   if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
396:                                                            "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
397:   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
398:     if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth);
399:     cellvertex = PETSC_TRUE;
400:   }

402:   /* Identify possible cells in the overlap */
403:   novl = 0;
404:   pown = NULL;
405:   if (size > 1) {
406:     IS             globalNum = NULL;
407:     const PetscInt *gNum;
408:     PetscBool      ovl  = PETSC_FALSE;

410:     PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);
411:     if (!globalNum) {
412:       if (view_ovl) {
413:         ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum);
414:       } else {
415:         DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);
416:       }
417:       PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);
418:       PetscObjectDereference((PetscObject)globalNum);
419:     }
420:     ISGetIndices(globalNum,&gNum);
421:     for (p=cStart; p<cEnd; p++) {
422:       if (gNum[p-cStart] < 0) {
423:         ovl = PETSC_TRUE;
424:         novl++;
425:       }
426:     }
427:     if (ovl) {
428:       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
429:          TODO: garbage collector? attach pown to dm?  */
430:       PetscBTCreate(cEnd-cStart,&pown);
431:       for (p=cStart; p<cEnd; p++) {
432:         if (gNum[p-cStart] < 0) continue;
433:         else {
434:           PetscBTSet(pown,p-cStart);
435:         }
436:       }
437:     }
438:     ISRestoreIndices(globalNum,&gNum);
439:   }

441:   /* return if this process is disabled */
442:   if (!enabled) {
443:     PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
444:     PetscViewerASCIIPrintf(viewer,"\ndimension\n");
445:     PetscViewerASCIIPrintf(viewer,"%D\n",dim);
446:     PetscViewerASCIIPrintf(viewer,"\nelements\n");
447:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
448:     PetscViewerASCIIPrintf(viewer,"\nboundary\n");
449:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
450:     PetscViewerASCIIPrintf(viewer,"\nvertices\n");
451:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
452:     PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
453:     PetscBTDestroy(&pown);
454:     return(0);
455:   }

457:   if (enable_mfem) {
458:     if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
459:       PetscInt    vpc = 0;
460:       char        fec[64];
461:       PetscInt    vids[8] = {0,1,2,3,4,5,6,7};
462:       PetscInt    hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3};
463:       PetscInt    quadv[8] = {0,1,3,2}, triv[3] = {0,1,2};
464:       PetscInt    *dof = NULL;
465:       PetscScalar *array,*ptr;

467:       PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);
468:       if (cEnd-cStart) {
469:         PetscInt fpc;

471:         DMPlexGetConeSize(dm,cStart,&fpc);
472:         switch(dim) {
473:           case 1:
474:             vpc = 2;
475:             dof = hexv;
476:             break;
477:           case 2:
478:             switch (fpc) {
479:               case 3:
480:                 vpc = 3;
481:                 dof = triv;
482:                 break;
483:               case 4:
484:                 vpc = 4;
485:                 dof = quadv;
486:                 break;
487:               default:
488:                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
489:                 break;
490:             }
491:             break;
492:           case 3:
493:             switch (fpc) {
494:               case 4: /* TODO: still need to understand L2 ordering for tets */
495:                 vpc = 4;
496:                 dof = tetv;
497:                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
498:                 break;
499:               case 6:
500:                 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc);
501:                 vpc = 8;
502:                 dof = hexv;
503:                 break;
504:               case 8:
505:                 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
506:                 vpc = 8;
507:                 dof = hexv;
508:                 break;
509:               default:
510:                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
511:                 break;
512:             }
513:             break;
514:           default:
515:             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
516:             break;
517:         }
518:         DMPlexReorderCell(dm,cStart,vids);
519:       }
520:       if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
521:       VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);
522:       PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);
523:       PetscObjectDereference((PetscObject)hovec);
524:       PetscObjectSetName((PetscObject)hovec,fec);
525:       VecGetArray(hovec,&array);
526:       ptr  = array;
527:       for (p=cStart;p<cEnd;p++) {
528:         PetscInt    csize,v,d;
529:         PetscScalar *vals = NULL;

531:         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
532:         DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);
533:         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);
534:         for (v=0;v<vpc;v++) {
535:           for (d=0;d<sdim;d++) {
536:             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
537:           }
538:         }
539:         ptr += vpc*sdim;
540:         DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);
541:       }
542:       VecRestoreArray(hovec,&array);
543:     }
544:   }
545:   /* if we have high-order coordinates in 3D, we need to specify the boundary */
546:   if (hovec && dim == 3) enable_boundary = PETSC_TRUE;

548:   /* header */
549:   PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");

551:   /* topological dimension */
552:   PetscViewerASCIIPrintf(viewer,"\ndimension\n");
553:   PetscViewerASCIIPrintf(viewer,"%D\n",dim);

555:   /* elements */
556:   minl = 1;
557:   label = NULL;
558:   if (enable_emark) {
559:     PetscInt lminl = PETSC_MAX_INT;

561:     DMGetLabel(dm,emark,&label);
562:     if (label) {
563:       IS       vals;
564:       PetscInt ldef;

566:       DMLabelGetDefaultValue(label,&ldef);
567:       DMLabelGetValueIS(label,&vals);
568:       ISGetMinMax(vals,&lminl,NULL);
569:       ISDestroy(&vals);
570:       lminl = PetscMin(ldef,lminl);
571:     }
572:     MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));
573:     if (minl == PETSC_MAX_INT) minl = 1;
574:   }
575:   PetscViewerASCIIPrintf(viewer,"\nelements\n");
576:   PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);
577:   for (p=cStart;p<cEnd;p++) {
578:     PetscInt       vids[8];
579:     PetscInt       i,nv = 0,cid = -1,mid = 1;

581:     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
582:     DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);
583:     DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);
584:     DMPlexReorderCell(dm,p,vids);
585:     PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
586:     for (i=0;i<nv;i++) {
587:       PetscViewerASCIIPrintf(viewer," %D",vids[i]);
588:     }
589:     PetscViewerASCIIPrintf(viewer,"\n");
590:   }

592:   /* boundary */
593:   PetscViewerASCIIPrintf(viewer,"\nboundary\n");
594:   if (!enable_boundary) {
595:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
596:   } else {
597:     DMLabel  perLabel;
598:     PetscBT  bfaces;
599:     PetscInt fStart,fEnd,*fcells;

601:     DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);
602:     PetscBTCreate(fEnd-fStart,&bfaces);
603:     DMPlexGetMaxSizes(dm,NULL,&p);
604:     PetscMalloc1(p,&fcells);
605:     DMGetLabel(dm,"glvis_periodic_cut",&perLabel);
606:     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
607:       DMCreateLabel(dm,"glvis_periodic_cut");
608:       DMGetLabel(dm,"glvis_periodic_cut",&perLabel);
609:       DMLabelSetDefaultValue(perLabel,1);
610:       for (p=cStart;p<cEnd;p++) {
611:         DMPolytopeType cellType;
612:         PetscInt       dof;

614:         DMPlexGetCellType(dm,p,&cellType);
615:         PetscSectionGetDof(coordSection,p,&dof);
616:         if (dof) {
617:           PetscInt    uvpc, v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
618:           PetscScalar *vals = NULL;

620:           uvpc = DMPolytopeTypeGetNumVertices(cellType);
621:           if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
622:           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);
623:           DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);
624:           DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);
625:           for (v=0;v<cellClosureSize;v++)
626:             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
627:               vidxs = cellClosure + 2*v;
628:               break;
629:             }
630:           if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
631:           for (v=0;v<uvpc;v++) {
632:             PetscInt s;

634:             for (s=0;s<sdim;s++) {
635:               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) {
636:                 DMLabelSetValue(perLabel,vidxs[2*v],2);
637:               }
638:             }
639:           }
640:           DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);
641:           DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);
642:         }
643:       }
644:       if (dim > 1) {
645:         PetscInt eEnd,eStart;

647:         DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);
648:         for (p=eStart;p<eEnd;p++) {
649:           const PetscInt *cone;
650:           PetscInt       coneSize,i;
651:           PetscBool      ispe = PETSC_TRUE;

653:           DMPlexGetCone(dm,p,&cone);
654:           DMPlexGetConeSize(dm,p,&coneSize);
655:           for (i=0;i<coneSize;i++) {
656:             PetscInt v;

658:             DMLabelGetValue(perLabel,cone[i],&v);
659:             ispe = (PetscBool)(ispe && (v==2));
660:           }
661:           if (ispe && coneSize) {
662:             PetscInt       ch, numChildren;
663:             const PetscInt *children;

665:             DMLabelSetValue(perLabel,p,2);
666:             DMPlexGetTreeChildren(dm,p,&numChildren,&children);
667:             for (ch = 0; ch < numChildren; ch++) {
668:               DMLabelSetValue(perLabel,children[ch],2);
669:             }
670:           }
671:         }
672:         if (dim > 2) {
673:           for (p=fStart;p<fEnd;p++) {
674:             const PetscInt *cone;
675:             PetscInt       coneSize,i;
676:             PetscBool      ispe = PETSC_TRUE;

678:             DMPlexGetCone(dm,p,&cone);
679:             DMPlexGetConeSize(dm,p,&coneSize);
680:             for (i=0;i<coneSize;i++) {
681:               PetscInt v;

683:               DMLabelGetValue(perLabel,cone[i],&v);
684:               ispe = (PetscBool)(ispe && (v==2));
685:             }
686:             if (ispe && coneSize) {
687:               PetscInt       ch, numChildren;
688:               const PetscInt *children;

690:               DMLabelSetValue(perLabel,p,2);
691:               DMPlexGetTreeChildren(dm,p,&numChildren,&children);
692:               for (ch = 0; ch < numChildren; ch++) {
693:                 DMLabelSetValue(perLabel,children[ch],2);
694:               }
695:             }
696:           }
697:         }
698:       }
699:     }
700:     for (p=fStart;p<fEnd;p++) {
701:       const PetscInt *support;
702:       PetscInt       supportSize;
703:       PetscBool      isbf = PETSC_FALSE;

705:       DMPlexGetSupportSize(dm,p,&supportSize);
706:       if (pown) {
707:         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
708:         PetscInt  i;

710:         DMPlexGetSupport(dm,p,&support);
711:         for (i=0;i<supportSize;i++) {
712:           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
713:           else has_ghost = PETSC_TRUE;
714:         }
715:         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
716:       } else {
717:         isbf = (PetscBool)(supportSize == 1);
718:       }
719:       if (!isbf && perLabel) {
720:         const PetscInt *cone;
721:         PetscInt       coneSize,i;

723:         DMPlexGetCone(dm,p,&cone);
724:         DMPlexGetConeSize(dm,p,&coneSize);
725:         isbf = PETSC_TRUE;
726:         for (i=0;i<coneSize;i++) {
727:           PetscInt v,d;

729:           DMLabelGetValue(perLabel,cone[i],&v);
730:           DMLabelGetDefaultValue(perLabel,&d);
731:           isbf = (PetscBool)(isbf && v != d);
732:         }
733:       }
734:       if (isbf) {
735:         PetscBTSet(bfaces,p-fStart);
736:       }
737:     }
738:     /* count boundary faces */
739:     for (p=fStart,bf=0;p<fEnd;p++) {
740:       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
741:         const PetscInt *support;
742:         PetscInt       supportSize,c;

744:         DMPlexGetSupportSize(dm,p,&supportSize);
745:         DMPlexGetSupport(dm,p,&support);
746:         for (c=0;c<supportSize;c++) {
747:           const    PetscInt *cone;
748:           PetscInt cell,cl,coneSize;

750:           cell = support[c];
751:           if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue;
752:           DMPlexGetCone(dm,cell,&cone);
753:           DMPlexGetConeSize(dm,cell,&coneSize);
754:           for (cl=0;cl<coneSize;cl++) {
755:             if (cone[cl] == p) {
756:               bf += 1;
757:               break;
758:             }
759:           }
760:         }
761:       }
762:     }
763:     minl = 1;
764:     label = NULL;
765:     if (enable_bmark) {
766:       PetscInt lminl = PETSC_MAX_INT;

768:       DMGetLabel(dm,bmark,&label);
769:       if (label) {
770:         IS       vals;
771:         PetscInt ldef;

773:         DMLabelGetDefaultValue(label,&ldef);
774:         DMLabelGetValueIS(label,&vals);
775:         ISGetMinMax(vals,&lminl,NULL);
776:         ISDestroy(&vals);
777:         lminl = PetscMin(ldef,lminl);
778:       }
779:       MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));
780:       if (minl == PETSC_MAX_INT) minl = 1;
781:     }
782:     PetscViewerASCIIPrintf(viewer,"%D\n",bf);
783:     for (p=fStart;p<fEnd;p++) {
784:       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
785:         const PetscInt *support;
786:         PetscInt       supportSize,c,nc = 0;

788:         DMPlexGetSupportSize(dm,p,&supportSize);
789:         DMPlexGetSupport(dm,p,&support);
790:         if (pown) {
791:           for (c=0;c<supportSize;c++) {
792:             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
793:               fcells[nc++] = support[c];
794:             }
795:           }
796:         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
797:         for (c=0;c<nc;c++) {
798:           const DMPolytopeType *faceTypes;
799:           DMPolytopeType       cellType;
800:           const PetscInt       *faceSizes,*cone;
801:           PetscInt             vids[8],*faces,st,i,coneSize,cell,cl,nv,cid = -1,mid = -1;

803:           cell = fcells[c];
804:           DMPlexGetCone(dm,cell,&cone);
805:           DMPlexGetConeSize(dm,cell,&coneSize);
806:           for (cl=0;cl<coneSize;cl++)
807:             if (cone[cl] == p)
808:               break;
809:           if (cl == coneSize) continue;

811:           /* face material id and type */
812:           DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);
813:           PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
814:           /* vertex ids */
815:           DMPlexGetCellType(dm,cell,&cellType);
816:           DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);
817:           DMPlexGetRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);
818:           st = 0;
819:           for (i=0;i<cl;i++) st += faceSizes[i];
820:           DMPlexInvertCell(faceTypes[cl],faces + st);
821:           for (i=0;i<faceSizes[cl];i++) {
822:             PetscViewerASCIIPrintf(viewer," %d",faces[st+i]);
823:           }
824:           PetscViewerASCIIPrintf(viewer,"\n");
825:           DMPlexRestoreRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);
826:           bf -= 1;
827:         }
828:       }
829:     }
830:     if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
831:     PetscBTDestroy(&bfaces);
832:     PetscFree(fcells);
833:   }

835:   /* mark owned vertices */
836:   vown = NULL;
837:   if (pown) {
838:     PetscBTCreate(vEnd-vStart,&vown);
839:     for (p=cStart;p<cEnd;p++) {
840:       PetscInt i,closureSize,*closure = NULL;

842:       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
843:       DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
844:       for (i=0;i<closureSize;i++) {
845:         const PetscInt pp = closure[2*i];

847:         if (pp >= vStart && pp < vEnd) {
848:           PetscBTSet(vown,pp-vStart);
849:         }
850:       }
851:       DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
852:     }
853:   }

855:   /* vertex_parents (Non-conforming meshes) */
856:   parentSection  = NULL;
857:   if (enable_ncmesh) {
858:     DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);
859:   }
860:   if (parentSection) {
861:     PetscInt vp,gvp;

863:     for (vp=0,p=vStart;p<vEnd;p++) {
864:       DMLabel  dlabel;
865:       PetscInt parent,depth;

867:       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
868:       DMPlexGetDepthLabel(dm,&dlabel);
869:       DMLabelGetValue(dlabel,p,&depth);
870:       DMPlexGetTreeParent(dm,p,&parent,NULL);
871:       if (parent != p) vp++;
872:     }
873:     MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
874:     if (gvp) {
875:       PetscInt  maxsupp;
876:       PetscBool *skip = NULL;

878:       PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");
879:       PetscViewerASCIIPrintf(viewer,"%D\n",vp);
880:       DMPlexGetMaxSizes(dm,NULL,&maxsupp);
881:       PetscMalloc1(maxsupp,&skip);
882:       for (p=vStart;p<vEnd;p++) {
883:         DMLabel  dlabel;
884:         PetscInt parent;

886:         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
887:         DMPlexGetDepthLabel(dm,&dlabel);
888:         DMPlexGetTreeParent(dm,p,&parent,NULL);
889:         if (parent != p) {
890:           PetscInt       vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
891:           PetscInt       i,nv,ssize,n,numChildren,depth = -1;
892:           const PetscInt *children;

894:           DMPlexGetConeSize(dm,parent,&ssize);
895:           switch (ssize) {
896:             case 2: /* edge */
897:               nv   = 0;
898:               DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);
899:               PetscViewerASCIIPrintf(viewer,"%D",p-vStart);
900:               for (i=0;i<nv;i++) {
901:                 PetscViewerASCIIPrintf(viewer," %D",vids[i]);
902:               }
903:               PetscViewerASCIIPrintf(viewer,"\n");
904:               vp--;
905:               break;
906:             case 4: /* face */
907:               DMPlexGetTreeChildren(dm,parent,&numChildren,&children);
908:               for (n=0;n<numChildren;n++) {
909:                 DMLabelGetValue(dlabel,children[n],&depth);
910:                 if (!depth) {
911:                   const PetscInt *hvsupp,*hesupp,*cone;
912:                   PetscInt       hvsuppSize,hesuppSize,coneSize;
913:                   PetscInt       hv = children[n],he = -1,f;

915:                   PetscArrayzero(skip,maxsupp);
916:                   DMPlexGetSupportSize(dm,hv,&hvsuppSize);
917:                   DMPlexGetSupport(dm,hv,&hvsupp);
918:                   for (i=0;i<hvsuppSize;i++) {
919:                     PetscInt ep;
920:                     DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);
921:                     if (ep != hvsupp[i]) {
922:                       he = hvsupp[i];
923:                     } else {
924:                       skip[i] = PETSC_TRUE;
925:                     }
926:                   }
927:                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
928:                   DMPlexGetCone(dm,he,&cone);
929:                   vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
930:                   DMPlexGetSupportSize(dm,he,&hesuppSize);
931:                   DMPlexGetSupport(dm,he,&hesupp);
932:                   for (f=0;f<hesuppSize;f++) {
933:                     PetscInt j;

935:                     DMPlexGetCone(dm,hesupp[f],&cone);
936:                     DMPlexGetConeSize(dm,hesupp[f],&coneSize);
937:                     for (j=0;j<coneSize;j++) {
938:                       PetscInt k;
939:                       for (k=0;k<hvsuppSize;k++) {
940:                         if (hvsupp[k] == cone[j]) {
941:                           skip[k] = PETSC_TRUE;
942:                           break;
943:                         }
944:                       }
945:                     }
946:                   }
947:                   for (i=0;i<hvsuppSize;i++) {
948:                     if (!skip[i]) {
949:                       DMPlexGetCone(dm,hvsupp[i],&cone);
950:                       vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
951:                     }
952:                   }
953:                   PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);
954:                   for (i=0;i<2;i++) {
955:                     PetscViewerASCIIPrintf(viewer," %D",vids[i]-vStart);
956:                   }
957:                   PetscViewerASCIIPrintf(viewer,"\n");
958:                   vp--;
959:                 }
960:               }
961:               break;
962:             default:
963:               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize);
964:           }
965:         }
966:       }
967:       PetscFree(skip);
968:     }
969:     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
970:   }
971:   PetscBTDestroy(&pown);
972:   PetscBTDestroy(&vown);

974:   /* vertices */
975:   if (hovec) { /* higher-order meshes */
976:     const char *fec;
977:     PetscInt   i,n,s;

979:     PetscViewerASCIIPrintf(viewer,"\nvertices\n");
980:     PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);
981:     PetscViewerASCIIPrintf(viewer,"nodes\n");
982:     PetscObjectGetName((PetscObject)hovec,&fec);
983:     PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");
984:     PetscViewerASCIIPrintf(viewer,"%s\n",fec);
985:     PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);
986:     PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n"); /*Ordering::byVDIM*/
987:     VecGetArrayRead(hovec,&array);
988:     VecGetLocalSize(hovec,&n);
989:     if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
990:     for (i=0;i<n/sdim;i++) {
991:       for (s=0;s<sdim;s++) {
992:         PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));
993:       }
994:       PetscViewerASCIIPrintf(viewer,"\n");
995:     }
996:     VecRestoreArrayRead(hovec,&array);
997:   } else {
998:     VecGetLocalSize(coordinates,&nvert);
999:     PetscViewerASCIIPrintf(viewer,"\nvertices\n");
1000:     PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);
1001:     PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
1002:     VecGetArrayRead(coordinates,&array);
1003:     for (p=0;p<nvert/sdim;p++) {
1004:       PetscInt s;
1005:       for (s=0;s<sdim;s++) {
1006:         PetscReal v = PetscRealPart(array[p*sdim+s]);

1008:         PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);
1009:       }
1010:       PetscViewerASCIIPrintf(viewer,"\n");
1011:     }
1012:     VecRestoreArrayRead(coordinates,&array);
1013:   }
1014:   return(0);
1015: }

1017: PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
1018: {
1021:   DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);
1022:   return(0);
1023: }