Actual source code: plexglvis.c

petsc-3.8.0 2017-09-26
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;
 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;

 61:   DMGetDimension(dm,&dim);
 62:   DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
 63:   DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
 64:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
 65:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
 66:   DMPlexGetCellNumbering(dm,&globalNum);
 67:   ISGetIndices(globalNum,&gNum);
 68:   PetscBTCreate(vEnd-vStart,&vown);
 69:   for (c = cStart, totc = 0; c < cEnd; c++) {
 70:     if (gNum[c-cStart] >= 0) {
 71:       PetscInt i,numPoints,*points = NULL;

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

 85:   DMCreateLocalVector(dm,&xlocal);
 86:   VecGetLocalSize(xlocal,&totdofs);
 87:   DMGetDefaultSection(dm,&s);
 88:   PetscSectionGetNumFields(s,&nfields);
 89:   for (f=0,maxfields=0;f<nfields;f++) {
 90:     PetscInt bs;

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

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

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

127:           PetscFEGetNumComponents(fem,&Nc);
128:           PetscFEGetDualSpace(fem,&sp);
129:           PetscDualSpaceGetType(sp,&spname);
130:           PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);
131:           if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
132:           PetscDualSpaceLagrangeGetContinuity(sp,&continuous);
133:           if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported");
134:           PetscDualSpaceGetOrder(sp,&order);
135:           if (continuous && order > 0) {
136:             PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%dD_P1",dim);
137:           } else {
138:             H1   = PETSC_FALSE;
139:             PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P%d",dim,order);
140:           }
141:           PetscStrallocpy(name,&fieldname[ctx->nf]);
142:           bs[ctx->nf]   = Nc;
143:           dims[ctx->nf] = dim;
144:           if (H1) {
145:             nlocal[ctx->nf] = Nc * Nv;
146:             PetscStrallocpy(fec,&fec_type[ctx->nf]);
147:             VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);
148:             for (i=0,cum=0;i<vEnd-vStart;i++) {
149:               PetscInt j,off;

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

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

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

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

208:   /* customize the viewer */
209:   PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fieldname,(const char**)fec_type,nlocal,bs,dims,DMPlexSampleGLVisFields_Private,ctx,DestroyGLVisViewerCtx_Private);
210:   for (f=0;f<ctx->nf;f++) {
211:     PetscFree(fieldname[f]);
212:     PetscFree(fec_type[f]);
213:   }
214:   PetscFree6(fieldname,nlocal,bs,dims,fec_type,idxs);
215:   return(0);
216: }

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

220: MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
221:                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
222:                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF, MFEM_UNDEF},
223:                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } };

225: static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid)
226: {
227:   DMLabel        dlabel;
228:   PetscInt       depth,csize;

232:   DMPlexGetDepthLabel(dm,&dlabel);
233:   DMLabelGetValue(dlabel,p,&depth);
234:   DMPlexGetConeSize(dm,p,&csize);
235:   if (label) {
236:     DMLabelGetValue(label,p,mid);
237:   } else {
238:     *mid = 1;
239:   }
240:   *cid = mfem_table_cid[depth][csize];
241:   return(0);
242: }

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

250:   DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
251:   DMGetDimension(dm,&dim);
252:   sdim = dim;
253:   if (csec) {
254:     DMGetCoordinateDim(dm,&sdim);
255:     PetscSectionGetOffset(csec,vStart,&off);
256:     off  = off/sdim;
257:     PetscSectionGetDof(csec,p,&dof);
258:   }
259:   if (!dof) {
260:     DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);
261:     for (i=0,q=0;i<numPoints*2;i+= 2)
262:       if ((points[i] >= vStart) && (points[i] < vEnd))
263:         vids[q++] = points[i]-vStart+off;
264:     DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);
265:   } else {
266:     PetscSectionGetOffset(csec,p,&off);
267:     PetscSectionGetDof(csec,p,&dof);
268:     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
269:   }
270:   DMPlexInvertCell(dim,q,vids);
271:   *nv  = q;
272:   return(0);
273: }

