xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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++) {
219566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&ctx->scctx[i]));
228135c375SStefano Zampini   }
239566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx->scctx));
249566063dSJacob Faibussowitsch   PetscCall(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++) {
359566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD));
369566063dSJacob Faibussowitsch     PetscCall(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;
589566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd));
609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum));
62b135d7daSStefano Zampini   if (!globalNum) {
639566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum));
649566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum));
659566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)globalNum));
66b135d7daSStefano Zampini   }
679566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(globalNum,&gNum));
689566063dSJacob Faibussowitsch   PetscCall(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++;
749566063dSJacob Faibussowitsch       PetscCall(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)) {
779566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(vown,points[i]-vStart));
788135c375SStefano Zampini         }
798135c375SStefano Zampini       }
809566063dSJacob Faibussowitsch       PetscCall(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 
859566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dm,&xlocal));
869566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xlocal,&totdofs));
879566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm,&s));
889566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s,&nfields));
898135c375SStefano Zampini   for (f=0,maxfields=0;f<nfields;f++) {
908135c375SStefano Zampini     PetscInt bs;
918135c375SStefano Zampini 
929566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetFieldComponents(s,f,&bs));
938135c375SStefano Zampini     maxfields += bs;
948135c375SStefano Zampini   }
959566063dSJacob Faibussowitsch   PetscCall(PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield));
969566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
979566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(maxfields,&ctx->scctx));
989566063dSJacob Faibussowitsch   PetscCall(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 
1069566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldName(s,f,&fname));
1079566063dSJacob Faibussowitsch       PetscCall(PetscStrlen(fname,&len));
1088135c375SStefano Zampini       if (len) {
1099566063dSJacob Faibussowitsch         PetscCall(PetscStrcpy(name,fname));
1108135c375SStefano Zampini       } else {
111*63a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(name,256,"Field%" PetscInt_FMT,f));
1128135c375SStefano Zampini       }
1139566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds,f,&disc));
1148135c375SStefano Zampini       if (disc) {
1158135c375SStefano Zampini         PetscClassId id;
1168135c375SStefano Zampini         PetscInt     Nc;
1178135c375SStefano Zampini         char         fec[64];
1188135c375SStefano Zampini 
1199566063dSJacob Faibussowitsch         PetscCall(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 
1279566063dSJacob Faibussowitsch           PetscCall(PetscFEGetNumComponents(fem,&Nc));
1289566063dSJacob Faibussowitsch           PetscCall(PetscFEGetDualSpace(fem,&sp));
1299566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetType(sp,&spname));
1309566063dSJacob Faibussowitsch           PetscCall(PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag));
13128b400f6SJacob Faibussowitsch           PetscCheck(islag,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
1329566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceLagrangeGetContinuity(sp,&continuous));
1339566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetOrder(sp,&order));
13428d58a37SPierre Jolivet           if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */
135*63a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%" PetscInt_FMT "D_P1",dim));
1368135c375SStefano Zampini           } else {
137*63a3b9bcSJacob Faibussowitsch             PetscCheck(continuous || !order,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %" PetscInt_FMT,order);
1388135c375SStefano Zampini             H1   = PETSC_FALSE;
139*63a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%" PetscInt_FMT "D_P%" PetscInt_FMT,dim,order));
1408135c375SStefano Zampini           }
1419566063dSJacob Faibussowitsch           PetscCall(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;
1469566063dSJacob Faibussowitsch             PetscCall(PetscStrallocpy(fec,&fec_type[ctx->nf]));
1479566063dSJacob Faibussowitsch             PetscCall(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;
1529566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetFieldOffset(s,i+vStart,f,&off));
1538135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1548135c375SStefano Zampini             }
1559566063dSJacob Faibussowitsch             PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield));
1568135c375SStefano Zampini           } else {
1578135c375SStefano Zampini             nlocal[ctx->nf] = Nc * totc;
1589566063dSJacob Faibussowitsch             PetscCall(PetscStrallocpy(fec,&fec_type[ctx->nf]));
1599566063dSJacob Faibussowitsch             PetscCall(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;
1649566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetFieldOffset(s,i+cStart,f,&off));
1658135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1668135c375SStefano Zampini             }
1679566063dSJacob Faibussowitsch             PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield));
1688135c375SStefano Zampini           }
1699566063dSJacob Faibussowitsch           PetscCall(VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]));
1709566063dSJacob Faibussowitsch           PetscCall(VecDestroy(&xfield));
1719566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&isfield));
1728135c375SStefano Zampini           ctx->nf++;
1738135c375SStefano Zampini         } else if (id == PETSCFV_CLASSID) {
1748135c375SStefano Zampini           PetscInt c;
1758135c375SStefano Zampini 
1769566063dSJacob Faibussowitsch           PetscCall(PetscFVGetNumComponents((PetscFV)disc,&Nc));
177*63a3b9bcSJacob Faibussowitsch           PetscCall(PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%" PetscInt_FMT "D_P0",dim));
1788135c375SStefano Zampini           for (c = 0; c < Nc; c++) {
1798135c375SStefano Zampini             char comp[256];
180*63a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(comp,256,"%s-Comp%" PetscInt_FMT,name,c));
1819566063dSJacob Faibussowitsch             PetscCall(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;
1859566063dSJacob Faibussowitsch             PetscCall(PetscStrallocpy(fec,&fec_type[ctx->nf]));
1869566063dSJacob Faibussowitsch             PetscCall(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;
1919566063dSJacob Faibussowitsch               PetscCall(PetscSectionGetFieldOffset(s,i+cStart,f,&off));
1928135c375SStefano Zampini               idxs[cum++] = off + c;
1938135c375SStefano Zampini             }
1949566063dSJacob Faibussowitsch             PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield));
1959566063dSJacob Faibussowitsch             PetscCall(VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]));
1969566063dSJacob Faibussowitsch             PetscCall(VecDestroy(&xfield));
1979566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&isfield));
1988135c375SStefano Zampini             ctx->nf++;
1998135c375SStefano Zampini           }
200*63a3b9bcSJacob Faibussowitsch         } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %" PetscInt_FMT,f);
201*63a3b9bcSJacob Faibussowitsch       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %" PetscInt_FMT,f);
2028135c375SStefano Zampini     }
2038135c375SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
2049566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&vown));
2059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xlocal));
2069566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(globalNum,&gNum));
2078135c375SStefano Zampini 
2084cac2994SStefano Zampini   /* create work vectors */
2094cac2994SStefano Zampini   for (f=0;f<ctx->nf;f++) {
2109566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]));
2119566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]));
2129566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(Ufield[f],bs[f]));
2139566063dSJacob Faibussowitsch     PetscCall(VecSetDM(Ufield[f],dm));
2144cac2994SStefano Zampini   }
2154cac2994SStefano Zampini 
2168135c375SStefano Zampini   /* customize the viewer */
2179566063dSJacob Faibussowitsch   PetscCall(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++) {
2199566063dSJacob Faibussowitsch     PetscCall(PetscFree(fieldname[f]));
2209566063dSJacob Faibussowitsch     PetscCall(PetscFree(fec_type[f]));
2219566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Ufield[f]));
2228135c375SStefano Zampini   }
2239566063dSJacob Faibussowitsch   PetscCall(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;
2459566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm,&dlabel));
2469566063dSJacob Faibussowitsch   PetscCall(DMLabelGetValue(dlabel,p,&pdepth));
2479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm,p,&csize));
2489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm,&depth));
2499566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
2508135c375SStefano Zampini   if (label) {
2519566063dSJacob Faibussowitsch     PetscCall(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 */
255*63a3b9bcSJacob Faibussowitsch     PetscCheckFalse(dim < 0 || dim > 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %" PetscInt_FMT,dim);
256*63a3b9bcSJacob Faibussowitsch     PetscCheck(csize <= 8,PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %" PetscInt_FMT " for point %" PetscInt_FMT,csize,p);
257*63a3b9bcSJacob Faibussowitsch     PetscCheck(depth == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %" PetscInt_FMT " for point %" PetscInt_FMT ". You should interpolate the mesh first",depth,p);
258044a5661SStefano Zampini     *cid = mfem_table_cid_unint[dim][csize];
259044a5661SStefano Zampini   } else {
260*63a3b9bcSJacob Faibussowitsch     PetscCheck(csize <= 6,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %" PetscInt_FMT " for point %" PetscInt_FMT,csize,p);
261*63a3b9bcSJacob Faibussowitsch     PetscCheckFalse(pdepth < 0 || pdepth > 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %" PetscInt_FMT " for point %" PetscInt_FMT,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;
2729566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd));
2739566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
274cc0d3ed7SStefano Zampini   sdim = dim;
275cc0d3ed7SStefano Zampini   if (csec) {
27684f354e3SLisandro Dalcin     PetscInt sStart,sEnd;
27784f354e3SLisandro Dalcin 
2789566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(dm,&sdim));
2799566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(csec,&sStart,&sEnd));
2809566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(csec,vStart,&off));
281cc0d3ed7SStefano Zampini     off  = off/sdim;
28284f354e3SLisandro Dalcin     if (p >= sStart && p < sEnd) {
2839566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(csec,p,&dof));
284cc0d3ed7SStefano Zampini     }
28584f354e3SLisandro Dalcin   }
286cc0d3ed7SStefano Zampini   if (!dof) {
2879566063dSJacob Faibussowitsch     PetscCall(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;
2919566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points));
292cc0d3ed7SStefano Zampini   } else {
2939566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(csec,p,&off));
2949566063dSJacob Faibussowitsch     PetscCall(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 
3010c2bc6bfSStefano Zampini static PetscErrorCode GLVisCreateFE(PetscFE femIn,char name[32],PetscFE *fem,IS *perm)
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;
3169566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)femIn, &comm));
3179566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(femIn,&P));
3189566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(femIn,&Q));
3199566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(Q,&K));
3209566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(K,&dim));
3219566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(P,&deg,NULL));
3229566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(P,&dof));
3239566063dSJacob Faibussowitsch   PetscCall(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;
3320c2bc6bfSStefano Zampini   if (isSimplex) deg = PetscMin(deg,3); /* Permutation not coded for degree higher than 3 */
333066ea43fSLisandro Dalcin   /* Create space */
3349566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreate(comm,&P));
3359566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL));
3369566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialSetTensor(P,isTensor));
3379566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumComponents(P,dof));
3389566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumVariables(P,dim));
3399566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetDegree(P,deg,PETSC_DETERMINE));
3409566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetUp(P));
341066ea43fSLisandro Dalcin   /* Create dual space */
3429566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceCreate(comm,&Q));
3439566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE));
3449566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTensor(Q,isTensor));
3459566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetContinuity(Q,continuity));
3469566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetNodeType(Q,nodeType,endpoint,0));
3479566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetNumComponents(Q,dof));
3489566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetOrder(Q,deg));
3499566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
3509566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(Q,K));
3519566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&K));
3529566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(Q));
353066ea43fSLisandro Dalcin   /* Create quadrature */
354066ea43fSLisandro Dalcin   if (isSimplex) {
3559566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(dim,  1,deg+1,-1,+1,&q));
3569566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(dim-1,1,deg+1,-1,+1,&fq));
357066ea43fSLisandro Dalcin   } else {
3589566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim,  1,deg+1,-1,+1,&q));
3599566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim-1,1,deg+1,-1,+1,&fq));
360066ea43fSLisandro Dalcin   }
361066ea43fSLisandro Dalcin   /* Create finite element */
3629566063dSJacob Faibussowitsch   PetscCall(PetscFECreate(comm,fem));
363*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name,32,"L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT,dim,deg));
3649566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*fem,name));
3659566063dSJacob Faibussowitsch   PetscCall(PetscFESetType(*fem,PETSCFEBASIC));
3669566063dSJacob Faibussowitsch   PetscCall(PetscFESetNumComponents(*fem,dof));
3679566063dSJacob Faibussowitsch   PetscCall(PetscFESetBasisSpace(*fem,P));
3689566063dSJacob Faibussowitsch   PetscCall(PetscFESetDualSpace(*fem,Q));
3699566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(*fem,q));
3709566063dSJacob Faibussowitsch   PetscCall(PetscFESetFaceQuadrature(*fem,fq));
3719566063dSJacob Faibussowitsch   PetscCall(PetscFESetUp(*fem));
3720c2bc6bfSStefano Zampini 
3730c2bc6bfSStefano Zampini   /* Both MFEM and PETSc are lexicographic, but PLEX stores the swapped cone */
3740c2bc6bfSStefano Zampini   *perm = NULL;
3750c2bc6bfSStefano Zampini   if (isSimplex && dim == 3) {
3760c2bc6bfSStefano Zampini     PetscInt celldofs,*pidx;
3770c2bc6bfSStefano Zampini 
3780c2bc6bfSStefano Zampini     PetscCall(PetscDualSpaceGetDimension(Q,&celldofs));
3790c2bc6bfSStefano Zampini     celldofs /= dof;
3800c2bc6bfSStefano Zampini     PetscCall(PetscMalloc1(celldofs,&pidx));
3810c2bc6bfSStefano Zampini     switch (celldofs) {
3820c2bc6bfSStefano Zampini     case 4:
3830c2bc6bfSStefano Zampini       pidx[0] = 2;
3840c2bc6bfSStefano Zampini       pidx[1] = 0;
3850c2bc6bfSStefano Zampini       pidx[2] = 1;
3860c2bc6bfSStefano Zampini       pidx[3] = 3;
3870c2bc6bfSStefano Zampini     break;
3880c2bc6bfSStefano Zampini     case 10:
3890c2bc6bfSStefano Zampini       pidx[0] = 5;
3900c2bc6bfSStefano Zampini       pidx[1] = 3;
3910c2bc6bfSStefano Zampini       pidx[2] = 0;
3920c2bc6bfSStefano Zampini       pidx[3] = 4;
3930c2bc6bfSStefano Zampini       pidx[4] = 1;
3940c2bc6bfSStefano Zampini       pidx[5] = 2;
3950c2bc6bfSStefano Zampini       pidx[6] = 8;
3960c2bc6bfSStefano Zampini       pidx[7] = 6;
3970c2bc6bfSStefano Zampini       pidx[8] = 7;
3980c2bc6bfSStefano Zampini       pidx[9] = 9;
3990c2bc6bfSStefano Zampini     break;
4000c2bc6bfSStefano Zampini     case 20:
4010c2bc6bfSStefano Zampini       pidx[ 0] = 9;
4020c2bc6bfSStefano Zampini       pidx[ 1] = 7;
4030c2bc6bfSStefano Zampini       pidx[ 2] = 4;
4040c2bc6bfSStefano Zampini       pidx[ 3] = 0;
4050c2bc6bfSStefano Zampini       pidx[ 4] = 8;
4060c2bc6bfSStefano Zampini       pidx[ 5] = 5;
4070c2bc6bfSStefano Zampini       pidx[ 6] = 1;
4080c2bc6bfSStefano Zampini       pidx[ 7] = 6;
4090c2bc6bfSStefano Zampini       pidx[ 8] = 2;
4100c2bc6bfSStefano Zampini       pidx[ 9] = 3;
4110c2bc6bfSStefano Zampini       pidx[10] = 15;
4120c2bc6bfSStefano Zampini       pidx[11] = 13;
4130c2bc6bfSStefano Zampini       pidx[12] = 10;
4140c2bc6bfSStefano Zampini       pidx[13] = 14;
4150c2bc6bfSStefano Zampini       pidx[14] = 11;
4160c2bc6bfSStefano Zampini       pidx[15] = 12;
4170c2bc6bfSStefano Zampini       pidx[16] = 18;
4180c2bc6bfSStefano Zampini       pidx[17] = 16;
4190c2bc6bfSStefano Zampini       pidx[18] = 17;
4200c2bc6bfSStefano Zampini       pidx[19] = 19;
4210c2bc6bfSStefano Zampini     break;
4220c2bc6bfSStefano Zampini     default:
423*63a3b9bcSJacob Faibussowitsch       SETERRQ(comm,PETSC_ERR_SUP,"Unhandled degree,dof pair %" PetscInt_FMT ",%" PetscInt_FMT,deg,celldofs);
4240c2bc6bfSStefano Zampini     break;
4250c2bc6bfSStefano Zampini     }
4260c2bc6bfSStefano Zampini     PetscCall(ISCreateBlock(PETSC_COMM_SELF,dof,celldofs,pidx,PETSC_OWN_POINTER,perm));
4270c2bc6bfSStefano Zampini   }
4280c2bc6bfSStefano Zampini 
429066ea43fSLisandro Dalcin   /* Cleanup */
4309566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&P));
4319566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&Q));
4329566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&q));
4339566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&fq));
434066ea43fSLisandro Dalcin   PetscFunctionReturn(0);
435066ea43fSLisandro Dalcin }
436066ea43fSLisandro Dalcin 
43777eacf09SStefano Zampini /*
43877eacf09SStefano Zampini    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
43977eacf09SStefano Zampini    Higher order meshes are also supported
44077eacf09SStefano Zampini */
4418135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
4428135c375SStefano Zampini {
4438135c375SStefano Zampini   DMLabel              label;
4440c2bc6bfSStefano Zampini   PetscSection         coordSection,parentSection,hoSection = NULL;
44577eacf09SStefano Zampini   Vec                  coordinates,hovec;
4468135c375SStefano Zampini   const PetscScalar    *array;
447f86f7544SStefano Zampini   PetscInt             bf,p,sdim,dim,depth,novl,minl;
448412e9a14SMatthew G. Knepley   PetscInt             cStart,cEnd,vStart,vEnd,nvert;
4493924b612SStefano Zampini   PetscMPIInt          size;
4503e6c54aaSStefano Zampini   PetscBool            localized,isascii;
45128d58a37SPierre Jolivet   PetscBool            enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE;
4523e6c54aaSStefano Zampini   PetscBT              pown,vown;
4538135c375SStefano Zampini   PetscContainer       glvis_container;
454044a5661SStefano Zampini   PetscBool            cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
455f86f7544SStefano Zampini   PetscBool            enable_emark,enable_bmark;
45677eacf09SStefano Zampini   const char           *fmt;
4577bf4dd16SStefano Zampini   char                 emark[64] = "",bmark[64] = "";
4588135c375SStefano Zampini 
4598135c375SStefano Zampini   PetscFunctionBegin;
4608135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4618135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
46328b400f6SJacob Faibussowitsch   PetscCheck(isascii,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
4649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size));
46508401ef6SPierre Jolivet   PetscCheck(size <= 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
4669566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
4670c2bc6bfSStefano Zampini   PetscCall(DMPlexGetDepth(dm,&depth));
4688135c375SStefano Zampini 
4698135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container));
47128b400f6SJacob Faibussowitsch   PetscCheck(glvis_container,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
4728135c375SStefano Zampini   {
4738135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
4749566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(glvis_container,(void**)&glvis_info));
4758135c375SStefano Zampini     enabled = glvis_info->enabled;
47677eacf09SStefano Zampini     fmt     = glvis_info->fmt;
4778135c375SStefano Zampini   }
47821414b21SStefano Zampini 
4790c2bc6bfSStefano Zampini   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh */
4809566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec));
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)hovec));
482066ea43fSLisandro Dalcin   if (!hovec) {
483066ea43fSLisandro Dalcin     DM           cdm;
484066ea43fSLisandro Dalcin     PetscFE      disc;
485066ea43fSLisandro Dalcin     PetscClassId classid;
486066ea43fSLisandro Dalcin 
4879566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm,&cdm));
4889566063dSJacob Faibussowitsch     PetscCall(DMGetField(cdm,0,NULL,(PetscObject*)&disc));
4899566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId((PetscObject)disc,&classid));
490066ea43fSLisandro Dalcin     if (classid == PETSCFE_CLASSID) {
491066ea43fSLisandro Dalcin       DM      hocdm;
492066ea43fSLisandro Dalcin       PetscFE hodisc;
493066ea43fSLisandro Dalcin       Vec     vec;
494066ea43fSLisandro Dalcin       Mat     mat;
495066ea43fSLisandro Dalcin       char    name[32],fec_type[64];
4960c2bc6bfSStefano Zampini       IS      perm = NULL;
497066ea43fSLisandro Dalcin 
4980c2bc6bfSStefano Zampini       PetscCall(GLVisCreateFE(disc,name,&hodisc,&perm));
4999566063dSJacob Faibussowitsch       PetscCall(DMClone(cdm,&hocdm));
5009566063dSJacob Faibussowitsch       PetscCall(DMSetField(hocdm,0,NULL,(PetscObject)hodisc));
5019566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&hodisc));
5029566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(hocdm));
503066ea43fSLisandro Dalcin 
5049566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(dm,&vec));
5059566063dSJacob Faibussowitsch       PetscCall(DMCreateGlobalVector(hocdm,&hovec));
5069566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation(cdm,hocdm,&mat,NULL));
5079566063dSJacob Faibussowitsch       PetscCall(MatInterpolate(mat,vec,hovec));
5089566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mat));
5090c2bc6bfSStefano Zampini       PetscCall(DMGetLocalSection(hocdm,&hoSection));
5100c2bc6bfSStefano Zampini       PetscCall(PetscSectionSetClosurePermutation(hoSection, (PetscObject)hocdm, depth, perm));
5110c2bc6bfSStefano Zampini       PetscCall(ISDestroy(&perm));
5129566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&hocdm));
5139566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name));
5149566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)hovec,fec_type));
515066ea43fSLisandro Dalcin     }
516066ea43fSLisandro Dalcin   }
51721414b21SStefano Zampini 
5189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
5199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGhostCellStratum(dm,&p,NULL));
520c3c203b2SStefano Zampini   if (p >= 0) cEnd = p;
5219566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd));
5229566063dSJacob Faibussowitsch   PetscCall(DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL));
5239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm,&localized));
5242c71b3e2SJacob Faibussowitsch   PetscCheckFalse(periodic && !localized && !hovec,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
5259566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm,&coordSection));
5269566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm,&sdim));
5279566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm,&coordinates));
52808401ef6SPierre Jolivet   PetscCheck(coordinates || hovec,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");
5298135c375SStefano Zampini 
5308135c375SStefano Zampini   /*
5318135c375SStefano Zampini      a couple of sections of the mesh specification are disabled
5323e96840aSStefano 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)
53377eacf09SStefano Zampini        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
5343e6c54aaSStefano Zampini                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
5358135c375SStefano Zampini   */
5363e96840aSStefano Zampini   enable_boundary = PETSC_FALSE;
5378135c375SStefano Zampini   enable_ncmesh   = PETSC_FALSE;
5383e6c54aaSStefano Zampini   enable_mfem     = PETSC_FALSE;
539f86f7544SStefano Zampini   enable_emark    = PETSC_FALSE;
540f86f7544SStefano Zampini   enable_bmark    = PETSC_FALSE;
5417bf4dd16SStefano Zampini   /* I'm tired of problems with negative values in the markers, disable them */
542d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");
5439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL));
5449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL));
5459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL));
5469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL));
5479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark));
5489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark));
549d0609cedSBarry Smith   PetscOptionsEnd();
550f86f7544SStefano Zampini   if (enable_bmark) enable_boundary = PETSC_TRUE;
551f86f7544SStefano Zampini 
5529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size));
5532c71b3e2SJacob Faibussowitsch   PetscCheckFalse(enable_ncmesh && size > 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
5542c71b3e2SJacob Faibussowitsch   PetscCheckFalse(enable_boundary && depth >= 0 && dim != depth,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
5557e1aca4eSStefano Zampini                                                              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
5562c71b3e2SJacob Faibussowitsch   PetscCheckFalse(enable_ncmesh && depth >= 0 && dim != depth,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
5577e1aca4eSStefano Zampini                                                            "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
558044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
559*63a3b9bcSJacob Faibussowitsch     PetscCheck(depth == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %" PetscInt_FMT ". You should interpolate the mesh first",depth);
560044a5661SStefano Zampini     cellvertex = PETSC_TRUE;
561044a5661SStefano Zampini   }
5628135c375SStefano Zampini 
5638135c375SStefano Zampini   /* Identify possible cells in the overlap */
5648135c375SStefano Zampini   novl = 0;
5658135c375SStefano Zampini   pown = NULL;
5663924b612SStefano Zampini   if (size > 1) {
5673e6c54aaSStefano Zampini     IS             globalNum = NULL;
5683e6c54aaSStefano Zampini     const PetscInt *gNum;
5693e6c54aaSStefano Zampini     PetscBool      ovl  = PETSC_FALSE;
5703e6c54aaSStefano Zampini 
5719566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum));
572b135d7daSStefano Zampini     if (!globalNum) {
57328d58a37SPierre Jolivet       if (view_ovl) {
5749566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum));
57528d58a37SPierre Jolivet       } else {
5769566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum));
57728d58a37SPierre Jolivet       }
5789566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum));
5799566063dSJacob Faibussowitsch       PetscCall(PetscObjectDereference((PetscObject)globalNum));
580b135d7daSStefano Zampini     }
5819566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(globalNum,&gNum));
5828135c375SStefano Zampini     for (p=cStart; p<cEnd; p++) {
5838135c375SStefano Zampini       if (gNum[p-cStart] < 0) {
5848135c375SStefano Zampini         ovl = PETSC_TRUE;
5858135c375SStefano Zampini         novl++;
5868135c375SStefano Zampini       }
5878135c375SStefano Zampini     }
5888135c375SStefano Zampini     if (ovl) {
5898135c375SStefano Zampini       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
5908135c375SStefano Zampini          TODO: garbage collector? attach pown to dm?  */
5919566063dSJacob Faibussowitsch       PetscCall(PetscBTCreate(cEnd-cStart,&pown));
5923e6c54aaSStefano Zampini       for (p=cStart; p<cEnd; p++) {
5933e6c54aaSStefano Zampini         if (gNum[p-cStart] < 0) continue;
5943e6c54aaSStefano Zampini         else {
5959566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(pown,p-cStart));
5968135c375SStefano Zampini         }
5978135c375SStefano Zampini       }
5983e6c54aaSStefano Zampini     }
5999566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(globalNum,&gNum));
6003e6c54aaSStefano Zampini   }
6018135c375SStefano Zampini 
6021d4815f0SStefano Zampini   /* vertex_parents (Non-conforming meshes) */
6031d4815f0SStefano Zampini   parentSection  = NULL;
6041d4815f0SStefano Zampini   if (enable_ncmesh) {
6059566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL));
6061d4815f0SStefano Zampini     enable_ncmesh = (PetscBool)(enable_ncmesh && parentSection);
6071d4815f0SStefano Zampini   }
6083e6c54aaSStefano Zampini   /* return if this process is disabled */
6098135c375SStefano Zampini   if (!enabled) {
6109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"MFEM mesh %s\n",enable_ncmesh ? "v1.1" : "v1.0"));
6119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\ndimension\n"));
612*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",dim));
6139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nelements\n"));
614*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
6159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nboundary\n"));
616*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
6179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
618*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
619*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",sdim));
6209566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&pown));
6219566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&hovec));
6228135c375SStefano Zampini     PetscFunctionReturn(0);
6238135c375SStefano Zampini   }
6248135c375SStefano Zampini 
6253e6c54aaSStefano Zampini   if (enable_mfem) {
6263e6c54aaSStefano Zampini     if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
6273e6c54aaSStefano Zampini       PetscInt    vpc = 0;
6283e6c54aaSStefano Zampini       char        fec[64];
62996ca5757SLisandro Dalcin       PetscInt    vids[8] = {0,1,2,3,4,5,6,7};
63096ca5757SLisandro Dalcin       PetscInt    hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3};
63196ca5757SLisandro Dalcin       PetscInt    quadv[8] = {0,1,3,2}, triv[3] = {0,1,2};
63296ca5757SLisandro Dalcin       PetscInt    *dof = NULL;
6333e6c54aaSStefano Zampini       PetscScalar *array,*ptr;
6343e6c54aaSStefano Zampini 
635*63a3b9bcSJacob Faibussowitsch       PetscCall(PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P1",dim));
6363e6c54aaSStefano Zampini       if (cEnd-cStart) {
6373e6c54aaSStefano Zampini         PetscInt fpc;
6383e6c54aaSStefano Zampini 
6399566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(dm,cStart,&fpc));
6403e6c54aaSStefano Zampini         switch(dim) {
6413e6c54aaSStefano Zampini           case 1:
6423e6c54aaSStefano Zampini             vpc = 2;
6433e6c54aaSStefano Zampini             dof = hexv;
6443e6c54aaSStefano Zampini             break;
6453e6c54aaSStefano Zampini           case 2:
6463e6c54aaSStefano Zampini             switch (fpc) {
6473e6c54aaSStefano Zampini               case 3:
6483e6c54aaSStefano Zampini                 vpc = 3;
649044a5661SStefano Zampini                 dof = triv;
6503e6c54aaSStefano Zampini                 break;
6513e6c54aaSStefano Zampini               case 4:
6523e6c54aaSStefano Zampini                 vpc = 4;
653044a5661SStefano Zampini                 dof = quadv;
6543e6c54aaSStefano Zampini                 break;
6553e6c54aaSStefano Zampini               default:
656*63a3b9bcSJacob Faibussowitsch                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %" PetscInt_FMT,fpc);
6573e6c54aaSStefano Zampini             }
6583e6c54aaSStefano Zampini             break;
6593e6c54aaSStefano Zampini           case 3:
6603e6c54aaSStefano Zampini             switch (fpc) {
661044a5661SStefano Zampini               case 4: /* TODO: still need to understand L2 ordering for tets */
6623e6c54aaSStefano Zampini                 vpc = 4;
663044a5661SStefano Zampini                 dof = tetv;
664044a5661SStefano Zampini                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
6653e6c54aaSStefano Zampini               case 6:
666*63a3b9bcSJacob Faibussowitsch                 PetscCheck(!cellvertex,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %" PetscInt_FMT,fpc);
667044a5661SStefano Zampini                 vpc = 8;
668044a5661SStefano Zampini                 dof = hexv;
669044a5661SStefano Zampini                 break;
670044a5661SStefano Zampini               case 8:
671*63a3b9bcSJacob Faibussowitsch                 PetscCheck(cellvertex,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %" PetscInt_FMT,fpc);
6723e6c54aaSStefano Zampini                 vpc = 8;
6733e6c54aaSStefano Zampini                 dof = hexv;
6743e6c54aaSStefano Zampini                 break;
6753e6c54aaSStefano Zampini               default:
676*63a3b9bcSJacob Faibussowitsch                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %" PetscInt_FMT,fpc);
6773e6c54aaSStefano Zampini             }
6783e6c54aaSStefano Zampini             break;
6793e6c54aaSStefano Zampini           default:
6803e6c54aaSStefano Zampini             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
6813e6c54aaSStefano Zampini         }
6829566063dSJacob Faibussowitsch         PetscCall(DMPlexReorderCell(dm,cStart,vids));
6833e6c54aaSStefano Zampini       }
68428b400f6SJacob Faibussowitsch       PetscCheck(dof,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
6859566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec));
6869566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)hovec,fec));
6879566063dSJacob Faibussowitsch       PetscCall(VecGetArray(hovec,&array));
6883e6c54aaSStefano Zampini       ptr  = array;
6893e6c54aaSStefano Zampini       for (p=cStart;p<cEnd;p++) {
6903e6c54aaSStefano Zampini         PetscInt    csize,v,d;
6913e6c54aaSStefano Zampini         PetscScalar *vals = NULL;
6923e6c54aaSStefano Zampini 
6933e6c54aaSStefano Zampini         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
6949566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals));
695*63a3b9bcSJacob Faibussowitsch         PetscCheckFalse(csize != vpc*sdim && csize != vpc*sdim*2,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported closure size %" PetscInt_FMT " (vpc %" PetscInt_FMT ", sdim %" PetscInt_FMT ")",csize,vpc,sdim);
6963e6c54aaSStefano Zampini         for (v=0;v<vpc;v++) {
6973e6c54aaSStefano Zampini           for (d=0;d<sdim;d++) {
6983e6c54aaSStefano Zampini             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
6993e6c54aaSStefano Zampini           }
7003e6c54aaSStefano Zampini         }
7013e6c54aaSStefano Zampini         ptr += vpc*sdim;
7029566063dSJacob Faibussowitsch         PetscCall(DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals));
7033e6c54aaSStefano Zampini       }
7049566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(hovec,&array));
7053e6c54aaSStefano Zampini     }
7063e6c54aaSStefano Zampini   }
7073e96840aSStefano Zampini   /* if we have high-order coordinates in 3D, we need to specify the boundary */
7083e96840aSStefano Zampini   if (hovec && dim == 3) enable_boundary = PETSC_TRUE;
7093e6c54aaSStefano Zampini 
7108135c375SStefano Zampini   /* header */
7119566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"MFEM mesh %s\n",enable_ncmesh ? "v1.1" : "v1.0"));
7128135c375SStefano Zampini 
7138135c375SStefano Zampini   /* topological dimension */
7149566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\ndimension\n"));
715*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",dim));
7168135c375SStefano Zampini 
7178135c375SStefano Zampini   /* elements */
718f86f7544SStefano Zampini   minl = 1;
719f86f7544SStefano Zampini   label = NULL;
720f86f7544SStefano Zampini   if (enable_emark) {
721f86f7544SStefano Zampini     PetscInt lminl = PETSC_MAX_INT;
722f86f7544SStefano Zampini 
7239566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm,emark,&label));
724f86f7544SStefano Zampini     if (label) {
725f86f7544SStefano Zampini       IS       vals;
726f86f7544SStefano Zampini       PetscInt ldef;
727f86f7544SStefano Zampini 
7289566063dSJacob Faibussowitsch       PetscCall(DMLabelGetDefaultValue(label,&ldef));
7299566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(label,&vals));
7309566063dSJacob Faibussowitsch       PetscCall(ISGetMinMax(vals,&lminl,NULL));
7319566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&vals));
732f86f7544SStefano Zampini       lminl = PetscMin(ldef,lminl);
733f86f7544SStefano Zampini     }
7341c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm)));
735f86f7544SStefano Zampini     if (minl == PETSC_MAX_INT) minl = 1;
736f86f7544SStefano Zampini   }
7379566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\nelements\n"));
738*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",cEnd-cStart-novl));
7398135c375SStefano Zampini   for (p=cStart;p<cEnd;p++) {
74096ca5757SLisandro Dalcin     PetscInt       vids[8];
74111a4995dSStefano Zampini     PetscInt       i,nv = 0,cid = -1,mid = 1;
7428135c375SStefano Zampini 
7433e6c54aaSStefano Zampini     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
7449566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid));
7459566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids));
7469566063dSJacob Faibussowitsch     PetscCall(DMPlexReorderCell(dm,p,vids));
747*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT,mid,cid));
7488135c375SStefano Zampini     for (i=0;i<nv;i++) {
749*63a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT,vids[i]));
7508135c375SStefano Zampini     }
7519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
7528135c375SStefano Zampini   }
7538135c375SStefano Zampini 
754cc0d3ed7SStefano Zampini   /* boundary */
7559566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\nboundary\n"));
756cc0d3ed7SStefano Zampini   if (!enable_boundary) {
757*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
758cc0d3ed7SStefano Zampini   } else {
75977eacf09SStefano Zampini     DMLabel  perLabel;
76077eacf09SStefano Zampini     PetscBT  bfaces;
761b135d7daSStefano Zampini     PetscInt fStart,fEnd,*fcells;
762cc0d3ed7SStefano Zampini 
7639566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm,1,&fStart,&fEnd));
7649566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(fEnd-fStart,&bfaces));
7659566063dSJacob Faibussowitsch     PetscCall(DMPlexGetMaxSizes(dm,NULL,&p));
7669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(p,&fcells));
7679566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm,"glvis_periodic_cut",&perLabel));
7680c2bc6bfSStefano Zampini     if (!perLabel && periodic) { /* this periodic cut can be moved up to DMPlex setup */
7699566063dSJacob Faibussowitsch       PetscCall(DMCreateLabel(dm,"glvis_periodic_cut"));
7709566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm,"glvis_periodic_cut",&perLabel));
7719566063dSJacob Faibussowitsch       PetscCall(DMLabelSetDefaultValue(perLabel,1));
77277eacf09SStefano Zampini       for (p=cStart;p<cEnd;p++) {
773c3c203b2SStefano Zampini         DMPolytopeType cellType;
774c3c203b2SStefano Zampini         PetscInt       dof;
775b135d7daSStefano Zampini 
7769566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(dm,p,&cellType));
7779566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(coordSection,p,&dof));
77877eacf09SStefano Zampini         if (dof) {
779c3c203b2SStefano Zampini           PetscInt    uvpc, v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
78077eacf09SStefano Zampini           PetscScalar *vals = NULL;
781c3c203b2SStefano Zampini 
782c3c203b2SStefano Zampini           uvpc = DMPolytopeTypeGetNumVertices(cellType);
783*63a3b9bcSJacob Faibussowitsch           PetscCheck(dof%sdim == 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %" PetscInt_FMT " and space dimension %" PetscInt_FMT,dof,sdim);
7849566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals));
7859566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure));
78677eacf09SStefano Zampini           for (v=0;v<cellClosureSize;v++)
78777eacf09SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
78877eacf09SStefano Zampini               vidxs = cellClosure + 2*v;
78977eacf09SStefano Zampini               break;
79077eacf09SStefano Zampini             }
79128b400f6SJacob Faibussowitsch           PetscCheck(vidxs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
792b135d7daSStefano Zampini           for (v=0;v<uvpc;v++) {
79377eacf09SStefano Zampini             PetscInt s;
794044a5661SStefano Zampini 
79577eacf09SStefano Zampini             for (s=0;s<sdim;s++) {
796b135d7daSStefano Zampini               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) {
7979566063dSJacob Faibussowitsch                 PetscCall(DMLabelSetValue(perLabel,vidxs[2*v],2));
79877eacf09SStefano Zampini               }
79977eacf09SStefano Zampini             }
80077eacf09SStefano Zampini           }
8019566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure));
8029566063dSJacob Faibussowitsch           PetscCall(DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals));
80377eacf09SStefano Zampini         }
80477eacf09SStefano Zampini       }
80577eacf09SStefano Zampini       if (dim > 1) {
806b135d7daSStefano Zampini         PetscInt eEnd,eStart;
807044a5661SStefano Zampini 
8089566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthStratum(dm,1,&eStart,&eEnd));
80977eacf09SStefano Zampini         for (p=eStart;p<eEnd;p++) {
81077eacf09SStefano Zampini           const PetscInt *cone;
81177eacf09SStefano Zampini           PetscInt       coneSize,i;
81277eacf09SStefano Zampini           PetscBool      ispe = PETSC_TRUE;
81377eacf09SStefano Zampini 
8149566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCone(dm,p,&cone));
8159566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm,p,&coneSize));
81677eacf09SStefano Zampini           for (i=0;i<coneSize;i++) {
81777eacf09SStefano Zampini             PetscInt v;
81877eacf09SStefano Zampini 
8199566063dSJacob Faibussowitsch             PetscCall(DMLabelGetValue(perLabel,cone[i],&v));
82077eacf09SStefano Zampini             ispe = (PetscBool)(ispe && (v==2));
82177eacf09SStefano Zampini           }
82277eacf09SStefano Zampini           if (ispe && coneSize) {
8233e96840aSStefano Zampini             PetscInt       ch, numChildren;
8243e96840aSStefano Zampini             const PetscInt *children;
8253e96840aSStefano Zampini 
8269566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(perLabel,p,2));
8279566063dSJacob Faibussowitsch             PetscCall(DMPlexGetTreeChildren(dm,p,&numChildren,&children));
8283e96840aSStefano Zampini             for (ch = 0; ch < numChildren; ch++) {
8299566063dSJacob Faibussowitsch               PetscCall(DMLabelSetValue(perLabel,children[ch],2));
8303e96840aSStefano Zampini             }
83177eacf09SStefano Zampini           }
83277eacf09SStefano Zampini         }
83377eacf09SStefano Zampini         if (dim > 2) {
83477eacf09SStefano Zampini           for (p=fStart;p<fEnd;p++) {
83577eacf09SStefano Zampini             const PetscInt *cone;
83677eacf09SStefano Zampini             PetscInt       coneSize,i;
83777eacf09SStefano Zampini             PetscBool      ispe = PETSC_TRUE;
83877eacf09SStefano Zampini 
8399566063dSJacob Faibussowitsch             PetscCall(DMPlexGetCone(dm,p,&cone));
8409566063dSJacob Faibussowitsch             PetscCall(DMPlexGetConeSize(dm,p,&coneSize));
84177eacf09SStefano Zampini             for (i=0;i<coneSize;i++) {
84277eacf09SStefano Zampini               PetscInt v;
84377eacf09SStefano Zampini 
8449566063dSJacob Faibussowitsch               PetscCall(DMLabelGetValue(perLabel,cone[i],&v));
84577eacf09SStefano Zampini               ispe = (PetscBool)(ispe && (v==2));
84677eacf09SStefano Zampini             }
84777eacf09SStefano Zampini             if (ispe && coneSize) {
8483e96840aSStefano Zampini               PetscInt       ch, numChildren;
8493e96840aSStefano Zampini               const PetscInt *children;
8503e96840aSStefano Zampini 
8519566063dSJacob Faibussowitsch               PetscCall(DMLabelSetValue(perLabel,p,2));
8529566063dSJacob Faibussowitsch               PetscCall(DMPlexGetTreeChildren(dm,p,&numChildren,&children));
8533e96840aSStefano Zampini               for (ch = 0; ch < numChildren; ch++) {
8549566063dSJacob Faibussowitsch                 PetscCall(DMLabelSetValue(perLabel,children[ch],2));
8553e96840aSStefano Zampini               }
85677eacf09SStefano Zampini             }
85777eacf09SStefano Zampini           }
85877eacf09SStefano Zampini         }
85977eacf09SStefano Zampini       }
86077eacf09SStefano Zampini     }
86177eacf09SStefano Zampini     for (p=fStart;p<fEnd;p++) {
86277eacf09SStefano Zampini       const PetscInt *support;
8638135c375SStefano Zampini       PetscInt       supportSize;
86477eacf09SStefano Zampini       PetscBool      isbf = PETSC_FALSE;
8658135c375SStefano Zampini 
8669566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(dm,p,&supportSize));
8673e6c54aaSStefano Zampini       if (pown) {
8688135c375SStefano Zampini         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
86977eacf09SStefano Zampini         PetscInt  i;
87077eacf09SStefano Zampini 
8719566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(dm,p,&support));
87277eacf09SStefano Zampini         for (i=0;i<supportSize;i++) {
87377eacf09SStefano Zampini           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
87477eacf09SStefano Zampini           else has_ghost = PETSC_TRUE;
87577eacf09SStefano Zampini         }
87677eacf09SStefano Zampini         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
87777eacf09SStefano Zampini       } else {
87877eacf09SStefano Zampini         isbf = (PetscBool)(supportSize == 1);
87977eacf09SStefano Zampini       }
88077eacf09SStefano Zampini       if (!isbf && perLabel) {
88177eacf09SStefano Zampini         const PetscInt *cone;
88277eacf09SStefano Zampini         PetscInt       coneSize,i;
88377eacf09SStefano Zampini 
8849566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(dm,p,&cone));
8859566063dSJacob Faibussowitsch         PetscCall(DMPlexGetConeSize(dm,p,&coneSize));
88677eacf09SStefano Zampini         isbf = PETSC_TRUE;
88777eacf09SStefano Zampini         for (i=0;i<coneSize;i++) {
88877eacf09SStefano Zampini           PetscInt v,d;
88977eacf09SStefano Zampini 
8909566063dSJacob Faibussowitsch           PetscCall(DMLabelGetValue(perLabel,cone[i],&v));
8919566063dSJacob Faibussowitsch           PetscCall(DMLabelGetDefaultValue(perLabel,&d));
89277eacf09SStefano Zampini           isbf = (PetscBool)(isbf && v != d);
89377eacf09SStefano Zampini         }
89477eacf09SStefano Zampini       }
89577eacf09SStefano Zampini       if (isbf) {
8969566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(bfaces,p-fStart));
89777eacf09SStefano Zampini       }
89877eacf09SStefano Zampini     }
89977eacf09SStefano Zampini     /* count boundary faces */
90077eacf09SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
90177eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
90277eacf09SStefano Zampini         const PetscInt *support;
90377eacf09SStefano Zampini         PetscInt       supportSize,c;
9048135c375SStefano Zampini 
9059566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(dm,p,&supportSize));
9069566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(dm,p,&support));
90777eacf09SStefano Zampini         for (c=0;c<supportSize;c++) {
9083e96840aSStefano Zampini           const    PetscInt *cone;
909b135d7daSStefano Zampini           PetscInt cell,cl,coneSize;
9103e96840aSStefano Zampini 
9113e96840aSStefano Zampini           cell = support[c];
9123e96840aSStefano Zampini           if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue;
9139566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCone(dm,cell,&cone));
9149566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm,cell,&coneSize));
915b135d7daSStefano Zampini           for (cl=0;cl<coneSize;cl++) {
9163e96840aSStefano Zampini             if (cone[cl] == p) {
9173e96840aSStefano Zampini               bf += 1;
9183e96840aSStefano Zampini               break;
9198135c375SStefano Zampini             }
92077eacf09SStefano Zampini           }
9213e96840aSStefano Zampini         }
9228135c375SStefano Zampini       }
9238135c375SStefano Zampini     }
924f86f7544SStefano Zampini     minl = 1;
925f86f7544SStefano Zampini     label = NULL;
926f86f7544SStefano Zampini     if (enable_bmark) {
927f86f7544SStefano Zampini       PetscInt lminl = PETSC_MAX_INT;
928f86f7544SStefano Zampini 
9299566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm,bmark,&label));
930f86f7544SStefano Zampini       if (label) {
931f86f7544SStefano Zampini         IS       vals;
932f86f7544SStefano Zampini         PetscInt ldef;
933f86f7544SStefano Zampini 
9349566063dSJacob Faibussowitsch         PetscCall(DMLabelGetDefaultValue(label,&ldef));
9359566063dSJacob Faibussowitsch         PetscCall(DMLabelGetValueIS(label,&vals));
9369566063dSJacob Faibussowitsch         PetscCall(ISGetMinMax(vals,&lminl,NULL));
9379566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&vals));
938f86f7544SStefano Zampini         lminl = PetscMin(ldef,lminl);
939f86f7544SStefano Zampini       }
9401c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm)));
941f86f7544SStefano Zampini       if (minl == PETSC_MAX_INT) minl = 1;
942f86f7544SStefano Zampini     }
943*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",bf));
9448135c375SStefano Zampini     for (p=fStart;p<fEnd;p++) {
94577eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
9468135c375SStefano Zampini         const PetscInt *support;
94777eacf09SStefano Zampini         PetscInt       supportSize,c,nc = 0;
9488135c375SStefano Zampini 
9499566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(dm,p,&supportSize));
9509566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(dm,p,&support));
9513e6c54aaSStefano Zampini         if (pown) {
95277eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
95377eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
95477eacf09SStefano Zampini               fcells[nc++] = support[c];
9558135c375SStefano Zampini             }
95677eacf09SStefano Zampini           }
95777eacf09SStefano Zampini         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
95877eacf09SStefano Zampini         for (c=0;c<nc;c++) {
959c3c203b2SStefano Zampini           const DMPolytopeType *faceTypes;
960c3c203b2SStefano Zampini           DMPolytopeType       cellType;
961c3c203b2SStefano Zampini           const PetscInt       *faceSizes,*cone;
962c3c203b2SStefano Zampini           PetscInt             vids[8],*faces,st,i,coneSize,cell,cl,nv,cid = -1,mid = -1;
9638135c375SStefano Zampini 
96477eacf09SStefano Zampini           cell = fcells[c];
9659566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCone(dm,cell,&cone));
9669566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm,cell,&coneSize));
967b135d7daSStefano Zampini           for (cl=0;cl<coneSize;cl++)
96877eacf09SStefano Zampini             if (cone[cl] == p)
9698135c375SStefano Zampini               break;
970b135d7daSStefano Zampini           if (cl == coneSize) continue;
9718135c375SStefano Zampini 
97277eacf09SStefano Zampini           /* face material id and type */
9739566063dSJacob Faibussowitsch           PetscCall(DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid));
974*63a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT,mid,cid));
97577eacf09SStefano Zampini           /* vertex ids */
9769566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(dm,cell,&cellType));
9779566063dSJacob Faibussowitsch           PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids));
9789566063dSJacob Faibussowitsch           PetscCall(DMPlexGetRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces));
979c3c203b2SStefano Zampini           st = 0;
980c3c203b2SStefano Zampini           for (i=0;i<cl;i++) st += faceSizes[i];
9819566063dSJacob Faibussowitsch           PetscCall(DMPlexInvertCell(faceTypes[cl],faces + st));
982c3c203b2SStefano Zampini           for (i=0;i<faceSizes[cl];i++) {
983*63a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT,faces[st+i]));
984b135d7daSStefano Zampini           }
9859566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
9869566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces));
9873e96840aSStefano Zampini           bf -= 1;
98877eacf09SStefano Zampini         }
9898135c375SStefano Zampini       }
9908135c375SStefano Zampini     }
991*63a3b9bcSJacob Faibussowitsch     PetscCheck(!bf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %" PetscInt_FMT,bf);
9929566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&bfaces));
9939566063dSJacob Faibussowitsch     PetscCall(PetscFree(fcells));
9948135c375SStefano Zampini   }
9958135c375SStefano Zampini 
9968135c375SStefano Zampini   /* mark owned vertices */
9973e6c54aaSStefano Zampini   vown = NULL;
9983e6c54aaSStefano Zampini   if (pown) {
9999566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(vEnd-vStart,&vown));
10008135c375SStefano Zampini     for (p=cStart;p<cEnd;p++) {
10018135c375SStefano Zampini       PetscInt i,closureSize,*closure = NULL;
10028135c375SStefano Zampini 
10033e6c54aaSStefano Zampini       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
10049566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure));
10058135c375SStefano Zampini       for (i=0;i<closureSize;i++) {
10068135c375SStefano Zampini         const PetscInt pp = closure[2*i];
10078135c375SStefano Zampini 
10088135c375SStefano Zampini         if (pp >= vStart && pp < vEnd) {
10099566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(vown,pp-vStart));
10108135c375SStefano Zampini         }
10118135c375SStefano Zampini       }
10129566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure));
10138135c375SStefano Zampini     }
10148135c375SStefano Zampini   }
10158135c375SStefano Zampini 
10168135c375SStefano Zampini   if (parentSection) {
10178135c375SStefano Zampini     PetscInt vp,gvp;
10188135c375SStefano Zampini 
10198135c375SStefano Zampini     for (vp=0,p=vStart;p<vEnd;p++) {
10208135c375SStefano Zampini       DMLabel  dlabel;
10218135c375SStefano Zampini       PetscInt parent,depth;
10228135c375SStefano Zampini 
10233e6c54aaSStefano Zampini       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
10249566063dSJacob Faibussowitsch       PetscCall(DMPlexGetDepthLabel(dm,&dlabel));
10259566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(dlabel,p,&depth));
10269566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTreeParent(dm,p,&parent,NULL));
10278135c375SStefano Zampini       if (parent != p) vp++;
10288135c375SStefano Zampini     }
10291c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
10308135c375SStefano Zampini     if (gvp) {
10318135c375SStefano Zampini       PetscInt  maxsupp;
10328135c375SStefano Zampini       PetscBool *skip = NULL;
10338135c375SStefano Zampini 
10349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n"));
1035*63a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",vp));
10369566063dSJacob Faibussowitsch       PetscCall(DMPlexGetMaxSizes(dm,NULL,&maxsupp));
10379566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(maxsupp,&skip));
10388135c375SStefano Zampini       for (p=vStart;p<vEnd;p++) {
10398135c375SStefano Zampini         DMLabel  dlabel;
10408135c375SStefano Zampini         PetscInt parent;
10418135c375SStefano Zampini 
10423e6c54aaSStefano Zampini         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
10439566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthLabel(dm,&dlabel));
10449566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTreeParent(dm,p,&parent,NULL));
10458135c375SStefano Zampini         if (parent != p) {
104696ca5757SLisandro Dalcin           PetscInt       vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
10473924b612SStefano Zampini           PetscInt       i,nv,ssize,n,numChildren,depth = -1;
10488135c375SStefano Zampini           const PetscInt *children;
10493924b612SStefano Zampini 
10509566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm,parent,&ssize));
10513924b612SStefano Zampini           switch (ssize) {
10528135c375SStefano Zampini             case 2: /* edge */
10538135c375SStefano Zampini               nv   = 0;
10549566063dSJacob Faibussowitsch               PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids));
1055*63a3b9bcSJacob Faibussowitsch               PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT,p-vStart));
10568135c375SStefano Zampini               for (i=0;i<nv;i++) {
1057*63a3b9bcSJacob Faibussowitsch                 PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT,vids[i]));
10588135c375SStefano Zampini               }
10599566063dSJacob Faibussowitsch               PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
10608135c375SStefano Zampini               vp--;
10618135c375SStefano Zampini               break;
10628135c375SStefano Zampini             case 4: /* face */
10639566063dSJacob Faibussowitsch               PetscCall(DMPlexGetTreeChildren(dm,parent,&numChildren,&children));
10648135c375SStefano Zampini               for (n=0;n<numChildren;n++) {
10659566063dSJacob Faibussowitsch                 PetscCall(DMLabelGetValue(dlabel,children[n],&depth));
10668135c375SStefano Zampini                 if (!depth) {
10678135c375SStefano Zampini                   const PetscInt *hvsupp,*hesupp,*cone;
10688135c375SStefano Zampini                   PetscInt       hvsuppSize,hesuppSize,coneSize;
1069451a39c7SStefano Zampini                   PetscInt       hv = children[n],he = -1,f;
10708135c375SStefano Zampini 
10719566063dSJacob Faibussowitsch                   PetscCall(PetscArrayzero(skip,maxsupp));
10729566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetSupportSize(dm,hv,&hvsuppSize));
10739566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetSupport(dm,hv,&hvsupp));
10748135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
10758135c375SStefano Zampini                     PetscInt ep;
10769566063dSJacob Faibussowitsch                     PetscCall(DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL));
10778135c375SStefano Zampini                     if (ep != hvsupp[i]) {
10788135c375SStefano Zampini                       he = hvsupp[i];
10798135c375SStefano Zampini                     } else {
10808135c375SStefano Zampini                       skip[i] = PETSC_TRUE;
10818135c375SStefano Zampini                     }
10828135c375SStefano Zampini                   }
1083*63a3b9bcSJacob Faibussowitsch                   PetscCheck(he != -1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %" PetscInt_FMT " support size %" PetscInt_FMT ": hanging edge not found",hv,hvsuppSize);
10849566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetCone(dm,he,&cone));
108596ca5757SLisandro Dalcin                   vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
10869566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetSupportSize(dm,he,&hesuppSize));
10879566063dSJacob Faibussowitsch                   PetscCall(DMPlexGetSupport(dm,he,&hesupp));
10888135c375SStefano Zampini                   for (f=0;f<hesuppSize;f++) {
10898135c375SStefano Zampini                     PetscInt j;
10908135c375SStefano Zampini 
10919566063dSJacob Faibussowitsch                     PetscCall(DMPlexGetCone(dm,hesupp[f],&cone));
10929566063dSJacob Faibussowitsch                     PetscCall(DMPlexGetConeSize(dm,hesupp[f],&coneSize));
10938135c375SStefano Zampini                     for (j=0;j<coneSize;j++) {
10948135c375SStefano Zampini                       PetscInt k;
10958135c375SStefano Zampini                       for (k=0;k<hvsuppSize;k++) {
10968135c375SStefano Zampini                         if (hvsupp[k] == cone[j]) {
10978135c375SStefano Zampini                           skip[k] = PETSC_TRUE;
10988135c375SStefano Zampini                           break;
10998135c375SStefano Zampini                         }
11008135c375SStefano Zampini                       }
11018135c375SStefano Zampini                     }
11028135c375SStefano Zampini                   }
11038135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
11048135c375SStefano Zampini                     if (!skip[i]) {
11059566063dSJacob Faibussowitsch                       PetscCall(DMPlexGetCone(dm,hvsupp[i],&cone));
110696ca5757SLisandro Dalcin                       vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
11078135c375SStefano Zampini                     }
11088135c375SStefano Zampini                   }
1109*63a3b9bcSJacob Faibussowitsch                   PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT,hv-vStart));
11108135c375SStefano Zampini                   for (i=0;i<2;i++) {
1111*63a3b9bcSJacob Faibussowitsch                     PetscCall(PetscViewerASCIIPrintf(viewer," %" PetscInt_FMT,vids[i]-vStart));
11128135c375SStefano Zampini                   }
11139566063dSJacob Faibussowitsch                   PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
11148135c375SStefano Zampini                   vp--;
11158135c375SStefano Zampini                 }
11168135c375SStefano Zampini               }
11178135c375SStefano Zampini               break;
11188135c375SStefano Zampini             default:
1119*63a3b9bcSJacob Faibussowitsch               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %" PetscInt_FMT,ssize);
11208135c375SStefano Zampini           }
11218135c375SStefano Zampini         }
11228135c375SStefano Zampini       }
11239566063dSJacob Faibussowitsch       PetscCall(PetscFree(skip));
11248135c375SStefano Zampini     }
1125*63a3b9bcSJacob Faibussowitsch     PetscCheck(!vp,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %" PetscInt_FMT " hanging vertices",vp);
11268135c375SStefano Zampini   }
11279566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&vown));
11288135c375SStefano Zampini 
11298135c375SStefano Zampini   /* vertices */
113077eacf09SStefano Zampini   if (hovec) { /* higher-order meshes */
113177eacf09SStefano Zampini     const char *fec;
11320286d493SLisandro Dalcin     PetscInt   i,n,s;
11339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
1134*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",vEnd-vStart));
11359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"nodes\n"));
11369566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)hovec,&fec));
11379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n"));
11389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",fec));
1139*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"VDim: %" PetscInt_FMT "\n",sdim));
11409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n")); /*Ordering::byVDIM*/
11410c2bc6bfSStefano Zampini     if (hoSection) {
11420c2bc6bfSStefano Zampini       DM cdm;
11430c2bc6bfSStefano Zampini 
11440c2bc6bfSStefano Zampini       PetscCall(VecGetDM(hovec,&cdm));
11450c2bc6bfSStefano Zampini       for (p=cStart;p<cEnd;p++) {
11460c2bc6bfSStefano Zampini         PetscScalar *vals = NULL;
11470c2bc6bfSStefano Zampini         PetscInt    csize;
11480c2bc6bfSStefano Zampini 
11490c2bc6bfSStefano Zampini         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
11500c2bc6bfSStefano Zampini         PetscCall(DMPlexVecGetClosure(cdm,hoSection,hovec,p,&csize,&vals));
1151*63a3b9bcSJacob Faibussowitsch         PetscCheck(csize%sdim == 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Size of closure %" PetscInt_FMT " incompatible with space dimension %" PetscInt_FMT,csize,sdim);
11520c2bc6bfSStefano Zampini         for (i=0;i<csize/sdim;i++) {
11530c2bc6bfSStefano Zampini           for (s=0;s<sdim;s++) {
11540c2bc6bfSStefano Zampini             PetscCall(PetscViewerASCIIPrintf(viewer,fmt,(double) PetscRealPart(vals[i*sdim+s])));
11550c2bc6bfSStefano Zampini           }
11560c2bc6bfSStefano Zampini           PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
11570c2bc6bfSStefano Zampini         }
11580c2bc6bfSStefano Zampini         PetscCall(DMPlexVecRestoreClosure(cdm,hoSection,hovec,p,&csize,&vals));
11590c2bc6bfSStefano Zampini       }
11600c2bc6bfSStefano Zampini     } else {
11619566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(hovec,&array));
11629566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(hovec,&n));
1163*63a3b9bcSJacob Faibussowitsch       PetscCheck(n%sdim == 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %" PetscInt_FMT " incompatible with space dimension %" PetscInt_FMT,n,sdim);
116477eacf09SStefano Zampini       for (i=0;i<n/sdim;i++) {
116577eacf09SStefano Zampini         for (s=0;s<sdim;s++) {
11669566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,fmt,(double) PetscRealPart(array[i*sdim+s])));
116777eacf09SStefano Zampini         }
11689566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
116977eacf09SStefano Zampini       }
11709566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(hovec,&array));
11710c2bc6bfSStefano Zampini     }
117277eacf09SStefano Zampini   } else {
11739566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coordinates,&nvert));
11749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
1175*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",nvert/sdim));
1176*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",sdim));
11779566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordinates,&array));
1178cc0d3ed7SStefano Zampini     for (p=0;p<nvert/sdim;p++) {
1179cc0d3ed7SStefano Zampini       PetscInt s;
1180cc0d3ed7SStefano Zampini       for (s=0;s<sdim;s++) {
11813e96840aSStefano Zampini         PetscReal v = PetscRealPart(array[p*sdim+s]);
11823e96840aSStefano Zampini 
11839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : (double) v));
11848135c375SStefano Zampini       }
11859566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
11868135c375SStefano Zampini     }
11879566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates,&array));
118877eacf09SStefano Zampini   }
11890c2bc6bfSStefano Zampini   PetscCall(PetscBTDestroy(&pown));
11909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&hovec));
11918135c375SStefano Zampini   PetscFunctionReturn(0);
11928135c375SStefano Zampini }
11938135c375SStefano Zampini 
11940286d493SLisandro Dalcin PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
11958135c375SStefano Zampini {
11968135c375SStefano Zampini   PetscFunctionBegin;
11979566063dSJacob Faibussowitsch   PetscCall(DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII));
11988135c375SStefano Zampini   PetscFunctionReturn(0);
11998135c375SStefano Zampini }
1200