xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
18135c375SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
28135c375SStefano Zampini #include <petsc/private/petscimpl.h>
38135c375SStefano Zampini #include <petsc/private/dmpleximpl.h>
48135c375SStefano Zampini #include <petscbt.h>
58135c375SStefano Zampini #include <petscdmplex.h>
68135c375SStefano Zampini #include <petscsf.h>
78135c375SStefano Zampini #include <petscds.h>
88135c375SStefano Zampini 
98135c375SStefano Zampini typedef struct {
108135c375SStefano Zampini   PetscInt   nf;
118135c375SStefano Zampini   VecScatter *scctx;
128135c375SStefano Zampini } GLVisViewerCtx;
138135c375SStefano Zampini 
148135c375SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
158135c375SStefano Zampini {
168135c375SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
178135c375SStefano Zampini   PetscInt       i;
188135c375SStefano Zampini 
198135c375SStefano Zampini   PetscFunctionBegin;
208135c375SStefano Zampini   for (i=0;i<ctx->nf;i++) {
21*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterDestroy(&ctx->scctx[i]));
228135c375SStefano Zampini   }
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ctx->scctx));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(vctx));
258135c375SStefano Zampini   PetscFunctionReturn(0);
268135c375SStefano Zampini }
278135c375SStefano Zampini 
288135c375SStefano Zampini static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
298135c375SStefano Zampini {
308135c375SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
318135c375SStefano Zampini   PetscInt       f;
328135c375SStefano Zampini 
338135c375SStefano Zampini   PetscFunctionBegin;
348135c375SStefano Zampini   for (f=0;f<nf;f++) {
35*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD));
36*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD));
378135c375SStefano Zampini   }
388135c375SStefano Zampini   PetscFunctionReturn(0);
398135c375SStefano Zampini }
408135c375SStefano Zampini 
418135c375SStefano Zampini /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */
428135c375SStefano Zampini PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer)
438135c375SStefano Zampini {
448135c375SStefano Zampini   DM             dm = (DM)odm;
454cac2994SStefano Zampini   Vec            xlocal,xfield,*Ufield;
468135c375SStefano Zampini   PetscDS        ds;
478135c375SStefano Zampini   IS             globalNum,isfield;
488135c375SStefano Zampini   PetscBT        vown;
498135c375SStefano Zampini   char           **fieldname = NULL,**fec_type = NULL;
508135c375SStefano Zampini   const PetscInt *gNum;
51bb77a09fSStefano Zampini   PetscInt       *nlocal,*bs,*idxs,*dims;
528135c375SStefano Zampini   PetscInt       f,maxfields,nfields,c,totc,totdofs,Nv,cum,i;
53b135d7daSStefano Zampini   PetscInt       dim,cStart,cEnd,vStart,vEnd;
548135c375SStefano Zampini   GLVisViewerCtx *ctx;
558135c375SStefano Zampini   PetscSection   s;
568135c375SStefano Zampini 
578135c375SStefano Zampini   PetscFunctionBegin;
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum));
62b135d7daSStefano Zampini   if (!globalNum) {
63*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum));
64*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum));
65*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectDereference((PetscObject)globalNum));
66b135d7daSStefano Zampini   }
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(globalNum,&gNum));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTCreate(vEnd-vStart,&vown));
698135c375SStefano Zampini   for (c = cStart, totc = 0; c < cEnd; c++) {
708135c375SStefano Zampini     if (gNum[c-cStart] >= 0) {
718135c375SStefano Zampini       PetscInt i,numPoints,*points = NULL;
728135c375SStefano Zampini 
738135c375SStefano Zampini       totc++;
74*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points));
758135c375SStefano Zampini       for (i=0;i<numPoints*2;i+= 2) {
768135c375SStefano Zampini         if ((points[i] >= vStart) && (points[i] < vEnd)) {
77*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscBTSet(vown,points[i]-vStart));
788135c375SStefano Zampini         }
798135c375SStefano Zampini       }
80*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points));
818135c375SStefano Zampini     }
828135c375SStefano Zampini   }
8377eacf09SStefano Zampini   for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++;
848135c375SStefano Zampini 
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(dm,&xlocal));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(xlocal,&totdofs));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm,&s));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetNumFields(s,&nfields));
898135c375SStefano Zampini   for (f=0,maxfields=0;f<nfields;f++) {
908135c375SStefano Zampini     PetscInt bs;
918135c375SStefano Zampini 
92*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetFieldComponents(s,f,&bs));
938135c375SStefano Zampini     maxfields += bs;
948135c375SStefano Zampini   }
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield));
96*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&ctx));
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(maxfields,&ctx->scctx));
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm,&ds));
998135c375SStefano Zampini   if (ds) {
1008135c375SStefano Zampini     for (f=0;f<nfields;f++) {
1018135c375SStefano Zampini       const char* fname;
1028135c375SStefano Zampini       char        name[256];
1038135c375SStefano Zampini       PetscObject disc;
1048135c375SStefano Zampini       size_t      len;
1058135c375SStefano Zampini 
106*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetFieldName(s,f,&fname));
107*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscStrlen(fname,&len));
1088135c375SStefano Zampini       if (len) {
109*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscStrcpy(name,fname));
1108135c375SStefano Zampini       } else {
111*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSNPrintf(name,256,"Field%D",f));
1128135c375SStefano Zampini       }
113*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSGetDiscretization(ds,f,&disc));
1148135c375SStefano Zampini       if (disc) {
1158135c375SStefano Zampini         PetscClassId id;
1168135c375SStefano Zampini         PetscInt     Nc;
1178135c375SStefano Zampini         char         fec[64];
1188135c375SStefano Zampini 
119*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscObjectGetClassId(disc, &id));
1208135c375SStefano Zampini         if (id == PETSCFE_CLASSID) {
1218135c375SStefano Zampini           PetscFE            fem = (PetscFE)disc;
1228135c375SStefano Zampini           PetscDualSpace     sp;
1238135c375SStefano Zampini           PetscDualSpaceType spname;
1248135c375SStefano Zampini           PetscInt           order;
1258135c375SStefano Zampini           PetscBool          islag,continuous,H1 = PETSC_TRUE;
1268135c375SStefano Zampini 
127*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscFEGetNumComponents(fem,&Nc));
128*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscFEGetDualSpace(fem,&sp));
129*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscDualSpaceGetType(sp,&spname));
130*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag));
1312c71b3e2SJacob Faibussowitsch           PetscCheckFalse(!islag,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
132*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscDualSpaceLagrangeGetContinuity(sp,&continuous));
133*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscDualSpaceGetOrder(sp,&order));
13428d58a37SPierre Jolivet           if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */
135*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim));
1368135c375SStefano Zampini           } else {
1372c71b3e2SJacob Faibussowitsch             PetscCheckFalse(!continuous && order,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %D",order);
1388135c375SStefano Zampini             H1   = PETSC_FALSE;
139*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order));
1408135c375SStefano Zampini           }
141*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscStrallocpy(name,&fieldname[ctx->nf]));
1428135c375SStefano Zampini           bs[ctx->nf]   = Nc;
143bb77a09fSStefano Zampini           dims[ctx->nf] = dim;
1448135c375SStefano Zampini           if (H1) {
1458135c375SStefano Zampini             nlocal[ctx->nf] = Nc * Nv;
146*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscStrallocpy(fec,&fec_type[ctx->nf]));
147*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield));
1488135c375SStefano Zampini             for (i=0,cum=0;i<vEnd-vStart;i++) {
1498135c375SStefano Zampini               PetscInt j,off;
1508135c375SStefano Zampini 
1518135c375SStefano Zampini               if (PetscUnlikely(!PetscBTLookup(vown,i))) continue;
152*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscSectionGetFieldOffset(s,i+vStart,f,&off));
1538135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1548135c375SStefano Zampini             }
155*5f80ce2aSJacob Faibussowitsch             CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield));
1568135c375SStefano Zampini           } else {
1578135c375SStefano Zampini             nlocal[ctx->nf] = Nc * totc;
158*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscStrallocpy(fec,&fec_type[ctx->nf]));
159*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield));
1608135c375SStefano Zampini             for (i=0,cum=0;i<cEnd-cStart;i++) {
1618135c375SStefano Zampini               PetscInt j,off;
1628135c375SStefano Zampini 
1638135c375SStefano Zampini               if (PetscUnlikely(gNum[i] < 0)) continue;
164*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscSectionGetFieldOffset(s,i+cStart,f,&off));
1658135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1668135c375SStefano Zampini             }
167*5f80ce2aSJacob Faibussowitsch             CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield));
1688135c375SStefano Zampini           }
169*5f80ce2aSJacob Faibussowitsch           CHKERRQ(VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]));
170*5f80ce2aSJacob Faibussowitsch           CHKERRQ(VecDestroy(&xfield));
171*5f80ce2aSJacob Faibussowitsch           CHKERRQ(ISDestroy(&isfield));
1728135c375SStefano Zampini           ctx->nf++;
1738135c375SStefano Zampini         } else if (id == PETSCFV_CLASSID) {
1748135c375SStefano Zampini           PetscInt c;
1758135c375SStefano Zampini 
176*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscFVGetNumComponents((PetscFV)disc,&Nc));
177*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim));
1788135c375SStefano Zampini           for (c = 0; c < Nc; c++) {
1798135c375SStefano Zampini             char comp[256];
180*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscSNPrintf(comp,256,"%s-Comp%D",name,c));
181*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscStrallocpy(comp,&fieldname[ctx->nf]));
1828135c375SStefano Zampini             bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */
1838135c375SStefano Zampini             nlocal[ctx->nf] = totc;
184bb77a09fSStefano Zampini             dims[ctx->nf] = dim;
185*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscStrallocpy(fec,&fec_type[ctx->nf]));
186*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,totc,&xfield));
1878135c375SStefano Zampini             for (i=0,cum=0;i<cEnd-cStart;i++) {
1888135c375SStefano Zampini               PetscInt off;
1898135c375SStefano Zampini 
1908135c375SStefano Zampini               if (PetscUnlikely(gNum[i])<0) continue;
191*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscSectionGetFieldOffset(s,i+cStart,f,&off));
1928135c375SStefano Zampini               idxs[cum++] = off + c;
1938135c375SStefano Zampini             }
194*5f80ce2aSJacob Faibussowitsch             CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield));
195*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]));
196*5f80ce2aSJacob Faibussowitsch             CHKERRQ(VecDestroy(&xfield));
197*5f80ce2aSJacob Faibussowitsch             CHKERRQ(ISDestroy(&isfield));
1988135c375SStefano Zampini             ctx->nf++;
1998135c375SStefano Zampini           }
20098921bdaSJacob Faibussowitsch         } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f);
20198921bdaSJacob Faibussowitsch       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f);
2028135c375SStefano Zampini     }
2038135c375SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&vown));
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xlocal));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(globalNum,&gNum));
2078135c375SStefano Zampini 
2084cac2994SStefano Zampini   /* create work vectors */
2094cac2994SStefano Zampini   for (f=0;f<ctx->nf;f++) {
210*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]));
211*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]));
212*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetBlockSize(Ufield[f],bs[f]));
213*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetDM(Ufield[f],dm));
2144cac2994SStefano Zampini   }
2154cac2994SStefano Zampini 
2168135c375SStefano Zampini   /* customize the viewer */
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private));
2188135c375SStefano Zampini   for (f=0;f<ctx->nf;f++) {
219*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(fieldname[f]));
220*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(fec_type[f]));
221*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&Ufield[f]));
2228135c375SStefano Zampini   }
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield));
2248135c375SStefano Zampini   PetscFunctionReturn(0);
2258135c375SStefano Zampini }
2268135c375SStefano Zampini 
227b135d7daSStefano Zampini typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid;
2288135c375SStefano Zampini 
2298135c375SStefano Zampini MFEM_cid mfem_table_cid[4][7]       = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_UNDEF},
2308135c375SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_UNDEF},
2318135c375SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF,MFEM_UNDEF},
232b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_PRISM,MFEM_CUBE } };
2338135c375SStefano Zampini 
234b135d7daSStefano Zampini 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},
235b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF},
236b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF},
237b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_CUBE } };
238044a5661SStefano Zampini 
239f86f7544SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid)
2408135c375SStefano Zampini {
2418135c375SStefano Zampini   DMLabel        dlabel;
242044a5661SStefano Zampini   PetscInt       depth,csize,pdepth,dim;
2438135c375SStefano Zampini 
2448135c375SStefano Zampini   PetscFunctionBegin;
245*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthLabel(dm,&dlabel));
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetValue(dlabel,p,&pdepth));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetConeSize(dm,p,&csize));
248*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepth(dm,&depth));
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
2508135c375SStefano Zampini   if (label) {
251*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetValue(label,p,mid));
252f86f7544SStefano Zampini     *mid = *mid - minl + 1; /* MFEM does not like negative markers */
25377eacf09SStefano Zampini   } else *mid = 1;
254044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
2552c71b3e2SJacob Faibussowitsch     PetscCheckFalse(dim < 0 || dim > 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim);
2562c71b3e2SJacob Faibussowitsch     PetscCheckFalse(csize > 8,PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p);
2572c71b3e2SJacob Faibussowitsch     PetscCheckFalse(depth != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p);
258044a5661SStefano Zampini     *cid = mfem_table_cid_unint[dim][csize];
259044a5661SStefano Zampini   } else {
2602c71b3e2SJacob Faibussowitsch     PetscCheckFalse(csize > 6,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p);
2612c71b3e2SJacob Faibussowitsch     PetscCheckFalse(pdepth < 0 || pdepth > 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p);
262044a5661SStefano Zampini     *cid = mfem_table_cid[pdepth][csize];
263044a5661SStefano Zampini   }
2648135c375SStefano Zampini   PetscFunctionReturn(0);
2658135c375SStefano Zampini }
2668135c375SStefano Zampini 
26796ca5757SLisandro Dalcin static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[])
2688135c375SStefano Zampini {
269cc0d3ed7SStefano Zampini   PetscInt       dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;
2708135c375SStefano Zampini 
2718135c375SStefano Zampini   PetscFunctionBegin;
272*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd));
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
274cc0d3ed7SStefano Zampini   sdim = dim;
275cc0d3ed7SStefano Zampini   if (csec) {
27684f354e3SLisandro Dalcin     PetscInt sStart,sEnd;
27784f354e3SLisandro Dalcin 
278*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinateDim(dm,&sdim));
279*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetChart(csec,&sStart,&sEnd));
280*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(csec,vStart,&off));
281cc0d3ed7SStefano Zampini     off  = off/sdim;
28284f354e3SLisandro Dalcin     if (p >= sStart && p < sEnd) {
283*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(csec,p,&dof));
284cc0d3ed7SStefano Zampini     }
28584f354e3SLisandro Dalcin   }
286cc0d3ed7SStefano Zampini   if (!dof) {
287*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points));
2888135c375SStefano Zampini     for (i=0,q=0;i<numPoints*2;i+= 2)
2898135c375SStefano Zampini       if ((points[i] >= vStart) && (points[i] < vEnd))
29096ca5757SLisandro Dalcin         vids[q++] = points[i]-vStart+off;
291*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points));
292cc0d3ed7SStefano Zampini   } else {
293*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(csec,p,&off));
294*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(csec,p,&dof));
29596ca5757SLisandro Dalcin     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
296cc0d3ed7SStefano Zampini   }
2978135c375SStefano Zampini   *nv = q;
2988135c375SStefano Zampini   PetscFunctionReturn(0);
2998135c375SStefano Zampini }
3008135c375SStefano Zampini 
301066ea43fSLisandro Dalcin static PetscErrorCode GLVisCreateFE(PetscFE femIn,char name[32],PetscFE *fem)
302066ea43fSLisandro Dalcin {
303066ea43fSLisandro Dalcin   DM              K;
304066ea43fSLisandro Dalcin   PetscSpace      P;
305066ea43fSLisandro Dalcin   PetscDualSpace  Q;
306066ea43fSLisandro Dalcin   PetscQuadrature q,fq;
307066ea43fSLisandro Dalcin   PetscInt        dim,deg,dof;
308066ea43fSLisandro Dalcin   DMPolytopeType  ptype;
309066ea43fSLisandro Dalcin   PetscBool       isSimplex,isTensor;
310066ea43fSLisandro Dalcin   PetscBool       continuity = PETSC_FALSE;
311066ea43fSLisandro Dalcin   PetscDTNodeType nodeType   = PETSCDTNODES_GAUSSJACOBI;
312066ea43fSLisandro Dalcin   PetscBool       endpoint   = PETSC_TRUE;
313066ea43fSLisandro Dalcin   MPI_Comm        comm;
314066ea43fSLisandro Dalcin 
315066ea43fSLisandro Dalcin   PetscFunctionBegin;
316*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)femIn, &comm));
317*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetBasisSpace(femIn,&P));
318*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetDualSpace(femIn,&Q));
319*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceGetDM(Q,&K));
320*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(K,&dim));
321*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceGetDegree(P,&deg,NULL));
322*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceGetNumComponents(P,&dof));
323*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(K,0,&ptype));
324066ea43fSLisandro Dalcin   switch (ptype) {
325066ea43fSLisandro Dalcin   case DM_POLYTOPE_QUADRILATERAL:
326066ea43fSLisandro Dalcin   case DM_POLYTOPE_HEXAHEDRON:
327066ea43fSLisandro Dalcin     isSimplex = PETSC_FALSE; break;
328066ea43fSLisandro Dalcin   default:
329066ea43fSLisandro Dalcin     isSimplex = PETSC_TRUE; break;
330066ea43fSLisandro Dalcin   }
331066ea43fSLisandro Dalcin   isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
332066ea43fSLisandro Dalcin   /* Create space */
333*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceCreate(comm,&P));
334*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL));
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpacePolynomialSetTensor(P,isTensor));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetNumComponents(P,dof));
337*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetNumVariables(P,dim));
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetDegree(P,deg,PETSC_DETERMINE));
339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceSetUp(P));
340066ea43fSLisandro Dalcin   /* Create dual space */
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceCreate(comm,&Q));
342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE));
343*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceLagrangeSetTensor(Q,isTensor));
344*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceLagrangeSetContinuity(Q,continuity));
345*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceLagrangeSetNodeType(Q,nodeType,endpoint,0));
346*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetNumComponents(Q,dof));
347*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetOrder(Q,deg));
348*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
349*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetDM(Q,K));
350*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&K));
351*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetUp(Q));
352066ea43fSLisandro Dalcin   /* Create quadrature */
353066ea43fSLisandro Dalcin   if (isSimplex) {
354*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTStroudConicalQuadrature(dim,  1,deg+1,-1,+1,&q));
355*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTStroudConicalQuadrature(dim-1,1,deg+1,-1,+1,&fq));
356066ea43fSLisandro Dalcin   } else {
357*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTGaussTensorQuadrature(dim,  1,deg+1,-1,+1,&q));
358*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTGaussTensorQuadrature(dim-1,1,deg+1,-1,+1,&fq));
359066ea43fSLisandro Dalcin   }
360066ea43fSLisandro Dalcin   /* Create finite element */
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreate(comm,fem));
362*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(name,32,"L2_T1_%DD_P%D",dim,deg));
363*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)*fem,name));
364*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetType(*fem,PETSCFEBASIC));
365*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetNumComponents(*fem,dof));
366*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetBasisSpace(*fem,P));
367*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetDualSpace(*fem,Q));
368*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetQuadrature(*fem,q));
369*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetFaceQuadrature(*fem,fq));
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFESetUp(*fem));
371066ea43fSLisandro Dalcin   /* Cleanup */
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceDestroy(&P));
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceDestroy(&Q));
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&q));
375*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&fq));
376066ea43fSLisandro Dalcin   PetscFunctionReturn(0);
377066ea43fSLisandro Dalcin }
378066ea43fSLisandro Dalcin 
37977eacf09SStefano Zampini /*
38077eacf09SStefano Zampini    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
38177eacf09SStefano Zampini    Higher order meshes are also supported
38277eacf09SStefano Zampini */
3838135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
3848135c375SStefano Zampini {
3858135c375SStefano Zampini   DMLabel              label;
3868135c375SStefano Zampini   PetscSection         coordSection,parentSection;
38777eacf09SStefano Zampini   Vec                  coordinates,hovec;
3888135c375SStefano Zampini   const PetscScalar    *array;
389f86f7544SStefano Zampini   PetscInt             bf,p,sdim,dim,depth,novl,minl;
390412e9a14SMatthew G. Knepley   PetscInt             cStart,cEnd,vStart,vEnd,nvert;
3913924b612SStefano Zampini   PetscMPIInt          size;
3923e6c54aaSStefano Zampini   PetscBool            localized,isascii;
39328d58a37SPierre Jolivet   PetscBool            enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE;
3943e6c54aaSStefano Zampini   PetscBT              pown,vown;
3958135c375SStefano Zampini   PetscErrorCode       ierr;
3968135c375SStefano Zampini   PetscContainer       glvis_container;
397044a5661SStefano Zampini   PetscBool            cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
398f86f7544SStefano Zampini   PetscBool            enable_emark,enable_bmark;
39977eacf09SStefano Zampini   const char           *fmt;
4007bf4dd16SStefano Zampini   char                 emark[64] = "",bmark[64] = "";
4018135c375SStefano Zampini 
4028135c375SStefano Zampini   PetscFunctionBegin;
4038135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4048135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
405*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
4062c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!isascii,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
407*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size));
4082c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm,&dim));
4108135c375SStefano Zampini 
4118135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
412*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container));
4132c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!glvis_container,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
4148135c375SStefano Zampini   {
4158135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
416*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscContainerGetPointer(glvis_container,(void**)&glvis_info));
4178135c375SStefano Zampini     enabled = glvis_info->enabled;
41877eacf09SStefano Zampini     fmt     = glvis_info->fmt;
4198135c375SStefano Zampini   }
42021414b21SStefano Zampini 
42121414b21SStefano Zampini   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
42221414b21SStefano Zampini      DMPlex does not currently support HO meshes, so there's no API for this */
423*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec));
424*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)hovec));
425066ea43fSLisandro Dalcin   if (!hovec) {
426066ea43fSLisandro Dalcin     DM           cdm;
427066ea43fSLisandro Dalcin     PetscFE      disc;
428066ea43fSLisandro Dalcin     PetscClassId classid;
429066ea43fSLisandro Dalcin 
430*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinateDM(dm,&cdm));
431*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetField(cdm,0,NULL,(PetscObject*)&disc));
432*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId((PetscObject)disc,&classid));
433066ea43fSLisandro Dalcin     if (classid == PETSCFE_CLASSID) {
434066ea43fSLisandro Dalcin       DM      hocdm;
435066ea43fSLisandro Dalcin       PetscFE hodisc;
436066ea43fSLisandro Dalcin       Vec     vec;
437066ea43fSLisandro Dalcin       Mat     mat;
438066ea43fSLisandro Dalcin       char    name[32],fec_type[64];
439066ea43fSLisandro Dalcin 
440*5f80ce2aSJacob Faibussowitsch       CHKERRQ(GLVisCreateFE(disc,name,&hodisc));
441*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMClone(cdm,&hocdm));
442*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetField(hocdm,0,NULL,(PetscObject)hodisc));
443*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFEDestroy(&hodisc));
444*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateDS(hocdm));
445066ea43fSLisandro Dalcin 
446*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetCoordinates(dm,&vec));
447*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateGlobalVector(hocdm,&hovec));
448*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateInterpolation(cdm,hocdm,&mat,NULL));
449*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatInterpolate(mat,vec,hovec));
450*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&mat));
451*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&hocdm));
452066ea43fSLisandro Dalcin 
453*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name));
454*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject)hovec,fec_type));
455066ea43fSLisandro Dalcin     }
456066ea43fSLisandro Dalcin   }
45721414b21SStefano Zampini 
458*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetGhostCellStratum(dm,&p,NULL));
460c3c203b2SStefano Zampini   if (p >= 0) cEnd = p;
461*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd));
462*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL));
463*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocalized(dm,&localized));
4642c71b3e2SJacob Faibussowitsch   PetscCheckFalse(periodic && !localized && !hovec,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
465*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dm,&coordSection));
466*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dm,&sdim));
467*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm,&coordinates));
4682c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!coordinates && !hovec,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");
4698135c375SStefano Zampini 
4708135c375SStefano Zampini   /*
4718135c375SStefano Zampini      a couple of sections of the mesh specification are disabled
4723e96840aSStefano Zampini        - 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)
47377eacf09SStefano Zampini        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
4743e6c54aaSStefano Zampini                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
4758135c375SStefano Zampini   */
4763e96840aSStefano Zampini   enable_boundary = PETSC_FALSE;
4778135c375SStefano Zampini   enable_ncmesh   = PETSC_FALSE;
4783e6c54aaSStefano Zampini   enable_mfem     = PETSC_FALSE;
479f86f7544SStefano Zampini   enable_emark    = PETSC_FALSE;
480f86f7544SStefano Zampini   enable_bmark    = PETSC_FALSE;
4817bf4dd16SStefano Zampini   /* I'm tired of problems with negative values in the markers, disable them */
4828135c375SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
483*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL));
484*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL));
485*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL));
486*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL));
487*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark));
488*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark));
4898135c375SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
490f86f7544SStefano Zampini   if (enable_bmark) enable_boundary = PETSC_TRUE;
491f86f7544SStefano Zampini 
492*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size));
4932c71b3e2SJacob Faibussowitsch   PetscCheckFalse(enable_ncmesh && size > 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
494*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepth(dm,&depth));
4952c71b3e2SJacob Faibussowitsch   PetscCheckFalse(enable_boundary && depth >= 0 && dim != depth,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
4967e1aca4eSStefano Zampini                                                              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
4972c71b3e2SJacob Faibussowitsch   PetscCheckFalse(enable_ncmesh && depth >= 0 && dim != depth,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
4987e1aca4eSStefano Zampini                                                            "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
499044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
5002c71b3e2SJacob Faibussowitsch     PetscCheckFalse(depth != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth);
501044a5661SStefano Zampini     cellvertex = PETSC_TRUE;
502044a5661SStefano Zampini   }
5038135c375SStefano Zampini 
5048135c375SStefano Zampini   /* Identify possible cells in the overlap */
5058135c375SStefano Zampini   novl = 0;
5068135c375SStefano Zampini   pown = NULL;
5073924b612SStefano Zampini   if (size > 1) {
5083e6c54aaSStefano Zampini     IS             globalNum = NULL;
5093e6c54aaSStefano Zampini     const PetscInt *gNum;
5103e6c54aaSStefano Zampini     PetscBool      ovl  = PETSC_FALSE;
5113e6c54aaSStefano Zampini 
512*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum));
513b135d7daSStefano Zampini     if (!globalNum) {
51428d58a37SPierre Jolivet       if (view_ovl) {
515*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum));
51628d58a37SPierre Jolivet       } else {
517*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum));
51828d58a37SPierre Jolivet       }
519*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum));
520*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectDereference((PetscObject)globalNum));
521b135d7daSStefano Zampini     }
522*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(globalNum,&gNum));
5238135c375SStefano Zampini     for (p=cStart; p<cEnd; p++) {
5248135c375SStefano Zampini       if (gNum[p-cStart] < 0) {
5258135c375SStefano Zampini         ovl = PETSC_TRUE;
5268135c375SStefano Zampini         novl++;
5278135c375SStefano Zampini       }
5288135c375SStefano Zampini     }
5298135c375SStefano Zampini     if (ovl) {
5308135c375SStefano Zampini       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
5318135c375SStefano Zampini          TODO: garbage collector? attach pown to dm?  */
532*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBTCreate(cEnd-cStart,&pown));
5333e6c54aaSStefano Zampini       for (p=cStart; p<cEnd; p++) {
5343e6c54aaSStefano Zampini         if (gNum[p-cStart] < 0) continue;
5353e6c54aaSStefano Zampini         else {
536*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscBTSet(pown,p-cStart));
5378135c375SStefano Zampini         }
5388135c375SStefano Zampini       }
5393e6c54aaSStefano Zampini     }
540*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(globalNum,&gNum));
5413e6c54aaSStefano Zampini   }
5428135c375SStefano Zampini 
5431d4815f0SStefano Zampini   /* vertex_parents (Non-conforming meshes) */
5441d4815f0SStefano Zampini   parentSection  = NULL;
5451d4815f0SStefano Zampini   if (enable_ncmesh) {
546*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL));
5471d4815f0SStefano Zampini     enable_ncmesh = (PetscBool)(enable_ncmesh && parentSection);
5481d4815f0SStefano Zampini   }
5493e6c54aaSStefano Zampini   /* return if this process is disabled */
5508135c375SStefano Zampini   if (!enabled) {
551*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"MFEM mesh %s\n",enable_ncmesh ? "v1.1" : "v1.0"));
552*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\ndimension\n"));
553*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",dim));
554*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nelements\n"));
555*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",0));
556*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nboundary\n"));
557*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",0));
558*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
559*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",0));
560*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",sdim));
561*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTDestroy(&pown));
562*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&hovec));
5638135c375SStefano Zampini     PetscFunctionReturn(0);
5648135c375SStefano Zampini   }
5658135c375SStefano Zampini 
5663e6c54aaSStefano Zampini   if (enable_mfem) {
5673e6c54aaSStefano Zampini     if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
5683e6c54aaSStefano Zampini       PetscInt    vpc = 0;
5693e6c54aaSStefano Zampini       char        fec[64];
57096ca5757SLisandro Dalcin       PetscInt    vids[8] = {0,1,2,3,4,5,6,7};
57196ca5757SLisandro Dalcin       PetscInt    hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3};
57296ca5757SLisandro Dalcin       PetscInt    quadv[8] = {0,1,3,2}, triv[3] = {0,1,2};
57396ca5757SLisandro Dalcin       PetscInt    *dof = NULL;
5743e6c54aaSStefano Zampini       PetscScalar *array,*ptr;
5753e6c54aaSStefano Zampini 
576*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim));
5773e6c54aaSStefano Zampini       if (cEnd-cStart) {
5783e6c54aaSStefano Zampini         PetscInt fpc;
5793e6c54aaSStefano Zampini 
580*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetConeSize(dm,cStart,&fpc));
5813e6c54aaSStefano Zampini         switch(dim) {
5823e6c54aaSStefano Zampini           case 1:
5833e6c54aaSStefano Zampini             vpc = 2;
5843e6c54aaSStefano Zampini             dof = hexv;
5853e6c54aaSStefano Zampini             break;
5863e6c54aaSStefano Zampini           case 2:
5873e6c54aaSStefano Zampini             switch (fpc) {
5883e6c54aaSStefano Zampini               case 3:
5893e6c54aaSStefano Zampini                 vpc = 3;
590044a5661SStefano Zampini                 dof = triv;
5913e6c54aaSStefano Zampini                 break;
5923e6c54aaSStefano Zampini               case 4:
5933e6c54aaSStefano Zampini                 vpc = 4;
594044a5661SStefano Zampini                 dof = quadv;
5953e6c54aaSStefano Zampini                 break;
5963e6c54aaSStefano Zampini               default:
59798921bdaSJacob Faibussowitsch                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
5983e6c54aaSStefano Zampini             }
5993e6c54aaSStefano Zampini             break;
6003e6c54aaSStefano Zampini           case 3:
6013e6c54aaSStefano Zampini             switch (fpc) {
602044a5661SStefano Zampini               case 4: /* TODO: still need to understand L2 ordering for tets */
6033e6c54aaSStefano Zampini                 vpc = 4;
604044a5661SStefano Zampini                 dof = tetv;
605044a5661SStefano Zampini                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
6063e6c54aaSStefano Zampini               case 6:
6072c71b3e2SJacob Faibussowitsch                 PetscCheckFalse(cellvertex,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc);
608044a5661SStefano Zampini                 vpc = 8;
609044a5661SStefano Zampini                 dof = hexv;
610044a5661SStefano Zampini                 break;
611044a5661SStefano Zampini               case 8:
6122c71b3e2SJacob Faibussowitsch                 PetscCheckFalse(!cellvertex,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
6133e6c54aaSStefano Zampini                 vpc = 8;
6143e6c54aaSStefano Zampini                 dof = hexv;
6153e6c54aaSStefano Zampini                 break;
6163e6c54aaSStefano Zampini               default:
61798921bdaSJacob Faibussowitsch                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
6183e6c54aaSStefano Zampini             }
6193e6c54aaSStefano Zampini             break;
6203e6c54aaSStefano Zampini           default:
6213e6c54aaSStefano Zampini             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
6223e6c54aaSStefano Zampini         }
623*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexReorderCell(dm,cStart,vids));
6243e6c54aaSStefano Zampini       }
6252c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!dof,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
626*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec));
627*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject)hovec,fec));
628*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArray(hovec,&array));
6293e6c54aaSStefano Zampini       ptr  = array;
6303e6c54aaSStefano Zampini       for (p=cStart;p<cEnd;p++) {
6313e6c54aaSStefano Zampini         PetscInt    csize,v,d;
6323e6c54aaSStefano Zampini         PetscScalar *vals = NULL;
6333e6c54aaSStefano Zampini 
6343e6c54aaSStefano Zampini         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
635*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals));
6362c71b3e2SJacob Faibussowitsch         PetscCheckFalse(csize != vpc*sdim && csize != vpc*sdim*2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported closure size %D (vpc %D, sdim %D)",csize,vpc,sdim);
6373e6c54aaSStefano Zampini         for (v=0;v<vpc;v++) {
6383e6c54aaSStefano Zampini           for (d=0;d<sdim;d++) {
6393e6c54aaSStefano Zampini             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
6403e6c54aaSStefano Zampini           }
6413e6c54aaSStefano Zampini         }
6423e6c54aaSStefano Zampini         ptr += vpc*sdim;
643*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals));
6443e6c54aaSStefano Zampini       }
645*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(hovec,&array));
6463e6c54aaSStefano Zampini     }
6473e6c54aaSStefano Zampini   }
6483e96840aSStefano Zampini   /* if we have high-order coordinates in 3D, we need to specify the boundary */
6493e96840aSStefano Zampini   if (hovec && dim == 3) enable_boundary = PETSC_TRUE;
6503e6c54aaSStefano Zampini 
6518135c375SStefano Zampini   /* header */
652*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"MFEM mesh %s\n",enable_ncmesh ? "v1.1" : "v1.0"));
6538135c375SStefano Zampini 
6548135c375SStefano Zampini   /* topological dimension */
655*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"\ndimension\n"));
656*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",dim));
6578135c375SStefano Zampini 
6588135c375SStefano Zampini   /* elements */
659f86f7544SStefano Zampini   minl = 1;
660f86f7544SStefano Zampini   label = NULL;
661f86f7544SStefano Zampini   if (enable_emark) {
662f86f7544SStefano Zampini     PetscInt lminl = PETSC_MAX_INT;
663f86f7544SStefano Zampini 
664*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dm,emark,&label));
665f86f7544SStefano Zampini     if (label) {
666f86f7544SStefano Zampini       IS       vals;
667f86f7544SStefano Zampini       PetscInt ldef;
668f86f7544SStefano Zampini 
669*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMLabelGetDefaultValue(label,&ldef));
670*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMLabelGetValueIS(label,&vals));
671*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISGetMinMax(vals,&lminl,NULL));
672*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&vals));
673f86f7544SStefano Zampini       lminl = PetscMin(ldef,lminl);
674f86f7544SStefano Zampini     }
675*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm)));
676f86f7544SStefano Zampini     if (minl == PETSC_MAX_INT) minl = 1;
677f86f7544SStefano Zampini   }
678*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nelements\n"));
679*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl));
6808135c375SStefano Zampini   for (p=cStart;p<cEnd;p++) {
68196ca5757SLisandro Dalcin     PetscInt       vids[8];
68211a4995dSStefano Zampini     PetscInt       i,nv = 0,cid = -1,mid = 1;
6838135c375SStefano Zampini 
6843e6c54aaSStefano Zampini     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
685*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid));
686*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids));
687*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexReorderCell(dm,p,vids));
688*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid));
6898135c375SStefano Zampini     for (i=0;i<nv;i++) {
690*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer," %D",vids[i]));
6918135c375SStefano Zampini     }
692*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n"));
6938135c375SStefano Zampini   }
6948135c375SStefano Zampini 
695cc0d3ed7SStefano Zampini   /* boundary */
696*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nboundary\n"));
697cc0d3ed7SStefano Zampini   if (!enable_boundary) {
698*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",0));
699cc0d3ed7SStefano Zampini   } else {
70077eacf09SStefano Zampini     DMLabel  perLabel;
70177eacf09SStefano Zampini     PetscBT  bfaces;
702b135d7daSStefano Zampini     PetscInt fStart,fEnd,*fcells;
703cc0d3ed7SStefano Zampini 
704*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(dm,1,&fStart,&fEnd));
705*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTCreate(fEnd-fStart,&bfaces));
706*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetMaxSizes(dm,NULL,&p));
707*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(p,&fcells));
708*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dm,"glvis_periodic_cut",&perLabel));
70977eacf09SStefano Zampini     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
710*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateLabel(dm,"glvis_periodic_cut"));
711*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetLabel(dm,"glvis_periodic_cut",&perLabel));
712*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMLabelSetDefaultValue(perLabel,1));
71377eacf09SStefano Zampini       for (p=cStart;p<cEnd;p++) {
714c3c203b2SStefano Zampini         DMPolytopeType cellType;
715c3c203b2SStefano Zampini         PetscInt       dof;
716b135d7daSStefano Zampini 
717*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetCellType(dm,p,&cellType));
718*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetDof(coordSection,p,&dof));
71977eacf09SStefano Zampini         if (dof) {
720c3c203b2SStefano Zampini           PetscInt    uvpc, v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
72177eacf09SStefano Zampini           PetscScalar *vals = NULL;
722c3c203b2SStefano Zampini 
723c3c203b2SStefano Zampini           uvpc = DMPolytopeTypeGetNumVertices(cellType);
7242c71b3e2SJacob Faibussowitsch           PetscCheckFalse(dof%sdim,PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
725*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals));
726*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure));
72777eacf09SStefano Zampini           for (v=0;v<cellClosureSize;v++)
72877eacf09SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
72977eacf09SStefano Zampini               vidxs = cellClosure + 2*v;
73077eacf09SStefano Zampini               break;
73177eacf09SStefano Zampini             }
7322c71b3e2SJacob Faibussowitsch           PetscCheckFalse(!vidxs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
733b135d7daSStefano Zampini           for (v=0;v<uvpc;v++) {
73477eacf09SStefano Zampini             PetscInt s;
735044a5661SStefano Zampini 
73677eacf09SStefano Zampini             for (s=0;s<sdim;s++) {
737b135d7daSStefano Zampini               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) {
738*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(DMLabelSetValue(perLabel,vidxs[2*v],2));
73977eacf09SStefano Zampini               }
74077eacf09SStefano Zampini             }
74177eacf09SStefano Zampini           }
742*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure));
743*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals));
74477eacf09SStefano Zampini         }
74577eacf09SStefano Zampini       }
74677eacf09SStefano Zampini       if (dim > 1) {
747b135d7daSStefano Zampini         PetscInt eEnd,eStart;
748044a5661SStefano Zampini 
749*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetDepthStratum(dm,1,&eStart,&eEnd));
75077eacf09SStefano Zampini         for (p=eStart;p<eEnd;p++) {
75177eacf09SStefano Zampini           const PetscInt *cone;
75277eacf09SStefano Zampini           PetscInt       coneSize,i;
75377eacf09SStefano Zampini           PetscBool      ispe = PETSC_TRUE;
75477eacf09SStefano Zampini 
755*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetCone(dm,p,&cone));
756*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetConeSize(dm,p,&coneSize));
75777eacf09SStefano Zampini           for (i=0;i<coneSize;i++) {
75877eacf09SStefano Zampini             PetscInt v;
75977eacf09SStefano Zampini 
760*5f80ce2aSJacob Faibussowitsch             CHKERRQ(DMLabelGetValue(perLabel,cone[i],&v));
76177eacf09SStefano Zampini             ispe = (PetscBool)(ispe && (v==2));
76277eacf09SStefano Zampini           }
76377eacf09SStefano Zampini           if (ispe && coneSize) {
7643e96840aSStefano Zampini             PetscInt       ch, numChildren;
7653e96840aSStefano Zampini             const PetscInt *children;
7663e96840aSStefano Zampini 
767*5f80ce2aSJacob Faibussowitsch             CHKERRQ(DMLabelSetValue(perLabel,p,2));
768*5f80ce2aSJacob Faibussowitsch             CHKERRQ(DMPlexGetTreeChildren(dm,p,&numChildren,&children));
7693e96840aSStefano Zampini             for (ch = 0; ch < numChildren; ch++) {
770*5f80ce2aSJacob Faibussowitsch               CHKERRQ(DMLabelSetValue(perLabel,children[ch],2));
7713e96840aSStefano Zampini             }
77277eacf09SStefano Zampini           }
77377eacf09SStefano Zampini         }
77477eacf09SStefano Zampini         if (dim > 2) {
77577eacf09SStefano Zampini           for (p=fStart;p<fEnd;p++) {
77677eacf09SStefano Zampini             const PetscInt *cone;
77777eacf09SStefano Zampini             PetscInt       coneSize,i;
77877eacf09SStefano Zampini             PetscBool      ispe = PETSC_TRUE;
77977eacf09SStefano Zampini 
780*5f80ce2aSJacob Faibussowitsch             CHKERRQ(DMPlexGetCone(dm,p,&cone));
781*5f80ce2aSJacob Faibussowitsch             CHKERRQ(DMPlexGetConeSize(dm,p,&coneSize));
78277eacf09SStefano Zampini             for (i=0;i<coneSize;i++) {
78377eacf09SStefano Zampini               PetscInt v;
78477eacf09SStefano Zampini 
785*5f80ce2aSJacob Faibussowitsch               CHKERRQ(DMLabelGetValue(perLabel,cone[i],&v));
78677eacf09SStefano Zampini               ispe = (PetscBool)(ispe && (v==2));
78777eacf09SStefano Zampini             }
78877eacf09SStefano Zampini             if (ispe && coneSize) {
7893e96840aSStefano Zampini               PetscInt       ch, numChildren;
7903e96840aSStefano Zampini               const PetscInt *children;
7913e96840aSStefano Zampini 
792*5f80ce2aSJacob Faibussowitsch               CHKERRQ(DMLabelSetValue(perLabel,p,2));
793*5f80ce2aSJacob Faibussowitsch               CHKERRQ(DMPlexGetTreeChildren(dm,p,&numChildren,&children));
7943e96840aSStefano Zampini               for (ch = 0; ch < numChildren; ch++) {
795*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(DMLabelSetValue(perLabel,children[ch],2));
7963e96840aSStefano Zampini               }
79777eacf09SStefano Zampini             }
79877eacf09SStefano Zampini           }
79977eacf09SStefano Zampini         }
80077eacf09SStefano Zampini       }
80177eacf09SStefano Zampini     }
80277eacf09SStefano Zampini     for (p=fStart;p<fEnd;p++) {
80377eacf09SStefano Zampini       const PetscInt *support;
8048135c375SStefano Zampini       PetscInt       supportSize;
80577eacf09SStefano Zampini       PetscBool      isbf = PETSC_FALSE;
8068135c375SStefano Zampini 
807*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetSupportSize(dm,p,&supportSize));
8083e6c54aaSStefano Zampini       if (pown) {
8098135c375SStefano Zampini         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
81077eacf09SStefano Zampini         PetscInt  i;
81177eacf09SStefano Zampini 
812*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetSupport(dm,p,&support));
81377eacf09SStefano Zampini         for (i=0;i<supportSize;i++) {
81477eacf09SStefano Zampini           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
81577eacf09SStefano Zampini           else has_ghost = PETSC_TRUE;
81677eacf09SStefano Zampini         }
81777eacf09SStefano Zampini         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
81877eacf09SStefano Zampini       } else {
81977eacf09SStefano Zampini         isbf = (PetscBool)(supportSize == 1);
82077eacf09SStefano Zampini       }
82177eacf09SStefano Zampini       if (!isbf && perLabel) {
82277eacf09SStefano Zampini         const PetscInt *cone;
82377eacf09SStefano Zampini         PetscInt       coneSize,i;
82477eacf09SStefano Zampini 
825*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetCone(dm,p,&cone));
826*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetConeSize(dm,p,&coneSize));
82777eacf09SStefano Zampini         isbf = PETSC_TRUE;
82877eacf09SStefano Zampini         for (i=0;i<coneSize;i++) {
82977eacf09SStefano Zampini           PetscInt v,d;
83077eacf09SStefano Zampini 
831*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMLabelGetValue(perLabel,cone[i],&v));
832*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMLabelGetDefaultValue(perLabel,&d));
83377eacf09SStefano Zampini           isbf = (PetscBool)(isbf && v != d);
83477eacf09SStefano Zampini         }
83577eacf09SStefano Zampini       }
83677eacf09SStefano Zampini       if (isbf) {
837*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscBTSet(bfaces,p-fStart));
83877eacf09SStefano Zampini       }
83977eacf09SStefano Zampini     }
84077eacf09SStefano Zampini     /* count boundary faces */
84177eacf09SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
84277eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
84377eacf09SStefano Zampini         const PetscInt *support;
84477eacf09SStefano Zampini         PetscInt       supportSize,c;
8458135c375SStefano Zampini 
846*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetSupportSize(dm,p,&supportSize));
847*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetSupport(dm,p,&support));
84877eacf09SStefano Zampini         for (c=0;c<supportSize;c++) {
8493e96840aSStefano Zampini           const    PetscInt *cone;
850b135d7daSStefano Zampini           PetscInt cell,cl,coneSize;
8513e96840aSStefano Zampini 
8523e96840aSStefano Zampini           cell = support[c];
8533e96840aSStefano Zampini           if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue;
854*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetCone(dm,cell,&cone));
855*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetConeSize(dm,cell,&coneSize));
856b135d7daSStefano Zampini           for (cl=0;cl<coneSize;cl++) {
8573e96840aSStefano Zampini             if (cone[cl] == p) {
8583e96840aSStefano Zampini               bf += 1;
8593e96840aSStefano Zampini               break;
8608135c375SStefano Zampini             }
86177eacf09SStefano Zampini           }
8623e96840aSStefano Zampini         }
8638135c375SStefano Zampini       }
8648135c375SStefano Zampini     }
865f86f7544SStefano Zampini     minl = 1;
866f86f7544SStefano Zampini     label = NULL;
867f86f7544SStefano Zampini     if (enable_bmark) {
868f86f7544SStefano Zampini       PetscInt lminl = PETSC_MAX_INT;
869f86f7544SStefano Zampini 
870*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetLabel(dm,bmark,&label));
871f86f7544SStefano Zampini       if (label) {
872f86f7544SStefano Zampini         IS       vals;
873f86f7544SStefano Zampini         PetscInt ldef;
874f86f7544SStefano Zampini 
875*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMLabelGetDefaultValue(label,&ldef));
876*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMLabelGetValueIS(label,&vals));
877*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetMinMax(vals,&lminl,NULL));
878*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&vals));
879f86f7544SStefano Zampini         lminl = PetscMin(ldef,lminl);
880f86f7544SStefano Zampini       }
881*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm)));
882f86f7544SStefano Zampini       if (minl == PETSC_MAX_INT) minl = 1;
883f86f7544SStefano Zampini     }
884*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",bf));
8858135c375SStefano Zampini     for (p=fStart;p<fEnd;p++) {
88677eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
8878135c375SStefano Zampini         const PetscInt *support;
88877eacf09SStefano Zampini         PetscInt       supportSize,c,nc = 0;
8898135c375SStefano Zampini 
890*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetSupportSize(dm,p,&supportSize));
891*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetSupport(dm,p,&support));
8923e6c54aaSStefano Zampini         if (pown) {
89377eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
89477eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
89577eacf09SStefano Zampini               fcells[nc++] = support[c];
8968135c375SStefano Zampini             }
89777eacf09SStefano Zampini           }
89877eacf09SStefano Zampini         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
89977eacf09SStefano Zampini         for (c=0;c<nc;c++) {
900c3c203b2SStefano Zampini           const DMPolytopeType *faceTypes;
901c3c203b2SStefano Zampini           DMPolytopeType       cellType;
902c3c203b2SStefano Zampini           const PetscInt       *faceSizes,*cone;
903c3c203b2SStefano Zampini           PetscInt             vids[8],*faces,st,i,coneSize,cell,cl,nv,cid = -1,mid = -1;
9048135c375SStefano Zampini 
90577eacf09SStefano Zampini           cell = fcells[c];
906*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetCone(dm,cell,&cone));
907*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetConeSize(dm,cell,&coneSize));
908b135d7daSStefano Zampini           for (cl=0;cl<coneSize;cl++)
90977eacf09SStefano Zampini             if (cone[cl] == p)
9108135c375SStefano Zampini               break;
911b135d7daSStefano Zampini           if (cl == coneSize) continue;
9128135c375SStefano Zampini 
91377eacf09SStefano Zampini           /* face material id and type */
914*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid));
915*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid));
91677eacf09SStefano Zampini           /* vertex ids */
917*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetCellType(dm,cell,&cellType));
918*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids));
919*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces));
920c3c203b2SStefano Zampini           st = 0;
921c3c203b2SStefano Zampini           for (i=0;i<cl;i++) st += faceSizes[i];
922*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexInvertCell(faceTypes[cl],faces + st));
923c3c203b2SStefano Zampini           for (i=0;i<faceSizes[cl];i++) {
924*5f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscViewerASCIIPrintf(viewer," %d",faces[st+i]));
925b135d7daSStefano Zampini           }
926*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n"));
927*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexRestoreRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces));
9283e96840aSStefano Zampini           bf -= 1;
92977eacf09SStefano Zampini         }
9308135c375SStefano Zampini       }
9318135c375SStefano Zampini     }
9322c71b3e2SJacob Faibussowitsch     PetscCheckFalse(bf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
933*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTDestroy(&bfaces));
934*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(fcells));
9358135c375SStefano Zampini   }
9368135c375SStefano Zampini 
9378135c375SStefano Zampini   /* mark owned vertices */
9383e6c54aaSStefano Zampini   vown = NULL;
9393e6c54aaSStefano Zampini   if (pown) {
940*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBTCreate(vEnd-vStart,&vown));
9418135c375SStefano Zampini     for (p=cStart;p<cEnd;p++) {
9428135c375SStefano Zampini       PetscInt i,closureSize,*closure = NULL;
9438135c375SStefano Zampini 
9443e6c54aaSStefano Zampini       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
945*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure));
9468135c375SStefano Zampini       for (i=0;i<closureSize;i++) {
9478135c375SStefano Zampini         const PetscInt pp = closure[2*i];
9488135c375SStefano Zampini 
9498135c375SStefano Zampini         if (pp >= vStart && pp < vEnd) {
950*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscBTSet(vown,pp-vStart));
9518135c375SStefano Zampini         }
9528135c375SStefano Zampini       }
953*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure));
9548135c375SStefano Zampini     }
9558135c375SStefano Zampini   }
9568135c375SStefano Zampini 
9578135c375SStefano Zampini   if (parentSection) {
9588135c375SStefano Zampini     PetscInt vp,gvp;
9598135c375SStefano Zampini 
9608135c375SStefano Zampini     for (vp=0,p=vStart;p<vEnd;p++) {
9618135c375SStefano Zampini       DMLabel  dlabel;
9628135c375SStefano Zampini       PetscInt parent,depth;
9638135c375SStefano Zampini 
9643e6c54aaSStefano Zampini       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
965*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetDepthLabel(dm,&dlabel));
966*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMLabelGetValue(dlabel,p,&depth));
967*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetTreeParent(dm,p,&parent,NULL));
9688135c375SStefano Zampini       if (parent != p) vp++;
9698135c375SStefano Zampini     }
970*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
9718135c375SStefano Zampini     if (gvp) {
9728135c375SStefano Zampini       PetscInt  maxsupp;
9738135c375SStefano Zampini       PetscBool *skip = NULL;
9748135c375SStefano Zampini 
975*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n"));
976*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",vp));
977*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetMaxSizes(dm,NULL,&maxsupp));
978*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(maxsupp,&skip));
9798135c375SStefano Zampini       for (p=vStart;p<vEnd;p++) {
9808135c375SStefano Zampini         DMLabel  dlabel;
9818135c375SStefano Zampini         PetscInt parent;
9828135c375SStefano Zampini 
9833e6c54aaSStefano Zampini         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
984*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetDepthLabel(dm,&dlabel));
985*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetTreeParent(dm,p,&parent,NULL));
9868135c375SStefano Zampini         if (parent != p) {
98796ca5757SLisandro Dalcin           PetscInt       vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
9883924b612SStefano Zampini           PetscInt       i,nv,ssize,n,numChildren,depth = -1;
9898135c375SStefano Zampini           const PetscInt *children;
9903924b612SStefano Zampini 
991*5f80ce2aSJacob Faibussowitsch           CHKERRQ(DMPlexGetConeSize(dm,parent,&ssize));
9923924b612SStefano Zampini           switch (ssize) {
9938135c375SStefano Zampini             case 2: /* edge */
9948135c375SStefano Zampini               nv   = 0;
995*5f80ce2aSJacob Faibussowitsch               CHKERRQ(DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids));
996*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D",p-vStart));
9978135c375SStefano Zampini               for (i=0;i<nv;i++) {
998*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(PetscViewerASCIIPrintf(viewer," %D",vids[i]));
9998135c375SStefano Zampini               }
1000*5f80ce2aSJacob Faibussowitsch               CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n"));
10018135c375SStefano Zampini               vp--;
10028135c375SStefano Zampini               break;
10038135c375SStefano Zampini             case 4: /* face */
1004*5f80ce2aSJacob Faibussowitsch               CHKERRQ(DMPlexGetTreeChildren(dm,parent,&numChildren,&children));
10058135c375SStefano Zampini               for (n=0;n<numChildren;n++) {
1006*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(DMLabelGetValue(dlabel,children[n],&depth));
10078135c375SStefano Zampini                 if (!depth) {
10088135c375SStefano Zampini                   const PetscInt *hvsupp,*hesupp,*cone;
10098135c375SStefano Zampini                   PetscInt       hvsuppSize,hesuppSize,coneSize;
1010451a39c7SStefano Zampini                   PetscInt       hv = children[n],he = -1,f;
10118135c375SStefano Zampini 
1012*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(PetscArrayzero(skip,maxsupp));
1013*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(DMPlexGetSupportSize(dm,hv,&hvsuppSize));
1014*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(DMPlexGetSupport(dm,hv,&hvsupp));
10158135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
10168135c375SStefano Zampini                     PetscInt ep;
1017*5f80ce2aSJacob Faibussowitsch                     CHKERRQ(DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL));
10188135c375SStefano Zampini                     if (ep != hvsupp[i]) {
10198135c375SStefano Zampini                       he = hvsupp[i];
10208135c375SStefano Zampini                     } else {
10218135c375SStefano Zampini                       skip[i] = PETSC_TRUE;
10228135c375SStefano Zampini                     }
10238135c375SStefano Zampini                   }
10242c71b3e2SJacob Faibussowitsch                   PetscCheckFalse(he == -1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
1025*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(DMPlexGetCone(dm,he,&cone));
102696ca5757SLisandro Dalcin                   vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
1027*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(DMPlexGetSupportSize(dm,he,&hesuppSize));
1028*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(DMPlexGetSupport(dm,he,&hesupp));
10298135c375SStefano Zampini                   for (f=0;f<hesuppSize;f++) {
10308135c375SStefano Zampini                     PetscInt j;
10318135c375SStefano Zampini 
1032*5f80ce2aSJacob Faibussowitsch                     CHKERRQ(DMPlexGetCone(dm,hesupp[f],&cone));
1033*5f80ce2aSJacob Faibussowitsch                     CHKERRQ(DMPlexGetConeSize(dm,hesupp[f],&coneSize));
10348135c375SStefano Zampini                     for (j=0;j<coneSize;j++) {
10358135c375SStefano Zampini                       PetscInt k;
10368135c375SStefano Zampini                       for (k=0;k<hvsuppSize;k++) {
10378135c375SStefano Zampini                         if (hvsupp[k] == cone[j]) {
10388135c375SStefano Zampini                           skip[k] = PETSC_TRUE;
10398135c375SStefano Zampini                           break;
10408135c375SStefano Zampini                         }
10418135c375SStefano Zampini                       }
10428135c375SStefano Zampini                     }
10438135c375SStefano Zampini                   }
10448135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
10458135c375SStefano Zampini                     if (!skip[i]) {
1046*5f80ce2aSJacob Faibussowitsch                       CHKERRQ(DMPlexGetCone(dm,hvsupp[i],&cone));
104796ca5757SLisandro Dalcin                       vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
10488135c375SStefano Zampini                     }
10498135c375SStefano Zampini                   }
1050*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D",hv-vStart));
10518135c375SStefano Zampini                   for (i=0;i<2;i++) {
1052*5f80ce2aSJacob Faibussowitsch                     CHKERRQ(PetscViewerASCIIPrintf(viewer," %D",vids[i]-vStart));
10538135c375SStefano Zampini                   }
1054*5f80ce2aSJacob Faibussowitsch                   CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n"));
10558135c375SStefano Zampini                   vp--;
10568135c375SStefano Zampini                 }
10578135c375SStefano Zampini               }
10588135c375SStefano Zampini               break;
10598135c375SStefano Zampini             default:
106098921bdaSJacob Faibussowitsch               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize);
10618135c375SStefano Zampini           }
10628135c375SStefano Zampini         }
10638135c375SStefano Zampini       }
1064*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(skip));
10658135c375SStefano Zampini     }
10662c71b3e2SJacob Faibussowitsch     PetscCheckFalse(vp,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
10678135c375SStefano Zampini   }
1068*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&pown));
1069*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBTDestroy(&vown));
10708135c375SStefano Zampini 
10718135c375SStefano Zampini   /* vertices */
107277eacf09SStefano Zampini   if (hovec) { /* higher-order meshes */
107377eacf09SStefano Zampini     const char *fec;
10740286d493SLisandro Dalcin     PetscInt   i,n,s;
107577eacf09SStefano Zampini 
1076*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
1077*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart));
1078*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"nodes\n"));
1079*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetName((PetscObject)hovec,&fec));
1080*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n"));
1081*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%s\n",fec));
1082*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim));
1083*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n")); /*Ordering::byVDIM*/
1084*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(hovec,&array));
1085*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(hovec,&n));
10862c71b3e2SJacob Faibussowitsch     PetscCheckFalse(n%sdim,PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
108777eacf09SStefano Zampini     for (i=0;i<n/sdim;i++) {
108877eacf09SStefano Zampini       for (s=0;s<sdim;s++) {
1089*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPrintf(viewer,fmt,(double) PetscRealPart(array[i*sdim+s])));
109077eacf09SStefano Zampini       }
1091*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n"));
109277eacf09SStefano Zampini     }
1093*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(hovec,&array));
109477eacf09SStefano Zampini   } else {
1095*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(coordinates,&nvert));
1096*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
1097*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim));
1098*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer,"%D\n",sdim));
1099*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(coordinates,&array));
1100cc0d3ed7SStefano Zampini     for (p=0;p<nvert/sdim;p++) {
1101cc0d3ed7SStefano Zampini       PetscInt s;
1102cc0d3ed7SStefano Zampini       for (s=0;s<sdim;s++) {
11033e96840aSStefano Zampini         PetscReal v = PetscRealPart(array[p*sdim+s]);
11043e96840aSStefano Zampini 
1105*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : (double) v));
11068135c375SStefano Zampini       }
1107*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n"));
11088135c375SStefano Zampini     }
1109*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(coordinates,&array));
111077eacf09SStefano Zampini   }
1111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&hovec));
11128135c375SStefano Zampini   PetscFunctionReturn(0);
11138135c375SStefano Zampini }
11148135c375SStefano Zampini 
11150286d493SLisandro Dalcin PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
11168135c375SStefano Zampini {
11178135c375SStefano Zampini   PetscFunctionBegin;
1118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII));
11198135c375SStefano Zampini   PetscFunctionReturn(0);
11208135c375SStefano Zampini }
1121