275: /* ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR */
276: static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
277: {
278:   DMLabel              label;
279:   PetscSection         coordSection,parentSection;
280:   Vec                  coordinates;
281:   IS                   globalNum = NULL;
282:   const PetscScalar    *array;
283:   const PetscInt       *gNum;
284:   PetscInt             bf,p,sdim,dim,depth,novl;
285:   PetscInt             cStart,cEnd,cEndInterior,vStart,vEnd,nvert;
286:   PetscMPIInt          commsize;
287:   PetscBool            localized,ovl,isascii,enable_boundary,enable_ncmesh;
288:   PetscBT              pown;
289:   PetscErrorCode       ierr;
290:   PetscContainer       glvis_container;
291:   PetscBool            enabled = PETSC_TRUE;

296:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
297:   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
298:   MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);
299:   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
300:   DMGetDimension(dm,&dim);
301:   DMPlexGetDepth(dm,&depth);
302:   if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");

304:   /* get container: determines if a process visualizes is portion of the data or not */
305:   PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
306:   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
307:   {
308:     PetscViewerGLVisInfo glvis_info;
309:     PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
310:     enabled = glvis_info->enabled;
311:   }
312:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
313:   DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
314:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
315:   DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
316:   DMGetCoordinatesLocalized(dm,&localized);
317:   DMGetCoordinateSection(dm,&coordSection);

319:   /*
320:      a couple of sections of the mesh specification are disabled
321:        - boundary: unless we want to visualize boundary attributes, the boundary is not needed for proper mesh visualization
322:        - 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
323:   */
324:   enable_boundary = PETSC_FALSE;
325:   enable_ncmesh   = PETSC_FALSE;
326:   PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");
327:   PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation; useful for debugging purposes",NULL,enable_boundary,&enable_boundary,NULL);
328:   PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation; useful for debugging purposes",NULL,enable_ncmesh,&enable_ncmesh,NULL);
329:   PetscOptionsEnd();

331:   /* Identify possible cells in the overlap */
332:   gNum = NULL;
333:   novl = 0;
334:   ovl  = PETSC_FALSE;
335:   pown = NULL;
336:   MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);
337:   if (commsize > 1) {
338:     DMPlexGetCellNumbering(dm,&globalNum);
339:     ISGetIndices(globalNum,&gNum);
340:     for (p=cStart; p<cEnd; p++) {
341:       if (gNum[p-cStart] < 0) {
342:         ovl = PETSC_TRUE;
343:         novl++;
344:       }
345:     }
346:     if (ovl) {
347:       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
348:          TODO: garbage collector? attach pown to dm?  */
349:       PetscBTCreate(PetscMax(cEnd-cStart,vEnd-vStart),&pown);
350:     }
351:   }

353:   if (!enabled) {
354:     PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
355:     PetscViewerASCIIPrintf(viewer,"\ndimension\n");
356:     PetscViewerASCIIPrintf(viewer,"%D\n",dim);
357:     PetscViewerASCIIPrintf(viewer,"\nelements\n");
358:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
359:     PetscViewerASCIIPrintf(viewer,"\nboundary\n");
360:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
361:     DMGetCoordinateDim(dm,&sdim);
362:     PetscViewerASCIIPrintf(viewer,"\nvertices\n");
363:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
364:     PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
365:     PetscBTDestroy(&pown);
366:     if (globalNum) {
367:       ISRestoreIndices(globalNum,&gNum);
368:     }
369:     return(0);
370:   }

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

375:   /* topological dimension */
376:   PetscViewerASCIIPrintf(viewer,"\ndimension\n");
377:   PetscViewerASCIIPrintf(viewer,"%D\n",dim);

379:   /* elements */
380:   /* TODO: is this the label we want to use for marking material IDs?
381:      We should decide to have a single marker for these stuff
382:      Proposal: DMSetMaterialLabel?
383:   */
384:   DMGetLabel(dm,"Cell Sets",&label);
385:   PetscViewerASCIIPrintf(viewer,"\nelements\n");
386:   PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);
387:   for (p=cStart;p<cEnd;p++) {
388:     int      vids[8];
389:     PetscInt i,nv = 0,cid = -1,mid = 1;

391:     if (ovl) {
392:       if (gNum[p-cStart] < 0) continue;
393:       else {
394:         PetscBTSet(pown,p-cStart);
395:       }
396:     }
397:     DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);
398:     DMPlexGetPointMFEMVertexIDs_Internal(dm,p,localized ? coordSection : NULL,&nv,vids);
399:     PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
400:     for (i=0;i<nv;i++) {
401:       PetscViewerASCIIPrintf(viewer," %d",vids[i]);
402:     }
403:     PetscViewerASCIIPrintf(viewer,"\n");
404:   }

