Actual source code: plexglvis.c
petsc-3.14.6 2021-03-30
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 (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim);
260: if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p);
261: if (depth != 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p);
262: *cid = mfem_table_cid_unint[dim][csize];
263: } else {
264: if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p);
265: if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p);
266: *cid = mfem_table_cid[pdepth][csize];
267: }
268: return(0);
269: }
271: static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[])
272: {
273: PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;
277: DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
278: DMGetDimension(dm,&dim);
279: sdim = dim;
280: if (csec) {
281: PetscInt sStart,sEnd;
283: DMGetCoordinateDim(dm,&sdim);
284: PetscSectionGetChart(csec,&sStart,&sEnd);
285: PetscSectionGetOffset(csec,vStart,&off);
286: off = off/sdim;
287: if (p >= sStart && p < sEnd) {
288: PetscSectionGetDof(csec,p,&dof);
289: }
290: }
291: if (!dof) {
292: DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);
293: for (i=0,q=0;i<numPoints*2;i+= 2)
294: if ((points[i] >= vStart) && (points[i] < vEnd))
295: vids[q++] = points[i]-vStart+off;
296: DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);
297: } else {
298: PetscSectionGetOffset(csec,p,&off);
299: PetscSectionGetDof(csec,p,&dof);
300: for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
301: }
302: *nv = q;
303: return(0);
304: }
306: static PetscErrorCode GLVisCreateFE(PetscFE femIn,char name[32],PetscFE *fem)
307: {
308: DM K;
309: PetscSpace P;
310: PetscDualSpace Q;
311: PetscQuadrature q,fq;
312: PetscInt dim,deg,dof;
313: DMPolytopeType ptype;
314: PetscBool isSimplex,isTensor;
315: PetscBool continuity = PETSC_FALSE;
316: PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI;
317: PetscBool endpoint = PETSC_TRUE;
318: MPI_Comm comm;
319: PetscErrorCode ierr;
322: comm = PetscObjectComm((PetscObject)femIn);
323: PetscFEGetBasisSpace(femIn,&P);
324: PetscFEGetDualSpace(femIn,&Q);
325: PetscDualSpaceGetDM(Q,&K);
326: DMGetDimension(K,&dim);
327: PetscSpaceGetDegree(P,°,NULL);
328: PetscSpaceGetNumComponents(P,&dof);
329: DMPlexGetCellType(K,0,&ptype);
330: switch (ptype) {
331: case DM_POLYTOPE_QUADRILATERAL:
332: case DM_POLYTOPE_HEXAHEDRON:
333: isSimplex = PETSC_FALSE; break;
334: default:
335: isSimplex = PETSC_TRUE; break;
336: }
337: isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
338: /* Create space */
339: PetscSpaceCreate(comm,&P);
340: PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);
341: PetscSpacePolynomialSetTensor(P,isTensor);
342: PetscSpaceSetNumComponents(P,dof);
343: PetscSpaceSetNumVariables(P,dim);
344: PetscSpaceSetDegree(P,deg,PETSC_DETERMINE);
345: PetscSpaceSetUp(P);
346: /* Create dual space */
347: PetscDualSpaceCreate(comm,&Q);
348: PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);
349: PetscDualSpaceLagrangeSetTensor(Q,isTensor);
350: PetscDualSpaceLagrangeSetContinuity(Q,continuity);
351: PetscDualSpaceLagrangeSetNodeType(Q,nodeType,endpoint,0);
352: PetscDualSpaceSetNumComponents(Q,dof);
353: PetscDualSpaceSetOrder(Q,deg);
354: PetscDualSpaceCreateReferenceCell(Q,dim,isSimplex,&K);
355: PetscDualSpaceSetDM(Q,K);
356: DMDestroy(&K);
357: PetscDualSpaceSetUp(Q);
358: /* Create quadrature */
359: if (isSimplex) {
360: PetscDTStroudConicalQuadrature(dim, 1,deg+1,-1,+1,&q);
361: PetscDTStroudConicalQuadrature(dim-1,1,deg+1,-1,+1,&fq);
362: } else {
363: PetscDTGaussTensorQuadrature(dim, 1,deg+1,-1,+1,&q);
364: PetscDTGaussTensorQuadrature(dim-1,1,deg+1,-1,+1,&fq);
365: }
366: /* Create finite element */
367: PetscFECreate(comm,fem);
368: PetscSNPrintf(name,32,"L2_T1_%DD_P%D",dim,deg);
369: PetscObjectSetName((PetscObject)*fem,name);
370: PetscFESetType(*fem,PETSCFEBASIC);
371: PetscFESetNumComponents(*fem,dof);
372: PetscFESetBasisSpace(*fem,P);
373: PetscFESetDualSpace(*fem,Q);
374: PetscFESetQuadrature(*fem,q);
375: PetscFESetFaceQuadrature(*fem,fq);
376: PetscFESetUp(*fem);
377: /* Cleanup */
378: PetscSpaceDestroy(&P);
379: PetscDualSpaceDestroy(&Q);
380: PetscQuadratureDestroy(&q);
381: PetscQuadratureDestroy(&fq);
382: return(0);
383: }
385: /*
386: ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
387: Higher order meshes are also supported
388: */
389: static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
390: {
391: DMLabel label;
392: PetscSection coordSection,parentSection;
393: Vec coordinates,hovec;
394: const PetscScalar *array;
395: PetscInt bf,p,sdim,dim,depth,novl,minl;
396: PetscInt cStart,cEnd,vStart,vEnd,nvert;
397: PetscMPIInt size;
398: PetscBool localized,isascii;
399: PetscBool enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE;
400: PetscBT pown,vown;
401: PetscErrorCode ierr;
402: PetscContainer glvis_container;
403: PetscBool cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
404: PetscBool enable_emark,enable_bmark;
405: const char *fmt;
406: char emark[64] = "",bmark[64] = "";
411: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
412: if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
413: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
414: if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
415: DMGetDimension(dm,&dim);
417: /* get container: determines if a process visualizes is portion of the data or not */
418: PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
419: if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
420: {
421: PetscViewerGLVisInfo glvis_info;
422: PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
423: enabled = glvis_info->enabled;
424: fmt = glvis_info->fmt;
425: }
427: /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
428: DMPlex does not currently support HO meshes, so there's no API for this */
429: PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);
430: PetscObjectReference((PetscObject)hovec);
431: if (!hovec) {
432: DM cdm;
433: PetscFE disc;
434: PetscClassId classid;
436: DMGetCoordinateDM(dm,&cdm);
437: DMGetField(cdm,0,NULL,(PetscObject*)&disc);
438: PetscObjectGetClassId((PetscObject)disc,&classid);
439: if (classid == PETSCFE_CLASSID) {
440: DM hocdm;
441: PetscFE hodisc;
442: Vec vec;
443: Mat mat;
444: char name[32],fec_type[64];
446: GLVisCreateFE(disc,name,&hodisc);
447: DMClone(cdm,&hocdm);
448: DMSetField(hocdm,0,NULL,(PetscObject)hodisc);
449: PetscFEDestroy(&hodisc);
450: DMCreateDS(hocdm);
452: DMGetCoordinates(dm,&vec);
453: DMCreateGlobalVector(hocdm,&hovec);
454: DMCreateInterpolation(cdm,hocdm,&mat,NULL);
455: MatInterpolate(mat,vec,hovec);
456: MatDestroy(&mat);
457: DMDestroy(&hocdm);
459: PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name);
460: PetscObjectSetName((PetscObject)hovec,fec_type);
461: }
462: }
464: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
465: DMPlexGetGhostCellStratum(dm,&p,NULL);
466: if (p >= 0) cEnd = p;
467: DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);
468: DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);
469: DMGetCoordinatesLocalized(dm,&localized);
470: if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
471: DMGetCoordinateSection(dm,&coordSection);
472: DMGetCoordinateDim(dm,&sdim);
473: DMGetCoordinatesLocal(dm,&coordinates);
474: if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");
476: /*
477: a couple of sections of the mesh specification are disabled
478: - 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)
479: - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
480: and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
481: */
482: enable_boundary = PETSC_FALSE;
483: enable_ncmesh = PETSC_FALSE;
484: enable_mfem = PETSC_FALSE;
485: enable_emark = PETSC_FALSE;
486: enable_bmark = PETSC_FALSE;
487: /* I'm tired of problems with negative values in the markers, disable them */
488: PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");
489: PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);
490: PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL);
491: PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL);
492: PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL);
493: PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);
494: PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);
495: PetscOptionsEnd();
496: if (enable_bmark) enable_boundary = PETSC_TRUE;
498: MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
499: if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
500: DMPlexGetDepth(dm,&depth);
501: if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
502: "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
503: if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
504: "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
505: if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
506: if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth);
507: cellvertex = PETSC_TRUE;
508: }
510: /* Identify possible cells in the overlap */
511: novl = 0;
512: pown = NULL;
513: if (size > 1) {
514: IS globalNum = NULL;
515: const PetscInt *gNum;
516: PetscBool ovl = PETSC_FALSE;
518: PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);
519: if (!globalNum) {
520: if (view_ovl) {
521: ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum);
522: } else {
523: DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);
524: }
525: PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);
526: PetscObjectDereference((PetscObject)globalNum);
527: }
528: ISGetIndices(globalNum,&gNum);
529: for (p=cStart; p<cEnd; p++) {
530: if (gNum[p-cStart] < 0) {
531: ovl = PETSC_TRUE;
532: novl++;
533: }
534: }
535: if (ovl) {
536: /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
537: TODO: garbage collector? attach pown to dm? */
538: PetscBTCreate(cEnd-cStart,&pown);
539: for (p=cStart; p<cEnd; p++) {
540: if (gNum[p-cStart] < 0) continue;
541: else {
542: PetscBTSet(pown,p-cStart);
543: }
544: }
545: }
546: ISRestoreIndices(globalNum,&gNum);
547: }
549: /* return if this process is disabled */
550: if (!enabled) {
551: PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
552: PetscViewerASCIIPrintf(viewer,"\ndimension\n");
553: PetscViewerASCIIPrintf(viewer,"%D\n",dim);
554: PetscViewerASCIIPrintf(viewer,"\nelements\n");
555: PetscViewerASCIIPrintf(viewer,"%D\n",0);
556: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
557: PetscViewerASCIIPrintf(viewer,"%D\n",0);
558: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
559: PetscViewerASCIIPrintf(viewer,"%D\n",0);
560: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
561: PetscBTDestroy(&pown);
562: VecDestroy(&hovec);
563: return(0);
564: }
566: if (enable_mfem) {
567: if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
568: PetscInt vpc = 0;
569: char fec[64];
570: PetscInt vids[8] = {0,1,2,3,4,5,6,7};
571: PetscInt hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3};
572: PetscInt quadv[8] = {0,1,3,2}, triv[3] = {0,1,2};
573: PetscInt *dof = NULL;
574: PetscScalar *array,*ptr;
576: PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);
577: if (cEnd-cStart) {
578: PetscInt fpc;
580: DMPlexGetConeSize(dm,cStart,&fpc);
581: switch(dim) {
582: case 1:
583: vpc = 2;
584: dof = hexv;
585: break;
586: case 2:
587: switch (fpc) {
588: case 3:
589: vpc = 3;
590: dof = triv;
591: break;
592: case 4:
593: vpc = 4;
594: dof = quadv;
595: break;
596: default:
597: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
598: }
599: break;
600: case 3:
601: switch (fpc) {
602: case 4: /* TODO: still need to understand L2 ordering for tets */
603: vpc = 4;
604: dof = tetv;
605: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
606: case 6:
607: if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc);
608: vpc = 8;
609: dof = hexv;
610: break;
611: case 8:
612: if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
613: vpc = 8;
614: dof = hexv;
615: break;
616: default:
617: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
618: }
619: break;
620: default:
621: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
622: }
623: DMPlexReorderCell(dm,cStart,vids);
624: }
625: if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
626: VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);
627: PetscObjectSetName((PetscObject)hovec,fec);
628: VecGetArray(hovec,&array);
629: ptr = array;
630: for (p=cStart;p<cEnd;p++) {
631: PetscInt csize,v,d;
632: PetscScalar *vals = NULL;
634: if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
635: DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);
636: 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);
637: for (v=0;v<vpc;v++) {
638: for (d=0;d<sdim;d++) {
639: ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
640: }
641: }
642: ptr += vpc*sdim;
643: DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);
644: }
645: VecRestoreArray(hovec,&array);
646: }
647: }
648: /* if we have high-order coordinates in 3D, we need to specify the boundary */
649: if (hovec && dim == 3) enable_boundary = PETSC_TRUE;
651: /* header */
652: PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
654: /* topological dimension */
655: PetscViewerASCIIPrintf(viewer,"\ndimension\n");
656: PetscViewerASCIIPrintf(viewer,"%D\n",dim);
658: /* elements */
659: minl = 1;
660: label = NULL;
661: if (enable_emark) {
662: PetscInt lminl = PETSC_MAX_INT;
664: DMGetLabel(dm,emark,&label);
665: if (label) {
666: IS vals;
667: PetscInt ldef;
669: DMLabelGetDefaultValue(label,&ldef);
670: DMLabelGetValueIS(label,&vals);
671: ISGetMinMax(vals,&lminl,NULL);
672: ISDestroy(&vals);
673: lminl = PetscMin(ldef,lminl);
674: }
675: MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));
676: if (minl == PETSC_MAX_INT) minl = 1;
677: }
678: PetscViewerASCIIPrintf(viewer,"\nelements\n");
679: PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);
680: for (p=cStart;p<cEnd;p++) {
681: PetscInt vids[8];
682: PetscInt i,nv = 0,cid = -1,mid = 1;
684: if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
685: DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);
686: DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);
687: DMPlexReorderCell(dm,p,vids);
688: PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
689: for (i=0;i<nv;i++) {
690: PetscViewerASCIIPrintf(viewer," %D",vids[i]);
691: }
692: PetscViewerASCIIPrintf(viewer,"\n");
693: }
695: /* boundary */
696: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
697: if (!enable_boundary) {
698: PetscViewerASCIIPrintf(viewer,"%D\n",0);
699: } else {
700: DMLabel perLabel;
701: PetscBT bfaces;
702: PetscInt fStart,fEnd,*fcells;
704: DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);
705: PetscBTCreate(fEnd-fStart,&bfaces);
706: DMPlexGetMaxSizes(dm,NULL,&p);
707: PetscMalloc1(p,&fcells);
708: DMGetLabel(dm,"glvis_periodic_cut",&perLabel);
709: if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
710: DMCreateLabel(dm,"glvis_periodic_cut");
711: DMGetLabel(dm,"glvis_periodic_cut",&perLabel);
712: DMLabelSetDefaultValue(perLabel,1);
713: for (p=cStart;p<cEnd;p++) {
714: DMPolytopeType cellType;
715: PetscInt dof;
717: DMPlexGetCellType(dm,p,&cellType);
718: PetscSectionGetDof(coordSection,p,&dof);
719: if (dof) {
720: PetscInt uvpc, v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
721: PetscScalar *vals = NULL;
723: uvpc = DMPolytopeTypeGetNumVertices(cellType);
724: if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
725: DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);
726: DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);
727: for (v=0;v<cellClosureSize;v++)
728: if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
729: vidxs = cellClosure + 2*v;
730: break;
731: }
732: if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
733: for (v=0;v<uvpc;v++) {
734: PetscInt s;
736: for (s=0;s<sdim;s++) {
737: if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) {
738: DMLabelSetValue(perLabel,vidxs[2*v],2);
739: }
740: }
741: }
742: DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);
743: DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);
744: }
745: }
746: if (dim > 1) {
747: PetscInt eEnd,eStart;
749: DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);
750: for (p=eStart;p<eEnd;p++) {
751: const PetscInt *cone;
752: PetscInt coneSize,i;
753: PetscBool ispe = PETSC_TRUE;
755: DMPlexGetCone(dm,p,&cone);
756: DMPlexGetConeSize(dm,p,&coneSize);
757: for (i=0;i<coneSize;i++) {
758: PetscInt v;
760: DMLabelGetValue(perLabel,cone[i],&v);
761: ispe = (PetscBool)(ispe && (v==2));
762: }
763: if (ispe && coneSize) {
764: PetscInt ch, numChildren;
765: const PetscInt *children;
767: DMLabelSetValue(perLabel,p,2);
768: DMPlexGetTreeChildren(dm,p,&numChildren,&children);
769: for (ch = 0; ch < numChildren; ch++) {
770: DMLabelSetValue(perLabel,children[ch],2);
771: }
772: }
773: }
774: if (dim > 2) {
775: for (p=fStart;p<fEnd;p++) {
776: const PetscInt *cone;
777: PetscInt coneSize,i;
778: PetscBool ispe = PETSC_TRUE;
780: DMPlexGetCone(dm,p,&cone);
781: DMPlexGetConeSize(dm,p,&coneSize);
782: for (i=0;i<coneSize;i++) {
783: PetscInt v;
785: DMLabelGetValue(perLabel,cone[i],&v);
786: ispe = (PetscBool)(ispe && (v==2));
787: }
788: if (ispe && coneSize) {
789: PetscInt ch, numChildren;
790: const PetscInt *children;
792: DMLabelSetValue(perLabel,p,2);
793: DMPlexGetTreeChildren(dm,p,&numChildren,&children);
794: for (ch = 0; ch < numChildren; ch++) {
795: DMLabelSetValue(perLabel,children[ch],2);
796: }
797: }
798: }
799: }
800: }
801: }
802: for (p=fStart;p<fEnd;p++) {
803: const PetscInt *support;
804: PetscInt supportSize;
805: PetscBool isbf = PETSC_FALSE;
807: DMPlexGetSupportSize(dm,p,&supportSize);
808: if (pown) {
809: PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
810: PetscInt i;
812: DMPlexGetSupport(dm,p,&support);
813: for (i=0;i<supportSize;i++) {
814: if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
815: else has_ghost = PETSC_TRUE;
816: }
817: isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
818: } else {
819: isbf = (PetscBool)(supportSize == 1);
820: }
821: if (!isbf && perLabel) {
822: const PetscInt *cone;
823: PetscInt coneSize,i;
825: DMPlexGetCone(dm,p,&cone);
826: DMPlexGetConeSize(dm,p,&coneSize);
827: isbf = PETSC_TRUE;
828: for (i=0;i<coneSize;i++) {
829: PetscInt v,d;
831: DMLabelGetValue(perLabel,cone[i],&v);
832: DMLabelGetDefaultValue(perLabel,&d);
833: isbf = (PetscBool)(isbf && v != d);
834: }
835: }
836: if (isbf) {
837: PetscBTSet(bfaces,p-fStart);
838: }
839: }
840: /* count boundary faces */
841: for (p=fStart,bf=0;p<fEnd;p++) {
842: if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
843: const PetscInt *support;
844: PetscInt supportSize,c;
846: DMPlexGetSupportSize(dm,p,&supportSize);
847: DMPlexGetSupport(dm,p,&support);
848: for (c=0;c<supportSize;c++) {
849: const PetscInt *cone;
850: PetscInt cell,cl,coneSize;
852: cell = support[c];
853: if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue;
854: DMPlexGetCone(dm,cell,&cone);
855: DMPlexGetConeSize(dm,cell,&coneSize);
856: for (cl=0;cl<coneSize;cl++) {
857: if (cone[cl] == p) {
858: bf += 1;
859: break;
860: }
861: }
862: }
863: }
864: }
865: minl = 1;
866: label = NULL;
867: if (enable_bmark) {
868: PetscInt lminl = PETSC_MAX_INT;
870: DMGetLabel(dm,bmark,&label);
871: if (label) {
872: IS vals;
873: PetscInt ldef;
875: DMLabelGetDefaultValue(label,&ldef);
876: DMLabelGetValueIS(label,&vals);
877: ISGetMinMax(vals,&lminl,NULL);
878: ISDestroy(&vals);
879: lminl = PetscMin(ldef,lminl);
880: }
881: MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));
882: if (minl == PETSC_MAX_INT) minl = 1;
883: }
884: PetscViewerASCIIPrintf(viewer,"%D\n",bf);
885: for (p=fStart;p<fEnd;p++) {
886: if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
887: const PetscInt *support;
888: PetscInt supportSize,c,nc = 0;
890: DMPlexGetSupportSize(dm,p,&supportSize);
891: DMPlexGetSupport(dm,p,&support);
892: if (pown) {
893: for (c=0;c<supportSize;c++) {
894: if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
895: fcells[nc++] = support[c];
896: }
897: }
898: } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
899: for (c=0;c<nc;c++) {
900: const DMPolytopeType *faceTypes;
901: DMPolytopeType cellType;
902: const PetscInt *faceSizes,*cone;
903: PetscInt vids[8],*faces,st,i,coneSize,cell,cl,nv,cid = -1,mid = -1;
905: cell = fcells[c];
906: DMPlexGetCone(dm,cell,&cone);
907: DMPlexGetConeSize(dm,cell,&coneSize);
908: for (cl=0;cl<coneSize;cl++)
909: if (cone[cl] == p)
910: break;
911: if (cl == coneSize) continue;
913: /* face material id and type */
914: DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);
915: PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);
916: /* vertex ids */
917: DMPlexGetCellType(dm,cell,&cellType);
918: DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);
919: DMPlexGetRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);
920: st = 0;
921: for (i=0;i<cl;i++) st += faceSizes[i];
922: DMPlexInvertCell(faceTypes[cl],faces + st);
923: for (i=0;i<faceSizes[cl];i++) {
924: PetscViewerASCIIPrintf(viewer," %d",faces[st+i]);
925: }
926: PetscViewerASCIIPrintf(viewer,"\n");
927: DMPlexRestoreRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);
928: bf -= 1;
929: }
930: }
931: }
932: if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
933: PetscBTDestroy(&bfaces);
934: PetscFree(fcells);
935: }
937: /* mark owned vertices */
938: vown = NULL;
939: if (pown) {
940: PetscBTCreate(vEnd-vStart,&vown);
941: for (p=cStart;p<cEnd;p++) {
942: PetscInt i,closureSize,*closure = NULL;
944: if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
945: DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
946: for (i=0;i<closureSize;i++) {
947: const PetscInt pp = closure[2*i];
949: if (pp >= vStart && pp < vEnd) {
950: PetscBTSet(vown,pp-vStart);
951: }
952: }
953: DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
954: }
955: }
957: /* vertex_parents (Non-conforming meshes) */
958: parentSection = NULL;
959: if (enable_ncmesh) {
960: DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);
961: }
962: if (parentSection) {
963: PetscInt vp,gvp;
965: for (vp=0,p=vStart;p<vEnd;p++) {
966: DMLabel dlabel;
967: PetscInt parent,depth;
969: if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
970: DMPlexGetDepthLabel(dm,&dlabel);
971: DMLabelGetValue(dlabel,p,&depth);
972: DMPlexGetTreeParent(dm,p,&parent,NULL);
973: if (parent != p) vp++;
974: }
975: MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
976: if (gvp) {
977: PetscInt maxsupp;
978: PetscBool *skip = NULL;
980: PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");
981: PetscViewerASCIIPrintf(viewer,"%D\n",vp);
982: DMPlexGetMaxSizes(dm,NULL,&maxsupp);
983: PetscMalloc1(maxsupp,&skip);
984: for (p=vStart;p<vEnd;p++) {
985: DMLabel dlabel;
986: PetscInt parent;
988: if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
989: DMPlexGetDepthLabel(dm,&dlabel);
990: DMPlexGetTreeParent(dm,p,&parent,NULL);
991: if (parent != p) {
992: PetscInt vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
993: PetscInt i,nv,ssize,n,numChildren,depth = -1;
994: const PetscInt *children;
996: DMPlexGetConeSize(dm,parent,&ssize);
997: switch (ssize) {
998: case 2: /* edge */
999: nv = 0;
1000: DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);
1001: PetscViewerASCIIPrintf(viewer,"%D",p-vStart);
1002: for (i=0;i<nv;i++) {
1003: PetscViewerASCIIPrintf(viewer," %D",vids[i]);
1004: }
1005: PetscViewerASCIIPrintf(viewer,"\n");
1006: vp--;
1007: break;
1008: case 4: /* face */
1009: DMPlexGetTreeChildren(dm,parent,&numChildren,&children);
1010: for (n=0;n<numChildren;n++) {
1011: DMLabelGetValue(dlabel,children[n],&depth);
1012: if (!depth) {
1013: const PetscInt *hvsupp,*hesupp,*cone;
1014: PetscInt hvsuppSize,hesuppSize,coneSize;
1015: PetscInt hv = children[n],he = -1,f;
1017: PetscArrayzero(skip,maxsupp);
1018: DMPlexGetSupportSize(dm,hv,&hvsuppSize);
1019: DMPlexGetSupport(dm,hv,&hvsupp);
1020: for (i=0;i<hvsuppSize;i++) {
1021: PetscInt ep;
1022: DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);
1023: if (ep != hvsupp[i]) {
1024: he = hvsupp[i];
1025: } else {
1026: skip[i] = PETSC_TRUE;
1027: }
1028: }
1029: if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
1030: DMPlexGetCone(dm,he,&cone);
1031: vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
1032: DMPlexGetSupportSize(dm,he,&hesuppSize);
1033: DMPlexGetSupport(dm,he,&hesupp);
1034: for (f=0;f<hesuppSize;f++) {
1035: PetscInt j;
1037: DMPlexGetCone(dm,hesupp[f],&cone);
1038: DMPlexGetConeSize(dm,hesupp[f],&coneSize);
1039: for (j=0;j<coneSize;j++) {
1040: PetscInt k;
1041: for (k=0;k<hvsuppSize;k++) {
1042: if (hvsupp[k] == cone[j]) {
1043: skip[k] = PETSC_TRUE;
1044: break;
1045: }
1046: }
1047: }
1048: }
1049: for (i=0;i<hvsuppSize;i++) {
1050: if (!skip[i]) {
1051: DMPlexGetCone(dm,hvsupp[i],&cone);
1052: vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
1053: }
1054: }
1055: PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);
1056: for (i=0;i<2;i++) {
1057: PetscViewerASCIIPrintf(viewer," %D",vids[i]-vStart);
1058: }
1059: PetscViewerASCIIPrintf(viewer,"\n");
1060: vp--;
1061: }
1062: }
1063: break;
1064: default:
1065: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize);
1066: }
1067: }
1068: }
1069: PetscFree(skip);
1070: }
1071: if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
1072: }
1073: PetscBTDestroy(&pown);
1074: PetscBTDestroy(&vown);
1076: /* vertices */
1077: if (hovec) { /* higher-order meshes */
1078: const char *fec;
1079: PetscInt i,n,s;
1081: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
1082: PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);
1083: PetscViewerASCIIPrintf(viewer,"nodes\n");
1084: PetscObjectGetName((PetscObject)hovec,&fec);
1085: PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");
1086: PetscViewerASCIIPrintf(viewer,"%s\n",fec);
1087: PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);
1088: PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n"); /*Ordering::byVDIM*/
1089: VecGetArrayRead(hovec,&array);
1090: VecGetLocalSize(hovec,&n);
1091: if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
1092: for (i=0;i<n/sdim;i++) {
1093: for (s=0;s<sdim;s++) {
1094: PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));
1095: }
1096: PetscViewerASCIIPrintf(viewer,"\n");
1097: }
1098: VecRestoreArrayRead(hovec,&array);
1099: } else {
1100: VecGetLocalSize(coordinates,&nvert);
1101: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
1102: PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);
1103: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
1104: VecGetArrayRead(coordinates,&array);
1105: for (p=0;p<nvert/sdim;p++) {
1106: PetscInt s;
1107: for (s=0;s<sdim;s++) {
1108: PetscReal v = PetscRealPart(array[p*sdim+s]);
1110: PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);
1111: }
1112: PetscViewerASCIIPrintf(viewer,"\n");
1113: }
1114: VecRestoreArrayRead(coordinates,&array);
1115: }
1116: VecDestroy(&hovec);
1117: return(0);
1118: }
1120: PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
1121: {
1124: DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);
1125: return(0);
1126: }