xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 066ea43f7f75752f012be6cd06b6107ebe84cc6d)
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   PetscErrorCode ierr;
198135c375SStefano Zampini 
208135c375SStefano Zampini   PetscFunctionBegin;
218135c375SStefano Zampini   for (i=0;i<ctx->nf;i++) {
228135c375SStefano Zampini     ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr);
238135c375SStefano Zampini   }
248135c375SStefano Zampini   ierr = PetscFree(ctx->scctx);CHKERRQ(ierr);
2551f51421SSatish Balay   ierr = PetscFree(vctx);CHKERRQ(ierr);
268135c375SStefano Zampini   PetscFunctionReturn(0);
278135c375SStefano Zampini }
288135c375SStefano Zampini 
298135c375SStefano Zampini static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
308135c375SStefano Zampini {
318135c375SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
328135c375SStefano Zampini   PetscInt       f;
338135c375SStefano Zampini   PetscErrorCode ierr;
348135c375SStefano Zampini 
358135c375SStefano Zampini   PetscFunctionBegin;
368135c375SStefano Zampini   for (f=0;f<nf;f++) {
378135c375SStefano Zampini     ierr = VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
388135c375SStefano Zampini     ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
398135c375SStefano Zampini   }
408135c375SStefano Zampini   PetscFunctionReturn(0);
418135c375SStefano Zampini }
428135c375SStefano Zampini 
438135c375SStefano Zampini /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */
448135c375SStefano Zampini PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer)
458135c375SStefano Zampini {
468135c375SStefano Zampini   DM             dm = (DM)odm;
474cac2994SStefano Zampini   Vec            xlocal,xfield,*Ufield;
488135c375SStefano Zampini   PetscDS        ds;
498135c375SStefano Zampini   IS             globalNum,isfield;
508135c375SStefano Zampini   PetscBT        vown;
518135c375SStefano Zampini   char           **fieldname = NULL,**fec_type = NULL;
528135c375SStefano Zampini   const PetscInt *gNum;
53bb77a09fSStefano Zampini   PetscInt       *nlocal,*bs,*idxs,*dims;
548135c375SStefano Zampini   PetscInt       f,maxfields,nfields,c,totc,totdofs,Nv,cum,i;
55b135d7daSStefano Zampini   PetscInt       dim,cStart,cEnd,vStart,vEnd;
568135c375SStefano Zampini   GLVisViewerCtx *ctx;
578135c375SStefano Zampini   PetscSection   s;
588135c375SStefano Zampini   PetscErrorCode ierr;
598135c375SStefano Zampini 
608135c375SStefano Zampini   PetscFunctionBegin;
618135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
628135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
638135c375SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
64b135d7daSStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr);
65b135d7daSStefano Zampini   if (!globalNum) {
66b135d7daSStefano Zampini     ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr);
67b135d7daSStefano Zampini     ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr);
68b135d7daSStefano Zampini     ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr);
69b135d7daSStefano Zampini   }
708135c375SStefano Zampini   ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
718135c375SStefano Zampini   ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
728135c375SStefano Zampini   for (c = cStart, totc = 0; c < cEnd; c++) {
738135c375SStefano Zampini     if (gNum[c-cStart] >= 0) {
748135c375SStefano Zampini       PetscInt i,numPoints,*points = NULL;
758135c375SStefano Zampini 
768135c375SStefano Zampini       totc++;
778135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
788135c375SStefano Zampini       for (i=0;i<numPoints*2;i+= 2) {
798135c375SStefano Zampini         if ((points[i] >= vStart) && (points[i] < vEnd)) {
808135c375SStefano Zampini           ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr);
818135c375SStefano Zampini         }
828135c375SStefano Zampini       }
838135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
848135c375SStefano Zampini     }
858135c375SStefano Zampini   }
8677eacf09SStefano Zampini   for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++;
878135c375SStefano Zampini 
888135c375SStefano Zampini   ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr);
898135c375SStefano Zampini   ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr);
9092fd8e1eSJed Brown   ierr = DMGetLocalSection(dm,&s);CHKERRQ(ierr);
918135c375SStefano Zampini   ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr);
928135c375SStefano Zampini   for (f=0,maxfields=0;f<nfields;f++) {
938135c375SStefano Zampini     PetscInt bs;
948135c375SStefano Zampini 
958135c375SStefano Zampini     ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr);
968135c375SStefano Zampini     maxfields += bs;
978135c375SStefano Zampini   }
984cac2994SStefano Zampini   ierr = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);CHKERRQ(ierr);
998135c375SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
1008135c375SStefano Zampini   ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr);
1018135c375SStefano Zampini   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
1028135c375SStefano Zampini   if (ds) {
1038135c375SStefano Zampini     for (f=0;f<nfields;f++) {
1048135c375SStefano Zampini       const char* fname;
1058135c375SStefano Zampini       char        name[256];
1068135c375SStefano Zampini       PetscObject disc;
1078135c375SStefano Zampini       size_t      len;
1088135c375SStefano Zampini 
1098135c375SStefano Zampini       ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr);
1108135c375SStefano Zampini       ierr = PetscStrlen(fname,&len);CHKERRQ(ierr);
1118135c375SStefano Zampini       if (len) {
1128135c375SStefano Zampini         ierr = PetscStrcpy(name,fname);CHKERRQ(ierr);
1138135c375SStefano Zampini       } else {
114bfadf1d8SStefano Zampini         ierr = PetscSNPrintf(name,256,"Field%D",f);CHKERRQ(ierr);
1158135c375SStefano Zampini       }
1168135c375SStefano Zampini       ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
1178135c375SStefano Zampini       if (disc) {
1188135c375SStefano Zampini         PetscClassId id;
1198135c375SStefano Zampini         PetscInt     Nc;
1208135c375SStefano Zampini         char         fec[64];
1218135c375SStefano Zampini 
1228135c375SStefano Zampini         ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
1238135c375SStefano Zampini         if (id == PETSCFE_CLASSID) {
1248135c375SStefano Zampini           PetscFE            fem = (PetscFE)disc;
1258135c375SStefano Zampini           PetscDualSpace     sp;
1268135c375SStefano Zampini           PetscDualSpaceType spname;
1278135c375SStefano Zampini           PetscInt           order;
1288135c375SStefano Zampini           PetscBool          islag,continuous,H1 = PETSC_TRUE;
1298135c375SStefano Zampini 
1308135c375SStefano Zampini           ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr);
1318135c375SStefano Zampini           ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr);
1328135c375SStefano Zampini           ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr);
1338135c375SStefano Zampini           ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr);
1348135c375SStefano Zampini           if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
1358135c375SStefano Zampini           ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr);
1368135c375SStefano Zampini           ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
13728d58a37SPierre Jolivet           if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */
138bfadf1d8SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr);
1398135c375SStefano Zampini           } else {
14028d58a37SPierre Jolivet             if (!continuous && order) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %D",order);
1418135c375SStefano Zampini             H1   = PETSC_FALSE;
142bfadf1d8SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order);CHKERRQ(ierr);
1438135c375SStefano Zampini           }
1448135c375SStefano Zampini           ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr);
1458135c375SStefano Zampini           bs[ctx->nf]   = Nc;
146bb77a09fSStefano Zampini           dims[ctx->nf] = dim;
1478135c375SStefano Zampini           if (H1) {
1488135c375SStefano Zampini             nlocal[ctx->nf] = Nc * Nv;
1498135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
1508135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr);
1518135c375SStefano Zampini             for (i=0,cum=0;i<vEnd-vStart;i++) {
1528135c375SStefano Zampini               PetscInt j,off;
1538135c375SStefano Zampini 
1548135c375SStefano Zampini               if (PetscUnlikely(!PetscBTLookup(vown,i))) continue;
1558135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr);
1568135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1578135c375SStefano Zampini             }
1588135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
1598135c375SStefano Zampini           } else {
1608135c375SStefano Zampini             nlocal[ctx->nf] = Nc * totc;
1618135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
1628135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr);
1638135c375SStefano Zampini             for (i=0,cum=0;i<cEnd-cStart;i++) {
1648135c375SStefano Zampini               PetscInt j,off;
1658135c375SStefano Zampini 
1668135c375SStefano Zampini               if (PetscUnlikely(gNum[i] < 0)) continue;
1678135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
1688135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1698135c375SStefano Zampini             }
1708135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
1718135c375SStefano Zampini           }
1729448b7f1SJunchao Zhang           ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
1738135c375SStefano Zampini           ierr = VecDestroy(&xfield);CHKERRQ(ierr);
1748135c375SStefano Zampini           ierr = ISDestroy(&isfield);CHKERRQ(ierr);
1758135c375SStefano Zampini           ctx->nf++;
1768135c375SStefano Zampini         } else if (id == PETSCFV_CLASSID) {
1778135c375SStefano Zampini           PetscInt c;
1788135c375SStefano Zampini 
1798135c375SStefano Zampini           ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr);
180bfadf1d8SStefano Zampini           ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim);CHKERRQ(ierr);
1818135c375SStefano Zampini           for (c = 0; c < Nc; c++) {
1828135c375SStefano Zampini             char comp[256];
183bfadf1d8SStefano Zampini             ierr = PetscSNPrintf(comp,256,"%s-Comp%D",name,c);CHKERRQ(ierr);
1848135c375SStefano Zampini             ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr);
1858135c375SStefano Zampini             bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */
1868135c375SStefano Zampini             nlocal[ctx->nf] = totc;
187bb77a09fSStefano Zampini             dims[ctx->nf] = dim;
1888135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
1898135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr);
1908135c375SStefano Zampini             for (i=0,cum=0;i<cEnd-cStart;i++) {
1918135c375SStefano Zampini               PetscInt off;
1928135c375SStefano Zampini 
1938135c375SStefano Zampini               if (PetscUnlikely(gNum[i])<0) continue;
1948135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
1958135c375SStefano Zampini               idxs[cum++] = off + c;
1968135c375SStefano Zampini             }
1978135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
1989448b7f1SJunchao Zhang             ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
1998135c375SStefano Zampini             ierr = VecDestroy(&xfield);CHKERRQ(ierr);
2008135c375SStefano Zampini             ierr = ISDestroy(&isfield);CHKERRQ(ierr);
2018135c375SStefano Zampini             ctx->nf++;
2028135c375SStefano Zampini           }
2038135c375SStefano Zampini         } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f);
2048135c375SStefano Zampini       } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f);
2058135c375SStefano Zampini     }
2068135c375SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
2078135c375SStefano Zampini   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
2088135c375SStefano Zampini   ierr = VecDestroy(&xlocal);CHKERRQ(ierr);
2098135c375SStefano Zampini   ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
2108135c375SStefano Zampini 
2114cac2994SStefano Zampini   /* create work vectors */
2124cac2994SStefano Zampini   for (f=0;f<ctx->nf;f++) {
2134cac2994SStefano Zampini     ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr);
2144cac2994SStefano Zampini     ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr);
2154cac2994SStefano Zampini     ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr);
2164cac2994SStefano Zampini     ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr);
2174cac2994SStefano Zampini   }
2184cac2994SStefano Zampini 
2198135c375SStefano Zampini   /* customize the viewer */
2204cac2994SStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
2218135c375SStefano Zampini   for (f=0;f<ctx->nf;f++) {
2228135c375SStefano Zampini     ierr = PetscFree(fieldname[f]);CHKERRQ(ierr);
2238135c375SStefano Zampini     ierr = PetscFree(fec_type[f]);CHKERRQ(ierr);
2244cac2994SStefano Zampini     ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr);
2258135c375SStefano Zampini   }
2264cac2994SStefano Zampini   ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr);
2278135c375SStefano Zampini   PetscFunctionReturn(0);
2288135c375SStefano Zampini }
2298135c375SStefano Zampini 
230b135d7daSStefano Zampini typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid;
2318135c375SStefano Zampini 
2328135c375SStefano Zampini MFEM_cid mfem_table_cid[4][7]       = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_UNDEF},
2338135c375SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_UNDEF},
2348135c375SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF,MFEM_UNDEF},
235b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_PRISM,MFEM_CUBE } };
2368135c375SStefano Zampini 
237b135d7daSStefano 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},
238b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF},
239b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF},
240b135d7daSStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_CUBE } };
241044a5661SStefano Zampini 
242f86f7544SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid)
2438135c375SStefano Zampini {
2448135c375SStefano Zampini   DMLabel        dlabel;
245044a5661SStefano Zampini   PetscInt       depth,csize,pdepth,dim;
2468135c375SStefano Zampini   PetscErrorCode ierr;
2478135c375SStefano Zampini 
2488135c375SStefano Zampini   PetscFunctionBegin;
2498135c375SStefano Zampini   ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
250044a5661SStefano Zampini   ierr = DMLabelGetValue(dlabel,p,&pdepth);CHKERRQ(ierr);
2518135c375SStefano Zampini   ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr);
252044a5661SStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
253044a5661SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
2548135c375SStefano Zampini   if (label) {
2558135c375SStefano Zampini     ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr);
256f86f7544SStefano Zampini     *mid = *mid - minl + 1; /* MFEM does not like negative markers */
25777eacf09SStefano Zampini   } else *mid = 1;
258044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
259044a5661SStefano Zampini     if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim);
260044a5661SStefano Zampini     if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p);
261044a5661SStefano Zampini     if (depth != 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p);
262044a5661SStefano Zampini     *cid = mfem_table_cid_unint[dim][csize];
263044a5661SStefano Zampini   } else {
264044a5661SStefano Zampini     if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p);
265044a5661SStefano Zampini     if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p);
266044a5661SStefano Zampini     *cid = mfem_table_cid[pdepth][csize];
267044a5661SStefano Zampini   }
2688135c375SStefano Zampini   PetscFunctionReturn(0);
2698135c375SStefano Zampini }
2708135c375SStefano Zampini 
27196ca5757SLisandro Dalcin static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[])
2728135c375SStefano Zampini {
273cc0d3ed7SStefano Zampini   PetscInt       dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;
2748135c375SStefano Zampini   PetscErrorCode ierr;
2758135c375SStefano Zampini 
2768135c375SStefano Zampini   PetscFunctionBegin;
2778135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
278cc0d3ed7SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
279cc0d3ed7SStefano Zampini   sdim = dim;
280cc0d3ed7SStefano Zampini   if (csec) {
28184f354e3SLisandro Dalcin     PetscInt sStart,sEnd;
28284f354e3SLisandro Dalcin 
283cc0d3ed7SStefano Zampini     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
28484f354e3SLisandro Dalcin     ierr = PetscSectionGetChart(csec,&sStart,&sEnd);CHKERRQ(ierr);
285cc0d3ed7SStefano Zampini     ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr);
286cc0d3ed7SStefano Zampini     off  = off/sdim;
28784f354e3SLisandro Dalcin     if (p >= sStart && p < sEnd) {
288cc0d3ed7SStefano Zampini       ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
289cc0d3ed7SStefano Zampini     }
29084f354e3SLisandro Dalcin   }
291cc0d3ed7SStefano Zampini   if (!dof) {
2928135c375SStefano Zampini     ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
2938135c375SStefano Zampini     for (i=0,q=0;i<numPoints*2;i+= 2)
2948135c375SStefano Zampini       if ((points[i] >= vStart) && (points[i] < vEnd))
29596ca5757SLisandro Dalcin         vids[q++] = points[i]-vStart+off;
2968135c375SStefano Zampini     ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
297cc0d3ed7SStefano Zampini   } else {
298cc0d3ed7SStefano Zampini     ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr);
299cc0d3ed7SStefano Zampini     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
30096ca5757SLisandro Dalcin     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
301cc0d3ed7SStefano Zampini   }
3028135c375SStefano Zampini   *nv = q;
3038135c375SStefano Zampini   PetscFunctionReturn(0);
3048135c375SStefano Zampini }
3058135c375SStefano Zampini 
306*066ea43fSLisandro Dalcin static PetscErrorCode GLVisCreateFE(PetscFE femIn,char name[32],PetscFE *fem)
307*066ea43fSLisandro Dalcin {
308*066ea43fSLisandro Dalcin   DM              K;
309*066ea43fSLisandro Dalcin   PetscSpace      P;
310*066ea43fSLisandro Dalcin   PetscDualSpace  Q;
311*066ea43fSLisandro Dalcin   PetscQuadrature q,fq;
312*066ea43fSLisandro Dalcin   PetscInt        dim,deg,dof;
313*066ea43fSLisandro Dalcin   DMPolytopeType  ptype;
314*066ea43fSLisandro Dalcin   PetscBool       isSimplex,isTensor;
315*066ea43fSLisandro Dalcin   PetscBool       continuity = PETSC_FALSE;
316*066ea43fSLisandro Dalcin   PetscDTNodeType nodeType   = PETSCDTNODES_GAUSSJACOBI;
317*066ea43fSLisandro Dalcin   PetscBool       endpoint   = PETSC_TRUE;
318*066ea43fSLisandro Dalcin   MPI_Comm        comm;
319*066ea43fSLisandro Dalcin   PetscErrorCode  ierr;
320*066ea43fSLisandro Dalcin 
321*066ea43fSLisandro Dalcin   PetscFunctionBegin;
322*066ea43fSLisandro Dalcin   comm = PetscObjectComm((PetscObject)femIn);
323*066ea43fSLisandro Dalcin   ierr = PetscFEGetBasisSpace(femIn,&P);CHKERRQ(ierr);
324*066ea43fSLisandro Dalcin   ierr = PetscFEGetDualSpace(femIn,&Q);CHKERRQ(ierr);
325*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceGetDM(Q,&K);CHKERRQ(ierr);
326*066ea43fSLisandro Dalcin   ierr = DMGetDimension(K,&dim);CHKERRQ(ierr);
327*066ea43fSLisandro Dalcin   ierr = PetscSpaceGetDegree(P,&deg,NULL);CHKERRQ(ierr);
328*066ea43fSLisandro Dalcin   ierr = PetscSpaceGetNumComponents(P,&dof);CHKERRQ(ierr);
329*066ea43fSLisandro Dalcin   ierr = DMPlexGetCellType(K,0,&ptype);CHKERRQ(ierr);
330*066ea43fSLisandro Dalcin   switch (ptype) {
331*066ea43fSLisandro Dalcin   case DM_POLYTOPE_QUADRILATERAL:
332*066ea43fSLisandro Dalcin   case DM_POLYTOPE_HEXAHEDRON:
333*066ea43fSLisandro Dalcin     isSimplex = PETSC_FALSE; break;
334*066ea43fSLisandro Dalcin   default:
335*066ea43fSLisandro Dalcin     isSimplex = PETSC_TRUE; break;
336*066ea43fSLisandro Dalcin   }
337*066ea43fSLisandro Dalcin   isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
338*066ea43fSLisandro Dalcin   /* Create space */
339*066ea43fSLisandro Dalcin   ierr = PetscSpaceCreate(comm,&P);CHKERRQ(ierr);
340*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
341*066ea43fSLisandro Dalcin   ierr = PetscSpacePolynomialSetTensor(P,isTensor);CHKERRQ(ierr);
342*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetNumComponents(P,dof);CHKERRQ(ierr);
343*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetNumVariables(P,dim);CHKERRQ(ierr);
344*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetDegree(P,deg,PETSC_DETERMINE);CHKERRQ(ierr);
345*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
346*066ea43fSLisandro Dalcin   /* Create dual space */
347*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceCreate(comm,&Q);CHKERRQ(ierr);
348*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
349*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceLagrangeSetTensor(Q,isTensor);CHKERRQ(ierr);
350*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceLagrangeSetContinuity(Q,continuity);CHKERRQ(ierr);
351*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceLagrangeSetNodeType(Q,nodeType,endpoint,0);CHKERRQ(ierr);
352*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetNumComponents(Q,dof);CHKERRQ(ierr);
353*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetOrder(Q,deg);CHKERRQ(ierr);
354*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceCreateReferenceCell(Q,dim,isSimplex,&K);CHKERRQ(ierr);
355*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetDM(Q,K);CHKERRQ(ierr);
356*066ea43fSLisandro Dalcin   ierr = DMDestroy(&K);CHKERRQ(ierr);
357*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
358*066ea43fSLisandro Dalcin   /* Create quadrature */
359*066ea43fSLisandro Dalcin   if (isSimplex) {
360*066ea43fSLisandro Dalcin     ierr = PetscDTStroudConicalQuadrature(dim,  1,deg+1,-1,+1,&q);CHKERRQ(ierr);
361*066ea43fSLisandro Dalcin     ierr = PetscDTStroudConicalQuadrature(dim-1,1,deg+1,-1,+1,&fq);CHKERRQ(ierr);
362*066ea43fSLisandro Dalcin   } else {
363*066ea43fSLisandro Dalcin     ierr = PetscDTGaussTensorQuadrature(dim,  1,deg+1,-1,+1,&q);CHKERRQ(ierr);
364*066ea43fSLisandro Dalcin     ierr = PetscDTGaussTensorQuadrature(dim-1,1,deg+1,-1,+1,&fq);CHKERRQ(ierr);
365*066ea43fSLisandro Dalcin   }
366*066ea43fSLisandro Dalcin   /* Create finite element */
367*066ea43fSLisandro Dalcin   ierr = PetscFECreate(comm,fem);CHKERRQ(ierr);
368*066ea43fSLisandro Dalcin   ierr = PetscSNPrintf(name,32,"L2_T1_%DD_P%D",dim,deg);CHKERRQ(ierr);
369*066ea43fSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject)*fem,name);CHKERRQ(ierr);
370*066ea43fSLisandro Dalcin   ierr = PetscFESetType(*fem,PETSCFEBASIC);CHKERRQ(ierr);
371*066ea43fSLisandro Dalcin   ierr = PetscFESetNumComponents(*fem,dof);CHKERRQ(ierr);
372*066ea43fSLisandro Dalcin   ierr = PetscFESetBasisSpace(*fem,P);CHKERRQ(ierr);
373*066ea43fSLisandro Dalcin   ierr = PetscFESetDualSpace(*fem,Q);CHKERRQ(ierr);
374*066ea43fSLisandro Dalcin   ierr = PetscFESetQuadrature(*fem,q);CHKERRQ(ierr);
375*066ea43fSLisandro Dalcin   ierr = PetscFESetFaceQuadrature(*fem,fq);CHKERRQ(ierr);
376*066ea43fSLisandro Dalcin   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
377*066ea43fSLisandro Dalcin   /* Cleanup */
378*066ea43fSLisandro Dalcin   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
379*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
380*066ea43fSLisandro Dalcin   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
381*066ea43fSLisandro Dalcin   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
382*066ea43fSLisandro Dalcin   PetscFunctionReturn(0);
383*066ea43fSLisandro Dalcin }
384*066ea43fSLisandro Dalcin 
38577eacf09SStefano Zampini /*
38677eacf09SStefano Zampini    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
38777eacf09SStefano Zampini    Higher order meshes are also supported
38877eacf09SStefano Zampini */
3898135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
3908135c375SStefano Zampini {
3918135c375SStefano Zampini   DMLabel              label;
3928135c375SStefano Zampini   PetscSection         coordSection,parentSection;
39377eacf09SStefano Zampini   Vec                  coordinates,hovec;
3948135c375SStefano Zampini   const PetscScalar    *array;
395f86f7544SStefano Zampini   PetscInt             bf,p,sdim,dim,depth,novl,minl;
396412e9a14SMatthew G. Knepley   PetscInt             cStart,cEnd,vStart,vEnd,nvert;
3973924b612SStefano Zampini   PetscMPIInt          size;
3983e6c54aaSStefano Zampini   PetscBool            localized,isascii;
39928d58a37SPierre Jolivet   PetscBool            enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE;
4003e6c54aaSStefano Zampini   PetscBT              pown,vown;
4018135c375SStefano Zampini   PetscErrorCode       ierr;
4028135c375SStefano Zampini   PetscContainer       glvis_container;
403044a5661SStefano Zampini   PetscBool            cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
404f86f7544SStefano Zampini   PetscBool            enable_emark,enable_bmark;
40577eacf09SStefano Zampini   const char           *fmt;
4067bf4dd16SStefano Zampini   char                 emark[64] = "",bmark[64] = "";
4078135c375SStefano Zampini 
4088135c375SStefano Zampini   PetscFunctionBegin;
4098135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4108135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
4118135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
4128135c375SStefano Zampini   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
4133924b612SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
4143924b612SStefano Zampini   if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
4158135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
4168135c375SStefano Zampini 
4178135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
4188135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
4198135c375SStefano Zampini   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
4208135c375SStefano Zampini   {
4218135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
4228135c375SStefano Zampini     ierr    = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
4238135c375SStefano Zampini     enabled = glvis_info->enabled;
42477eacf09SStefano Zampini     fmt     = glvis_info->fmt;
4258135c375SStefano Zampini   }
42621414b21SStefano Zampini 
42721414b21SStefano Zampini   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
42821414b21SStefano Zampini      DMPlex does not currently support HO meshes, so there's no API for this */
42921414b21SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr);
430*066ea43fSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)hovec);CHKERRQ(ierr);
431*066ea43fSLisandro Dalcin   if (!hovec) {
432*066ea43fSLisandro Dalcin     DM           cdm;
433*066ea43fSLisandro Dalcin     PetscFE      disc;
434*066ea43fSLisandro Dalcin     PetscClassId classid;
435*066ea43fSLisandro Dalcin 
436*066ea43fSLisandro Dalcin     ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr);
437*066ea43fSLisandro Dalcin     ierr = DMGetField(cdm,0,NULL,(PetscObject*)&disc);CHKERRQ(ierr);
438*066ea43fSLisandro Dalcin     ierr = PetscObjectGetClassId((PetscObject)disc,&classid);CHKERRQ(ierr);
439*066ea43fSLisandro Dalcin     if (classid == PETSCFE_CLASSID) {
440*066ea43fSLisandro Dalcin       DM      hocdm;
441*066ea43fSLisandro Dalcin       PetscFE hodisc;
442*066ea43fSLisandro Dalcin       Vec     vec;
443*066ea43fSLisandro Dalcin       Mat     mat;
444*066ea43fSLisandro Dalcin       char    name[32],fec_type[64];
445*066ea43fSLisandro Dalcin 
446*066ea43fSLisandro Dalcin       ierr = GLVisCreateFE(disc,name,&hodisc);CHKERRQ(ierr);
447*066ea43fSLisandro Dalcin       ierr = DMClone(cdm,&hocdm);CHKERRQ(ierr);
448*066ea43fSLisandro Dalcin       ierr = DMSetField(hocdm,0,NULL,(PetscObject)hodisc);CHKERRQ(ierr);
449*066ea43fSLisandro Dalcin       ierr = PetscFEDestroy(&hodisc);CHKERRQ(ierr);
450*066ea43fSLisandro Dalcin       ierr = DMCreateDS(hocdm);CHKERRQ(ierr);
451*066ea43fSLisandro Dalcin 
452*066ea43fSLisandro Dalcin       ierr = DMGetCoordinates(dm,&vec);CHKERRQ(ierr);
453*066ea43fSLisandro Dalcin       ierr = DMCreateGlobalVector(hocdm,&hovec);CHKERRQ(ierr);
454*066ea43fSLisandro Dalcin       ierr = DMCreateInterpolation(cdm,hocdm,&mat,NULL);CHKERRQ(ierr);
455*066ea43fSLisandro Dalcin       ierr = MatInterpolate(mat,vec,hovec);CHKERRQ(ierr);
456*066ea43fSLisandro Dalcin       ierr = MatDestroy(&mat);CHKERRQ(ierr);
457*066ea43fSLisandro Dalcin       ierr = DMDestroy(&hocdm);CHKERRQ(ierr);
458*066ea43fSLisandro Dalcin 
459*066ea43fSLisandro Dalcin       ierr = PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name);CHKERRQ(ierr);
460*066ea43fSLisandro Dalcin       ierr = PetscObjectSetName((PetscObject)hovec,fec_type);CHKERRQ(ierr);
461*066ea43fSLisandro Dalcin     }
462*066ea43fSLisandro Dalcin   }
46321414b21SStefano Zampini 
464c3c203b2SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
465c3c203b2SStefano Zampini   ierr = DMPlexGetGhostCellStratum(dm,&p,NULL);CHKERRQ(ierr);
466c3c203b2SStefano Zampini   if (p >= 0) cEnd = p;
4678135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
46877eacf09SStefano Zampini   ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr);
469cc0d3ed7SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
47021414b21SStefano Zampini   if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
471cc0d3ed7SStefano Zampini   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
47277eacf09SStefano Zampini   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
47377eacf09SStefano Zampini   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
47421414b21SStefano Zampini   if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");
4758135c375SStefano Zampini 
4768135c375SStefano Zampini   /*
4778135c375SStefano Zampini      a couple of sections of the mesh specification are disabled
4783e96840aSStefano 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)
47977eacf09SStefano Zampini        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
4803e6c54aaSStefano Zampini                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
4818135c375SStefano Zampini   */
4823e96840aSStefano Zampini   enable_boundary = PETSC_FALSE;
4838135c375SStefano Zampini   enable_ncmesh   = PETSC_FALSE;
4843e6c54aaSStefano Zampini   enable_mfem     = PETSC_FALSE;
485f86f7544SStefano Zampini   enable_emark    = PETSC_FALSE;
486f86f7544SStefano Zampini   enable_bmark    = PETSC_FALSE;
4877bf4dd16SStefano Zampini   /* I'm tired of problems with negative values in the markers, disable them */
4888135c375SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
4893e6c54aaSStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr);
4903e6c54aaSStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL);CHKERRQ(ierr);
4913e6c54aaSStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL);CHKERRQ(ierr);
49228d58a37SPierre Jolivet   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL);CHKERRQ(ierr);
493f86f7544SStefano Zampini   ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);CHKERRQ(ierr);
494f86f7544SStefano Zampini   ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);CHKERRQ(ierr);
4958135c375SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
496f86f7544SStefano Zampini   if (enable_bmark) enable_boundary = PETSC_TRUE;
497f86f7544SStefano Zampini 
4983924b612SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
4993924b612SStefano Zampini   if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
5007e1aca4eSStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
5017e1aca4eSStefano Zampini   if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
5027e1aca4eSStefano Zampini                                                              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
5037e1aca4eSStefano Zampini   if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
5047e1aca4eSStefano Zampini                                                            "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
505044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
506044a5661SStefano Zampini     if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth);
507044a5661SStefano Zampini     cellvertex = PETSC_TRUE;
508044a5661SStefano Zampini   }
5098135c375SStefano Zampini 
5108135c375SStefano Zampini   /* Identify possible cells in the overlap */
5118135c375SStefano Zampini   novl = 0;
5128135c375SStefano Zampini   pown = NULL;
5133924b612SStefano Zampini   if (size > 1) {
5143e6c54aaSStefano Zampini     IS             globalNum = NULL;
5153e6c54aaSStefano Zampini     const PetscInt *gNum;
5163e6c54aaSStefano Zampini     PetscBool      ovl  = PETSC_FALSE;
5173e6c54aaSStefano Zampini 
518b135d7daSStefano Zampini     ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr);
519b135d7daSStefano Zampini     if (!globalNum) {
52028d58a37SPierre Jolivet       if (view_ovl) {
52128d58a37SPierre Jolivet         ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum);CHKERRQ(ierr);
52228d58a37SPierre Jolivet       } else {
523b135d7daSStefano Zampini         ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr);
52428d58a37SPierre Jolivet       }
525b135d7daSStefano Zampini       ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr);
526b135d7daSStefano Zampini       ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr);
527b135d7daSStefano Zampini     }
5288135c375SStefano Zampini     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
5298135c375SStefano Zampini     for (p=cStart; p<cEnd; p++) {
5308135c375SStefano Zampini       if (gNum[p-cStart] < 0) {
5318135c375SStefano Zampini         ovl = PETSC_TRUE;
5328135c375SStefano Zampini         novl++;
5338135c375SStefano Zampini       }
5348135c375SStefano Zampini     }
5358135c375SStefano Zampini     if (ovl) {
5368135c375SStefano Zampini       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
5378135c375SStefano Zampini          TODO: garbage collector? attach pown to dm?  */
5383e6c54aaSStefano Zampini       ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr);
5393e6c54aaSStefano Zampini       for (p=cStart; p<cEnd; p++) {
5403e6c54aaSStefano Zampini         if (gNum[p-cStart] < 0) continue;
5413e6c54aaSStefano Zampini         else {
5423e6c54aaSStefano Zampini           ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
5438135c375SStefano Zampini         }
5448135c375SStefano Zampini       }
5453e6c54aaSStefano Zampini     }
5463e6c54aaSStefano Zampini     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
5473e6c54aaSStefano Zampini   }
5488135c375SStefano Zampini 
5493e6c54aaSStefano Zampini   /* return if this process is disabled */
5508135c375SStefano Zampini   if (!enabled) {
5518135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
5528135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
5538135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
5548135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
5558135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
5568135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
5578135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
5588135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
5598135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
5608135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
5618135c375SStefano Zampini     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
562*066ea43fSLisandro Dalcin     ierr = VecDestroy(&hovec);CHKERRQ(ierr);
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 
576bfadf1d8SStefano Zampini       ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr);
5773e6c54aaSStefano Zampini       if (cEnd-cStart) {
5783e6c54aaSStefano Zampini         PetscInt fpc;
5793e6c54aaSStefano Zampini 
5803e6c54aaSStefano Zampini         ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
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:
5973e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
5983e6c54aaSStefano Zampini                 break;
5993e6c54aaSStefano Zampini             }
6003e6c54aaSStefano Zampini             break;
6013e6c54aaSStefano Zampini           case 3:
6023e6c54aaSStefano Zampini             switch (fpc) {
603044a5661SStefano Zampini               case 4: /* TODO: still need to understand L2 ordering for tets */
6043e6c54aaSStefano Zampini                 vpc = 4;
605044a5661SStefano Zampini                 dof = tetv;
606044a5661SStefano Zampini                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
6073e6c54aaSStefano Zampini                 break;
6083e6c54aaSStefano Zampini               case 6:
609044a5661SStefano Zampini                 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc);
610044a5661SStefano Zampini                 vpc = 8;
611044a5661SStefano Zampini                 dof = hexv;
612044a5661SStefano Zampini                 break;
613044a5661SStefano Zampini               case 8:
614044a5661SStefano Zampini                 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
6153e6c54aaSStefano Zampini                 vpc = 8;
6163e6c54aaSStefano Zampini                 dof = hexv;
6173e6c54aaSStefano Zampini                 break;
6183e6c54aaSStefano Zampini               default:
6193e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
6203e6c54aaSStefano Zampini                 break;
6213e6c54aaSStefano Zampini             }
6223e6c54aaSStefano Zampini             break;
6233e6c54aaSStefano Zampini           default:
6243e6c54aaSStefano Zampini             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
6253e6c54aaSStefano Zampini             break;
6263e6c54aaSStefano Zampini         }
62796ca5757SLisandro Dalcin         ierr = DMPlexReorderCell(dm,cStart,vids);CHKERRQ(ierr);
6283e6c54aaSStefano Zampini       }
629044a5661SStefano Zampini       if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
6303e6c54aaSStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr);
6313e6c54aaSStefano Zampini       ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr);
6323e6c54aaSStefano Zampini       ierr = VecGetArray(hovec,&array);CHKERRQ(ierr);
6333e6c54aaSStefano Zampini       ptr  = array;
6343e6c54aaSStefano Zampini       for (p=cStart;p<cEnd;p++) {
6353e6c54aaSStefano Zampini         PetscInt    csize,v,d;
6363e6c54aaSStefano Zampini         PetscScalar *vals = NULL;
6373e6c54aaSStefano Zampini 
6383e6c54aaSStefano Zampini         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
6393e6c54aaSStefano Zampini         ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
6403e6c54aaSStefano Zampini         if (csize != vpc*sdim && csize != vpc*sdim*2) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported closure size %D (vpc %D, sdim %D)",csize,vpc,sdim);
6413e6c54aaSStefano Zampini         for (v=0;v<vpc;v++) {
6423e6c54aaSStefano Zampini           for (d=0;d<sdim;d++) {
6433e6c54aaSStefano Zampini             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
6443e6c54aaSStefano Zampini           }
6453e6c54aaSStefano Zampini         }
6463e6c54aaSStefano Zampini         ptr += vpc*sdim;
6473e6c54aaSStefano Zampini         ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
6483e6c54aaSStefano Zampini       }
6493e6c54aaSStefano Zampini       ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr);
6503e6c54aaSStefano Zampini     }
6513e6c54aaSStefano Zampini   }
6523e96840aSStefano Zampini   /* if we have high-order coordinates in 3D, we need to specify the boundary */
6533e96840aSStefano Zampini   if (hovec && dim == 3) enable_boundary = PETSC_TRUE;
6543e6c54aaSStefano Zampini 
6558135c375SStefano Zampini   /* header */
6568135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
6578135c375SStefano Zampini 
6588135c375SStefano Zampini   /* topological dimension */
6598135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
6608135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
6618135c375SStefano Zampini 
6628135c375SStefano Zampini   /* elements */
663f86f7544SStefano Zampini   minl = 1;
664f86f7544SStefano Zampini   label = NULL;
665f86f7544SStefano Zampini   if (enable_emark) {
666f86f7544SStefano Zampini     PetscInt lminl = PETSC_MAX_INT;
667f86f7544SStefano Zampini 
6687bf4dd16SStefano Zampini     ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr);
669f86f7544SStefano Zampini     if (label) {
670f86f7544SStefano Zampini       IS       vals;
671f86f7544SStefano Zampini       PetscInt ldef;
672f86f7544SStefano Zampini 
673f86f7544SStefano Zampini       ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr);
674f86f7544SStefano Zampini       ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr);
675f86f7544SStefano Zampini       ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr);
676f86f7544SStefano Zampini       ierr = ISDestroy(&vals);CHKERRQ(ierr);
677f86f7544SStefano Zampini       lminl = PetscMin(ldef,lminl);
678f86f7544SStefano Zampini     }
679f86f7544SStefano Zampini     ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
680f86f7544SStefano Zampini     if (minl == PETSC_MAX_INT) minl = 1;
681f86f7544SStefano Zampini   }
6828135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
6838135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
6848135c375SStefano Zampini   for (p=cStart;p<cEnd;p++) {
68596ca5757SLisandro Dalcin     PetscInt       vids[8];
68611a4995dSStefano Zampini     PetscInt       i,nv = 0,cid = -1,mid = 1;
6878135c375SStefano Zampini 
6883e6c54aaSStefano Zampini     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
689f86f7544SStefano Zampini     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr);
69077eacf09SStefano Zampini     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
69196ca5757SLisandro Dalcin     ierr = DMPlexReorderCell(dm,p,vids);CHKERRQ(ierr);
6928135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
6938135c375SStefano Zampini     for (i=0;i<nv;i++) {
69496ca5757SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]);CHKERRQ(ierr);
6958135c375SStefano Zampini     }
6968135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6978135c375SStefano Zampini   }
6988135c375SStefano Zampini 
699cc0d3ed7SStefano Zampini   /* boundary */
700cc0d3ed7SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
701cc0d3ed7SStefano Zampini   if (!enable_boundary) {
702cc0d3ed7SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
703cc0d3ed7SStefano Zampini   } else {
70477eacf09SStefano Zampini     DMLabel  perLabel;
70577eacf09SStefano Zampini     PetscBT  bfaces;
706b135d7daSStefano Zampini     PetscInt fStart,fEnd,*fcells;
707cc0d3ed7SStefano Zampini 
7088135c375SStefano Zampini     ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
70977eacf09SStefano Zampini     ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr);
71077eacf09SStefano Zampini     ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr);
71177eacf09SStefano Zampini     ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr);
71277eacf09SStefano Zampini     ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
71377eacf09SStefano Zampini     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
71477eacf09SStefano Zampini       ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr);
71577eacf09SStefano Zampini       ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
71677eacf09SStefano Zampini       ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr);
71777eacf09SStefano Zampini       for (p=cStart;p<cEnd;p++) {
718c3c203b2SStefano Zampini         DMPolytopeType cellType;
719c3c203b2SStefano Zampini         PetscInt       dof;
720b135d7daSStefano Zampini 
721c3c203b2SStefano Zampini         ierr = DMPlexGetCellType(dm,p,&cellType);CHKERRQ(ierr);
72277eacf09SStefano Zampini         ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
72377eacf09SStefano Zampini         if (dof) {
724c3c203b2SStefano Zampini           PetscInt    uvpc, v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
72577eacf09SStefano Zampini           PetscScalar *vals = NULL;
726c3c203b2SStefano Zampini 
727c3c203b2SStefano Zampini           uvpc = DMPolytopeTypeGetNumVertices(cellType);
72877eacf09SStefano Zampini           if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
72977eacf09SStefano Zampini           ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
73077eacf09SStefano Zampini           ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
73177eacf09SStefano Zampini           for (v=0;v<cellClosureSize;v++)
73277eacf09SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
73377eacf09SStefano Zampini               vidxs = cellClosure + 2*v;
73477eacf09SStefano Zampini               break;
73577eacf09SStefano Zampini             }
73677eacf09SStefano Zampini           if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
737b135d7daSStefano Zampini           for (v=0;v<uvpc;v++) {
73877eacf09SStefano Zampini             PetscInt s;
739044a5661SStefano Zampini 
74077eacf09SStefano Zampini             for (s=0;s<sdim;s++) {
741b135d7daSStefano Zampini               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) {
74277eacf09SStefano Zampini                 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr);
74377eacf09SStefano Zampini               }
74477eacf09SStefano Zampini             }
74577eacf09SStefano Zampini           }
74677eacf09SStefano Zampini           ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
74777eacf09SStefano Zampini           ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
74877eacf09SStefano Zampini         }
74977eacf09SStefano Zampini       }
75077eacf09SStefano Zampini       if (dim > 1) {
751b135d7daSStefano Zampini         PetscInt eEnd,eStart;
752044a5661SStefano Zampini 
75377eacf09SStefano Zampini         ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr);
75477eacf09SStefano Zampini         for (p=eStart;p<eEnd;p++) {
75577eacf09SStefano Zampini           const PetscInt *cone;
75677eacf09SStefano Zampini           PetscInt       coneSize,i;
75777eacf09SStefano Zampini           PetscBool      ispe = PETSC_TRUE;
75877eacf09SStefano Zampini 
75977eacf09SStefano Zampini           ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
76077eacf09SStefano Zampini           ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
76177eacf09SStefano Zampini           for (i=0;i<coneSize;i++) {
76277eacf09SStefano Zampini             PetscInt v;
76377eacf09SStefano Zampini 
76477eacf09SStefano Zampini             ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
76577eacf09SStefano Zampini             ispe = (PetscBool)(ispe && (v==2));
76677eacf09SStefano Zampini           }
76777eacf09SStefano Zampini           if (ispe && coneSize) {
7683e96840aSStefano Zampini             PetscInt       ch, numChildren;
7693e96840aSStefano Zampini             const PetscInt *children;
7703e96840aSStefano Zampini 
77177eacf09SStefano Zampini             ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
7723e96840aSStefano Zampini             ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr);
7733e96840aSStefano Zampini             for (ch = 0; ch < numChildren; ch++) {
7743e96840aSStefano Zampini               ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr);
7753e96840aSStefano Zampini             }
77677eacf09SStefano Zampini           }
77777eacf09SStefano Zampini         }
77877eacf09SStefano Zampini         if (dim > 2) {
77977eacf09SStefano Zampini           for (p=fStart;p<fEnd;p++) {
78077eacf09SStefano Zampini             const PetscInt *cone;
78177eacf09SStefano Zampini             PetscInt       coneSize,i;
78277eacf09SStefano Zampini             PetscBool      ispe = PETSC_TRUE;
78377eacf09SStefano Zampini 
78477eacf09SStefano Zampini             ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
78577eacf09SStefano Zampini             ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
78677eacf09SStefano Zampini             for (i=0;i<coneSize;i++) {
78777eacf09SStefano Zampini               PetscInt v;
78877eacf09SStefano Zampini 
78977eacf09SStefano Zampini               ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
79077eacf09SStefano Zampini               ispe = (PetscBool)(ispe && (v==2));
79177eacf09SStefano Zampini             }
79277eacf09SStefano Zampini             if (ispe && coneSize) {
7933e96840aSStefano Zampini               PetscInt       ch, numChildren;
7943e96840aSStefano Zampini               const PetscInt *children;
7953e96840aSStefano Zampini 
79677eacf09SStefano Zampini               ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
7973e96840aSStefano Zampini               ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr);
7983e96840aSStefano Zampini               for (ch = 0; ch < numChildren; ch++) {
7993e96840aSStefano Zampini                 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr);
8003e96840aSStefano Zampini               }
80177eacf09SStefano Zampini             }
80277eacf09SStefano Zampini           }
80377eacf09SStefano Zampini         }
80477eacf09SStefano Zampini       }
80577eacf09SStefano Zampini     }
80677eacf09SStefano Zampini     for (p=fStart;p<fEnd;p++) {
80777eacf09SStefano Zampini       const PetscInt *support;
8088135c375SStefano Zampini       PetscInt       supportSize;
80977eacf09SStefano Zampini       PetscBool      isbf = PETSC_FALSE;
8108135c375SStefano Zampini 
8118135c375SStefano Zampini       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
8123e6c54aaSStefano Zampini       if (pown) {
8138135c375SStefano Zampini         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
81477eacf09SStefano Zampini         PetscInt  i;
81577eacf09SStefano Zampini 
81677eacf09SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
81777eacf09SStefano Zampini         for (i=0;i<supportSize;i++) {
81877eacf09SStefano Zampini           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
81977eacf09SStefano Zampini           else has_ghost = PETSC_TRUE;
82077eacf09SStefano Zampini         }
82177eacf09SStefano Zampini         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
82277eacf09SStefano Zampini       } else {
82377eacf09SStefano Zampini         isbf = (PetscBool)(supportSize == 1);
82477eacf09SStefano Zampini       }
82577eacf09SStefano Zampini       if (!isbf && perLabel) {
82677eacf09SStefano Zampini         const PetscInt *cone;
82777eacf09SStefano Zampini         PetscInt       coneSize,i;
82877eacf09SStefano Zampini 
82977eacf09SStefano Zampini         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
83077eacf09SStefano Zampini         ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
83177eacf09SStefano Zampini         isbf = PETSC_TRUE;
83277eacf09SStefano Zampini         for (i=0;i<coneSize;i++) {
83377eacf09SStefano Zampini           PetscInt v,d;
83477eacf09SStefano Zampini 
83577eacf09SStefano Zampini           ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
83677eacf09SStefano Zampini           ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr);
83777eacf09SStefano Zampini           isbf = (PetscBool)(isbf && v != d);
83877eacf09SStefano Zampini         }
83977eacf09SStefano Zampini       }
84077eacf09SStefano Zampini       if (isbf) {
84177eacf09SStefano Zampini         ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr);
84277eacf09SStefano Zampini       }
84377eacf09SStefano Zampini     }
84477eacf09SStefano Zampini     /* count boundary faces */
84577eacf09SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
84677eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
84777eacf09SStefano Zampini         const PetscInt *support;
84877eacf09SStefano Zampini         PetscInt       supportSize,c;
8498135c375SStefano Zampini 
8508135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
8518135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
85277eacf09SStefano Zampini         for (c=0;c<supportSize;c++) {
8533e96840aSStefano Zampini           const    PetscInt *cone;
854b135d7daSStefano Zampini           PetscInt cell,cl,coneSize;
8553e96840aSStefano Zampini 
8563e96840aSStefano Zampini           cell = support[c];
8573e96840aSStefano Zampini           if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue;
8583e96840aSStefano Zampini           ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr);
859b135d7daSStefano Zampini           ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
860b135d7daSStefano Zampini           for (cl=0;cl<coneSize;cl++) {
8613e96840aSStefano Zampini             if (cone[cl] == p) {
8623e96840aSStefano Zampini               bf += 1;
8633e96840aSStefano Zampini               break;
8648135c375SStefano Zampini             }
86577eacf09SStefano Zampini           }
8663e96840aSStefano Zampini         }
8678135c375SStefano Zampini       }
8688135c375SStefano Zampini     }
869f86f7544SStefano Zampini     minl = 1;
870f86f7544SStefano Zampini     label = NULL;
871f86f7544SStefano Zampini     if (enable_bmark) {
872f86f7544SStefano Zampini       PetscInt lminl = PETSC_MAX_INT;
873f86f7544SStefano Zampini 
8747bf4dd16SStefano Zampini       ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr);
875f86f7544SStefano Zampini       if (label) {
876f86f7544SStefano Zampini         IS       vals;
877f86f7544SStefano Zampini         PetscInt ldef;
878f86f7544SStefano Zampini 
879f86f7544SStefano Zampini         ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr);
880f86f7544SStefano Zampini         ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr);
881f86f7544SStefano Zampini         ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr);
882f86f7544SStefano Zampini         ierr = ISDestroy(&vals);CHKERRQ(ierr);
883f86f7544SStefano Zampini         lminl = PetscMin(ldef,lminl);
884f86f7544SStefano Zampini       }
885f86f7544SStefano Zampini       ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
886f86f7544SStefano Zampini       if (minl == PETSC_MAX_INT) minl = 1;
887f86f7544SStefano Zampini     }
8888135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
8898135c375SStefano Zampini     for (p=fStart;p<fEnd;p++) {
89077eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
8918135c375SStefano Zampini         const PetscInt *support;
89277eacf09SStefano Zampini         PetscInt       supportSize,c,nc = 0;
8938135c375SStefano Zampini 
8948135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
8958135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
8963e6c54aaSStefano Zampini         if (pown) {
89777eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
89877eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
89977eacf09SStefano Zampini               fcells[nc++] = support[c];
9008135c375SStefano Zampini             }
90177eacf09SStefano Zampini           }
90277eacf09SStefano Zampini         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
90377eacf09SStefano Zampini         for (c=0;c<nc;c++) {
904c3c203b2SStefano Zampini           const DMPolytopeType *faceTypes;
905c3c203b2SStefano Zampini           DMPolytopeType       cellType;
906c3c203b2SStefano Zampini           const PetscInt       *faceSizes,*cone;
907c3c203b2SStefano Zampini           PetscInt             vids[8],*faces,st,i,coneSize,cell,cl,nv,cid = -1,mid = -1;
9088135c375SStefano Zampini 
90977eacf09SStefano Zampini           cell = fcells[c];
91077eacf09SStefano Zampini           ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr);
911b135d7daSStefano Zampini           ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr);
912b135d7daSStefano Zampini           for (cl=0;cl<coneSize;cl++)
91377eacf09SStefano Zampini             if (cone[cl] == p)
9148135c375SStefano Zampini               break;
915b135d7daSStefano Zampini           if (cl == coneSize) continue;
9168135c375SStefano Zampini 
91777eacf09SStefano Zampini           /* face material id and type */
918f86f7544SStefano Zampini           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr);
9198135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
92077eacf09SStefano Zampini           /* vertex ids */
921c3c203b2SStefano Zampini           ierr = DMPlexGetCellType(dm,cell,&cellType);CHKERRQ(ierr);
92277eacf09SStefano Zampini           ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
923c3c203b2SStefano Zampini           ierr = DMPlexGetRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);CHKERRQ(ierr);
924c3c203b2SStefano Zampini           st = 0;
925c3c203b2SStefano Zampini           for (i=0;i<cl;i++) st += faceSizes[i];
926c3c203b2SStefano Zampini           ierr = DMPlexInvertCell(faceTypes[cl],faces + st);CHKERRQ(ierr);
927c3c203b2SStefano Zampini           for (i=0;i<faceSizes[cl];i++) {
928c3c203b2SStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer," %d",faces[st+i]);CHKERRQ(ierr);
929b135d7daSStefano Zampini           }
9308135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
931c3c203b2SStefano Zampini           ierr = DMPlexRestoreRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);CHKERRQ(ierr);
9323e96840aSStefano Zampini           bf -= 1;
93377eacf09SStefano Zampini         }
9348135c375SStefano Zampini       }
9358135c375SStefano Zampini     }
93677eacf09SStefano Zampini     if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
93777eacf09SStefano Zampini     ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr);
93877eacf09SStefano Zampini     ierr = PetscFree(fcells);CHKERRQ(ierr);
9398135c375SStefano Zampini   }
9408135c375SStefano Zampini 
9418135c375SStefano Zampini   /* mark owned vertices */
9423e6c54aaSStefano Zampini   vown = NULL;
9433e6c54aaSStefano Zampini   if (pown) {
9443e6c54aaSStefano Zampini     ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
9458135c375SStefano Zampini     for (p=cStart;p<cEnd;p++) {
9468135c375SStefano Zampini       PetscInt i,closureSize,*closure = NULL;
9478135c375SStefano Zampini 
9483e6c54aaSStefano Zampini       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
9498135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
9508135c375SStefano Zampini       for (i=0;i<closureSize;i++) {
9518135c375SStefano Zampini         const PetscInt pp = closure[2*i];
9528135c375SStefano Zampini 
9538135c375SStefano Zampini         if (pp >= vStart && pp < vEnd) {
9543e6c54aaSStefano Zampini           ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr);
9558135c375SStefano Zampini         }
9568135c375SStefano Zampini       }
9578135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
9588135c375SStefano Zampini     }
9598135c375SStefano Zampini   }
9608135c375SStefano Zampini 
9618135c375SStefano Zampini   /* vertex_parents (Non-conforming meshes) */
9628135c375SStefano Zampini   parentSection  = NULL;
9638135c375SStefano Zampini   if (enable_ncmesh) {
9648135c375SStefano Zampini     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
9658135c375SStefano Zampini   }
9668135c375SStefano Zampini   if (parentSection) {
9678135c375SStefano Zampini     PetscInt vp,gvp;
9688135c375SStefano Zampini 
9698135c375SStefano Zampini     for (vp=0,p=vStart;p<vEnd;p++) {
9708135c375SStefano Zampini       DMLabel  dlabel;
9718135c375SStefano Zampini       PetscInt parent,depth;
9728135c375SStefano Zampini 
9733e6c54aaSStefano Zampini       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
9748135c375SStefano Zampini       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
9758135c375SStefano Zampini       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
9768135c375SStefano Zampini       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
9778135c375SStefano Zampini       if (parent != p) vp++;
9788135c375SStefano Zampini     }
9798135c375SStefano Zampini     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
9808135c375SStefano Zampini     if (gvp) {
9818135c375SStefano Zampini       PetscInt  maxsupp;
9828135c375SStefano Zampini       PetscBool *skip = NULL;
9838135c375SStefano Zampini 
9848135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
9858135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
9868135c375SStefano Zampini       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
9878135c375SStefano Zampini       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
9888135c375SStefano Zampini       for (p=vStart;p<vEnd;p++) {
9898135c375SStefano Zampini         DMLabel  dlabel;
9908135c375SStefano Zampini         PetscInt parent;
9918135c375SStefano Zampini 
9923e6c54aaSStefano Zampini         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
9938135c375SStefano Zampini         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
9948135c375SStefano Zampini         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
9958135c375SStefano Zampini         if (parent != p) {
99696ca5757SLisandro Dalcin           PetscInt       vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
9973924b612SStefano Zampini           PetscInt       i,nv,ssize,n,numChildren,depth = -1;
9988135c375SStefano Zampini           const PetscInt *children;
9993924b612SStefano Zampini 
10003924b612SStefano Zampini           ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr);
10013924b612SStefano Zampini           switch (ssize) {
10028135c375SStefano Zampini             case 2: /* edge */
10038135c375SStefano Zampini               nv   = 0;
1004cc0d3ed7SStefano Zampini               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
10058135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
10068135c375SStefano Zampini               for (i=0;i<nv;i++) {
100796ca5757SLisandro Dalcin                 ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]);CHKERRQ(ierr);
10088135c375SStefano Zampini               }
10098135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
10108135c375SStefano Zampini               vp--;
10118135c375SStefano Zampini               break;
10128135c375SStefano Zampini             case 4: /* face */
10138135c375SStefano Zampini               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
10148135c375SStefano Zampini               for (n=0;n<numChildren;n++) {
10158135c375SStefano Zampini                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
10168135c375SStefano Zampini                 if (!depth) {
10178135c375SStefano Zampini                   const PetscInt *hvsupp,*hesupp,*cone;
10188135c375SStefano Zampini                   PetscInt       hvsuppSize,hesuppSize,coneSize;
1019451a39c7SStefano Zampini                   PetscInt       hv = children[n],he = -1,f;
10208135c375SStefano Zampini 
1021580bdb30SBarry Smith                   ierr = PetscArrayzero(skip,maxsupp);CHKERRQ(ierr);
10228135c375SStefano Zampini                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
10238135c375SStefano Zampini                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
10248135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
10258135c375SStefano Zampini                     PetscInt ep;
10268135c375SStefano Zampini                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
10278135c375SStefano Zampini                     if (ep != hvsupp[i]) {
10288135c375SStefano Zampini                       he = hvsupp[i];
10298135c375SStefano Zampini                     } else {
10308135c375SStefano Zampini                       skip[i] = PETSC_TRUE;
10318135c375SStefano Zampini                     }
10328135c375SStefano Zampini                   }
1033451a39c7SStefano Zampini                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
10348135c375SStefano Zampini                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
103596ca5757SLisandro Dalcin                   vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
10368135c375SStefano Zampini                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
10378135c375SStefano Zampini                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
10388135c375SStefano Zampini                   for (f=0;f<hesuppSize;f++) {
10398135c375SStefano Zampini                     PetscInt j;
10408135c375SStefano Zampini 
10418135c375SStefano Zampini                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
10428135c375SStefano Zampini                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
10438135c375SStefano Zampini                     for (j=0;j<coneSize;j++) {
10448135c375SStefano Zampini                       PetscInt k;
10458135c375SStefano Zampini                       for (k=0;k<hvsuppSize;k++) {
10468135c375SStefano Zampini                         if (hvsupp[k] == cone[j]) {
10478135c375SStefano Zampini                           skip[k] = PETSC_TRUE;
10488135c375SStefano Zampini                           break;
10498135c375SStefano Zampini                         }
10508135c375SStefano Zampini                       }
10518135c375SStefano Zampini                     }
10528135c375SStefano Zampini                   }
10538135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
10548135c375SStefano Zampini                     if (!skip[i]) {
10558135c375SStefano Zampini                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
105696ca5757SLisandro Dalcin                       vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
10578135c375SStefano Zampini                     }
10588135c375SStefano Zampini                   }
10598135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
10608135c375SStefano Zampini                   for (i=0;i<2;i++) {
106196ca5757SLisandro Dalcin                     ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]-vStart);CHKERRQ(ierr);
10628135c375SStefano Zampini                   }
10638135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
10648135c375SStefano Zampini                   vp--;
10658135c375SStefano Zampini                 }
10668135c375SStefano Zampini               }
10678135c375SStefano Zampini               break;
10688135c375SStefano Zampini             default:
10693924b612SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize);
10708135c375SStefano Zampini           }
10718135c375SStefano Zampini         }
10728135c375SStefano Zampini       }
10738135c375SStefano Zampini       ierr = PetscFree(skip);CHKERRQ(ierr);
10748135c375SStefano Zampini     }
10758135c375SStefano Zampini     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
10768135c375SStefano Zampini   }
10778135c375SStefano Zampini   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
10783e6c54aaSStefano Zampini   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
10798135c375SStefano Zampini 
10808135c375SStefano Zampini   /* vertices */
108177eacf09SStefano Zampini   if (hovec) { /* higher-order meshes */
108277eacf09SStefano Zampini     const char *fec;
10830286d493SLisandro Dalcin     PetscInt   i,n,s;
108477eacf09SStefano Zampini 
108577eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
108677eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
108777eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
108877eacf09SStefano Zampini     ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr);
108977eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
10903e7cab22SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr);
109177eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
109277eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
109377eacf09SStefano Zampini     ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr);
109477eacf09SStefano Zampini     ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr);
109577eacf09SStefano Zampini     if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
109677eacf09SStefano Zampini     for (i=0;i<n/sdim;i++) {
109777eacf09SStefano Zampini       for (s=0;s<sdim;s++) {
109877eacf09SStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr);
109977eacf09SStefano Zampini       }
110077eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
110177eacf09SStefano Zampini     }
110277eacf09SStefano Zampini     ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr);
110377eacf09SStefano Zampini   } else {
1104cc0d3ed7SStefano Zampini     ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr);
11058135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
1106cc0d3ed7SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr);
11078135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
11088135c375SStefano Zampini     ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
1109cc0d3ed7SStefano Zampini     for (p=0;p<nvert/sdim;p++) {
1110cc0d3ed7SStefano Zampini       PetscInt s;
1111cc0d3ed7SStefano Zampini       for (s=0;s<sdim;s++) {
11123e96840aSStefano Zampini         PetscReal v = PetscRealPart(array[p*sdim+s]);
11133e96840aSStefano Zampini 
11143e96840aSStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);CHKERRQ(ierr);
11158135c375SStefano Zampini       }
1116cc0d3ed7SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
11178135c375SStefano Zampini     }
11188135c375SStefano Zampini     ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
111977eacf09SStefano Zampini   }
1120*066ea43fSLisandro Dalcin   ierr = VecDestroy(&hovec);CHKERRQ(ierr);
11218135c375SStefano Zampini   PetscFunctionReturn(0);
11228135c375SStefano Zampini }
11238135c375SStefano Zampini 
11240286d493SLisandro Dalcin PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
11258135c375SStefano Zampini {
11268135c375SStefano Zampini   PetscErrorCode ierr;
11278135c375SStefano Zampini   PetscFunctionBegin;
11280286d493SLisandro Dalcin   ierr = DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);CHKERRQ(ierr);
11298135c375SStefano Zampini   PetscFunctionReturn(0);
11308135c375SStefano Zampini }
1131