406:   /* boundary */
407:   PetscViewerASCIIPrintf(viewer,"\nboundary\n");
408:   if (!enable_boundary) {
409:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
410:   } else {
411:     PetscInt fStart,fEnd,fEndInterior;
412:     PetscInt *faces = NULL,fpc = 0,vpf = 0;
413:     PetscInt fv1[]     = {0,1},
414:              fv2tri[]  = {0,1,
415:                           1,2,
416:                           2,0},
417:              fv2quad[] = {0,1,
418:                           1,2,
419:                           2,3,
420:                           3,0},
421:              fv3tet[]  = {0,1,2,
422:                           0,3,1,
423:                           0,2,3,
424:                           2,1,3},
425:              fv3hex[]  = {0,1,2,3,
426:                        4,5,6,7,
427:                        0,3,5,4,
428:                        2,1,7,6,
429:                        3,2,6,5,
430:                        0,4,7,1};

432:     if (localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Boundary specification with periodic mesh");
433:     /* determine orientation of boundary mesh */
434:     if (cEnd-cStart) {
435:       DMPlexGetConeSize(dm,cStart,&fpc);
436:       switch(dim) {
437:         case 1:
438:           if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
439:           faces = fv1;
440:           vpf = 1;
441:           break;
442:         case 2:
443:           switch (fpc) {
444:             case 3:
445:               faces = fv2tri;
446:               vpf   = 2;
447:               break;
448:             case 4:
449:               faces = fv2quad;
450:               vpf   = 2;
451:               break;
452:             default:
453:               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
454:               break;
455:           }
456:           break;
457:         case 3:
458:           switch (fpc) {
459:             case 4:
460:               faces = fv3tet;
461:               vpf   = 3;
462:               break;
463:             case 6:
464:               faces = fv3hex;
465:               vpf   = 4;
466:               break;
467:             default:
468:               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
469:               break;
470:           }
471:           break;
472:         default:
473:           SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
474:           break;
475:       }
476:     }

478:     DMPlexGetHybridBounds(dm, NULL, &fEndInterior, NULL, NULL);
479:     DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);
480:     fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
481:     if (!ovl) {
482:       for (p=fStart,bf=0;p<fEnd;p++) {
483:         PetscInt supportSize;

485:         DMPlexGetSupportSize(dm,p,&supportSize);
486:         if (supportSize == 1) bf++;
487:       }
488:     } else {
489:       for (p=fStart,bf=0;p<fEnd;p++) {
490:         const PetscInt *support;
491:         PetscInt       i,supportSize;
492:         PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;

494:         DMPlexGetSupportSize(dm,p,&supportSize);
495:         DMPlexGetSupport(dm,p,&support);
496:         for (i=0;i<supportSize;i++) {
497:           if (PetscBTLookup(pown,support[i]-cStart)) has_owned = PETSC_TRUE;
498:           else has_ghost = PETSC_TRUE;
499:         }
500:         if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) bf++;
501:       }
502:     }
503:     /* TODO: is this the label we want to use for marking boundary facets?
504:        We should decide to have a single marker for these stuff
505:        Proposal: DMSetBoundaryLabel?
506:     */
507:     DMGetLabel(dm,"marker",&label);
508:     PetscViewerASCIIPrintf(viewer,"%D\n",bf);
509:     if (!ovl) {
510:       for (p=fStart;p<fEnd;p++) {
511:         PetscInt supportSize;

513:         DMPlexGetSupportSize(dm,p,&supportSize);
514:         if (supportSize == 1) {
515:           const    PetscInt *support, *cone;
516:           PetscInt i,c,v,cid = -1,mid = -1;
517:           PetscInt cellClosureSize,*cellClosure = NULL;

519:           DMPlexGetSupport(dm,p,&support);
520:           DMPlexGetCone(dm,support[0],&cone);
521:           for (c=0;c<fpc;c++)
522:             if (cone[c] == p)
523:               break;

525:           DMPlexGetTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);
526:           for (v=0;v<cellClosureSize;v++)
527:             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
528:               break;
529:           DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);
530:           PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
531:           for (i=0;i<vpf;i++) {
532:             PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);
533:           }
534:           PetscViewerASCIIPrintf(viewer,"\n");
535:           DMPlexRestoreTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);
536:         }
537:       }
538:     } else {
539:       for (p=fStart;p<fEnd;p++) {
540:         const PetscInt *support;
541:         PetscInt       i,supportSize,validcell = -1;
542:         PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;

544:         DMPlexGetSupportSize(dm,p,&supportSize);
545:         DMPlexGetSupport(dm,p,&support);
546:         for (i=0;i<supportSize;i++) {
547:           if (PetscBTLookup(pown,support[i]-cStart)) {
548:             has_owned = PETSC_TRUE;
549:             validcell = support[i];
550:           } else has_ghost = PETSC_TRUE;
551:         }
552:         if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) {
553:           const    PetscInt *support, *cone;
554:           PetscInt i,c,v,cid = -1,mid = -1;
555:           PetscInt cellClosureSize,*cellClosure = NULL;

557:           DMPlexGetSupport(dm,p,&support);
558:           DMPlexGetCone(dm,validcell,&cone);
559:           for (c=0;c<fpc;c++)
560:             if (cone[c] == p)
561:               break;

563:           DMPlexGetTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);
564:           for (v=0;v<cellClosureSize;v++)
565:             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
566:               break;
567:           DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);
568:           PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
569:           for (i=0;i<vpf;i++) {
570:             PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);
571:           }
572:           PetscViewerASCIIPrintf(viewer,"\n");
573:           DMPlexRestoreTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);
574:         }
575:       }
576:     }
577:   }

579:   /* mark owned vertices */
580:   if (ovl) {
581:     PetscBTMemzero(vEnd-vStart,pown);
582:     for (p=cStart;p<cEnd;p++) {
583:       PetscInt i,closureSize,*closure = NULL;

585:       if (gNum[p-cStart] < 0) continue;
586:       DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
587:       for (i=0;i<closureSize;i++) {
588:         const PetscInt pp = closure[2*i];

590:         if (pp >= vStart && pp < vEnd) {
591:           PetscBTSet(pown,pp-vStart);
592:         }
593:       }
594:       DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
595:     }
596:   }
597:   if (globalNum) {
598:     ISRestoreIndices(globalNum,&gNum);
599:   }

601:   /* vertex_parents (Non-conforming meshes) */
602:   parentSection  = NULL;
603:   if (enable_ncmesh) {
604:     DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);
605:   }
606:   if (parentSection) {
607:     PetscInt vp,gvp;

609:     for (vp=0,p=vStart;p<vEnd;p++) {
610:       DMLabel  dlabel;
611:       PetscInt parent,depth;

613:       if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
614:       DMPlexGetDepthLabel(dm,&dlabel);
615:       DMLabelGetValue(dlabel,p,&depth);
616:       DMPlexGetTreeParent(dm,p,&parent,NULL);
617:       if (parent != p) vp++;
618:     }
619:     MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
620:     if (gvp) {
621:       PetscInt  maxsupp;
622:       PetscBool *skip = NULL;

624:       PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");
625:       PetscViewerASCIIPrintf(viewer,"%D\n",vp);
626:       DMPlexGetMaxSizes(dm,NULL,&maxsupp);
627:       PetscMalloc1(maxsupp,&skip);
628:       for (p=vStart;p<vEnd;p++) {
629:         DMLabel  dlabel;
630:         PetscInt parent;

632:         if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
633:         DMPlexGetDepthLabel(dm,&dlabel);
634:         DMPlexGetTreeParent(dm,p,&parent,NULL);
635:         if (parent != p) {
636:           int            vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
637:           PetscInt       i,nv,size,n,numChildren,depth = -1;
638:           const PetscInt *children;
639:           DMPlexGetConeSize(dm,parent,&size);
640:           switch (size) {
641:             case 2: /* edge */
642:               nv   = 0;
643:               DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);
644:               PetscViewerASCIIPrintf(viewer,"%D",p-vStart);
645:               for (i=0;i<nv;i++) {
646:                 PetscViewerASCIIPrintf(viewer," %d",vids[i]);
647:               }
648:               PetscViewerASCIIPrintf(viewer,"\n");
649:               vp--;
650:               break;
651:             case 4: /* face */
652:               DMPlexGetTreeChildren(dm,parent,&numChildren,&children);
653:               for (n=0;n<numChildren;n++) {
654:                 DMLabelGetValue(dlabel,children[n],&depth);
655:                 if (!depth) {
656:                   const PetscInt *hvsupp,*hesupp,*cone;
657:                   PetscInt       hvsuppSize,hesuppSize,coneSize;
658:                   PetscInt       hv = children[n],he = -1,f;

660:                   PetscMemzero(skip,maxsupp*sizeof(PetscBool));
661:                   DMPlexGetSupportSize(dm,hv,&hvsuppSize);
662:                   DMPlexGetSupport(dm,hv,&hvsupp);
663:                   for (i=0;i<hvsuppSize;i++) {
664:                     PetscInt ep;
665:                     DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);
666:                     if (ep != hvsupp[i]) {
667:                       he = hvsupp[i];
668:                     } else {
669:                       skip[i] = PETSC_TRUE;
670:                     }
671:                   }
672:                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
673:                   DMPlexGetCone(dm,he,&cone);
674:                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
675:                   DMPlexGetSupportSize(dm,he,&hesuppSize);
676:                   DMPlexGetSupport(dm,he,&hesupp);
677:                   for (f=0;f<hesuppSize;f++) {
678:                     PetscInt j;

680:                     DMPlexGetCone(dm,hesupp[f],&cone);
681:                     DMPlexGetConeSize(dm,hesupp[f],&coneSize);
682:                     for (j=0;j<coneSize;j++) {
683:                       PetscInt k;
684:                       for (k=0;k<hvsuppSize;k++) {
685:                         if (hvsupp[k] == cone[j]) {
686:                           skip[k] = PETSC_TRUE;
687:                           break;
688:                         }
689:                       }
690:                     }
691:                   }
692:                   for (i=0;i<hvsuppSize;i++) {
693:                     if (!skip[i]) {
694:                       DMPlexGetCone(dm,hvsupp[i],&cone);
695:                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
696:                     }
697:                   }
698:                   PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);
699:                   for (i=0;i<2;i++) {
700:                     PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);
701:                   }
702:                   PetscViewerASCIIPrintf(viewer,"\n");
703:                   vp--;
704:                 }
705:               }
706:               break;
707:             default:
708:               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size);
709:           }
710:         }
711:       }
712:       PetscFree(skip);
713:     }
714:     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
715:   }
716:   PetscBTDestroy(&pown);

718:   /* vertices */
719:   DMGetCoordinateDim(dm,&sdim);
720:   DMGetCoordinatesLocal(dm,&coordinates);
721:   VecGetLocalSize(coordinates,&nvert);
722:   PetscViewerASCIIPrintf(viewer,"\nvertices\n");
723:   PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);
724:   PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
725:   VecGetArrayRead(coordinates,&array);
726:   for (p=0;p<nvert/sdim;p++) {
727:     PetscInt s;
728:     for (s=0;s<sdim;s++) {
729:       PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[p*sdim+s]));
730:     }
731:     PetscViewerASCIIPrintf(viewer,"\n");
732:   }
733:   VecRestoreArrayRead(coordinates,&array);
734:   return(0);
735: }

737: /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
738: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
739: {
741:   PetscBool      isglvis,isascii;

746:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);
747:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
748:   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
749:   if (isglvis) {
750:     PetscViewer          view;
751:     PetscViewerGLVisType type;

753:     PetscViewerGLVisGetType_Private(viewer,&type);
754:     PetscViewerGLVisGetDMWindow_Private(viewer,&view);
755:     if (view) { /* in the socket case, it may happen that the connection failed */
756:       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
757:         PetscMPIInt size,rank;
758:         MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
759:         MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
760:         PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);
761:       }
762:       DMPlexView_GLVis_ASCII(dm,view);
763:       PetscViewerFlush(view);
764:       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
765:         PetscInt    dim;
766:         const char* name;

768:         PetscObjectGetName((PetscObject)dm,&name);
769:         DMGetDimension(dm,&dim);
770:         PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);
771:         PetscBarrier((PetscObject)dm);
772:       }
773:     }
774:   } else {
775:     DMPlexView_GLVis_ASCII(dm,viewer);
776:   }
777:   return(0);
778: }