xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 7bf4dd16b704174a5c32122f284e36aa4ab43524) !
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;
5511a4995dSStefano Zampini   PetscInt       dim,cStart,cEnd,cEndInterior,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);
648135c375SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
658135c375SStefano Zampini   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
668135c375SStefano Zampini   ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
678135c375SStefano Zampini   ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
688135c375SStefano Zampini   ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
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++;
748135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
758135c375SStefano Zampini       for (i=0;i<numPoints*2;i+= 2) {
768135c375SStefano Zampini         if ((points[i] >= vStart) && (points[i] < vEnd)) {
778135c375SStefano Zampini           ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr);
788135c375SStefano Zampini         }
798135c375SStefano Zampini       }
808135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
818135c375SStefano Zampini     }
828135c375SStefano Zampini   }
8377eacf09SStefano Zampini   for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++;
848135c375SStefano Zampini 
858135c375SStefano Zampini   ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr);
868135c375SStefano Zampini   ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr);
878135c375SStefano Zampini   ierr = DMGetDefaultSection(dm,&s);CHKERRQ(ierr);
888135c375SStefano Zampini   ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr);
898135c375SStefano Zampini   for (f=0,maxfields=0;f<nfields;f++) {
908135c375SStefano Zampini     PetscInt bs;
918135c375SStefano Zampini 
928135c375SStefano Zampini     ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr);
938135c375SStefano Zampini     maxfields += bs;
948135c375SStefano Zampini   }
954cac2994SStefano Zampini   ierr = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);CHKERRQ(ierr);
968135c375SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
978135c375SStefano Zampini   ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr);
988135c375SStefano Zampini   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
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 
1068135c375SStefano Zampini       ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr);
1078135c375SStefano Zampini       ierr = PetscStrlen(fname,&len);CHKERRQ(ierr);
1088135c375SStefano Zampini       if (len) {
1098135c375SStefano Zampini         ierr = PetscStrcpy(name,fname);CHKERRQ(ierr);
1108135c375SStefano Zampini       } else {
1118135c375SStefano Zampini         ierr = PetscSNPrintf(name,256,"Field%d",f);CHKERRQ(ierr);
1128135c375SStefano Zampini       }
1138135c375SStefano Zampini       ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
1148135c375SStefano Zampini       if (disc) {
1158135c375SStefano Zampini         PetscClassId id;
1168135c375SStefano Zampini         PetscInt     Nc;
1178135c375SStefano Zampini         char         fec[64];
1188135c375SStefano Zampini 
1198135c375SStefano Zampini         ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
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 
1278135c375SStefano Zampini           ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr);
1288135c375SStefano Zampini           ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr);
1298135c375SStefano Zampini           ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr);
1308135c375SStefano Zampini           ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr);
1318135c375SStefano Zampini           if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
1328135c375SStefano Zampini           ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr);
1338135c375SStefano Zampini           if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported");
1348135c375SStefano Zampini           ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
1358135c375SStefano Zampini           if (continuous && order > 0) {
1368135c375SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%dD_P1",dim);CHKERRQ(ierr);
1378135c375SStefano Zampini           } else {
1388135c375SStefano Zampini             H1   = PETSC_FALSE;
1398135c375SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P%d",dim,order);CHKERRQ(ierr);
1408135c375SStefano Zampini           }
1418135c375SStefano Zampini           ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr);
1428135c375SStefano Zampini           bs[ctx->nf]   = Nc;
143bb77a09fSStefano Zampini           dims[ctx->nf] = dim;
1448135c375SStefano Zampini           if (H1) {
1458135c375SStefano Zampini             nlocal[ctx->nf] = Nc * Nv;
1468135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
1478135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr);
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;
1528135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr);
1538135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1548135c375SStefano Zampini             }
1558135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
1568135c375SStefano Zampini           } else {
1578135c375SStefano Zampini             nlocal[ctx->nf] = Nc * totc;
1588135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
1598135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr);
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;
1648135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
1658135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
1668135c375SStefano Zampini             }
1678135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
1688135c375SStefano Zampini           }
1698135c375SStefano Zampini           ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
1708135c375SStefano Zampini           ierr = VecDestroy(&xfield);CHKERRQ(ierr);
1718135c375SStefano Zampini           ierr = ISDestroy(&isfield);CHKERRQ(ierr);
1728135c375SStefano Zampini           ctx->nf++;
1738135c375SStefano Zampini         } else if (id == PETSCFV_CLASSID) {
1748135c375SStefano Zampini           PetscInt c;
1758135c375SStefano Zampini 
1768135c375SStefano Zampini           ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr);
1778135c375SStefano Zampini           ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P0",dim);CHKERRQ(ierr);
1788135c375SStefano Zampini           for (c = 0; c < Nc; c++) {
1798135c375SStefano Zampini             char comp[256];
1808135c375SStefano Zampini             ierr = PetscSNPrintf(comp,256,"%s-Comp%d",name,c);CHKERRQ(ierr);
1818135c375SStefano Zampini             ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr);
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;
1858135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
1868135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr);
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;
1918135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
1928135c375SStefano Zampini               idxs[cum++] = off + c;
1938135c375SStefano Zampini             }
1948135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
1958135c375SStefano Zampini             ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
1968135c375SStefano Zampini             ierr = VecDestroy(&xfield);CHKERRQ(ierr);
1978135c375SStefano Zampini             ierr = ISDestroy(&isfield);CHKERRQ(ierr);
1988135c375SStefano Zampini             ctx->nf++;
1998135c375SStefano Zampini           }
2008135c375SStefano Zampini         } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f);
2018135c375SStefano Zampini       } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f);
2028135c375SStefano Zampini     }
2038135c375SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
2048135c375SStefano Zampini   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
2058135c375SStefano Zampini   ierr = VecDestroy(&xlocal);CHKERRQ(ierr);
2068135c375SStefano Zampini   ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
2078135c375SStefano Zampini 
2084cac2994SStefano Zampini   /* create work vectors */
2094cac2994SStefano Zampini   for (f=0;f<ctx->nf;f++) {
2104cac2994SStefano Zampini     ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr);
2114cac2994SStefano Zampini     ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr);
2124cac2994SStefano Zampini     ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr);
2134cac2994SStefano Zampini     ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr);
2144cac2994SStefano Zampini   }
2154cac2994SStefano Zampini 
2168135c375SStefano Zampini   /* customize the viewer */
2174cac2994SStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
2188135c375SStefano Zampini   for (f=0;f<ctx->nf;f++) {
2198135c375SStefano Zampini     ierr = PetscFree(fieldname[f]);CHKERRQ(ierr);
2208135c375SStefano Zampini     ierr = PetscFree(fec_type[f]);CHKERRQ(ierr);
2214cac2994SStefano Zampini     ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr);
2228135c375SStefano Zampini   }
2234cac2994SStefano Zampini   ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr);
2248135c375SStefano Zampini   PetscFunctionReturn(0);
2258135c375SStefano Zampini }
2268135c375SStefano Zampini 
2278135c375SStefano Zampini typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,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},
2328135c375SStefano Zampini                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } };
2338135c375SStefano Zampini 
234cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid)
2358135c375SStefano Zampini {
2368135c375SStefano Zampini   DMLabel        dlabel;
2378135c375SStefano Zampini   PetscInt       depth,csize;
2388135c375SStefano Zampini   PetscErrorCode ierr;
2398135c375SStefano Zampini 
2408135c375SStefano Zampini   PetscFunctionBegin;
2418135c375SStefano Zampini   ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
2428135c375SStefano Zampini   ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
2438135c375SStefano Zampini   ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr);
2448135c375SStefano Zampini   if (label) {
2458135c375SStefano Zampini     ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr);
24677eacf09SStefano Zampini   } else *mid = 1;
2478135c375SStefano Zampini   *cid = mfem_table_cid[depth][csize];
2488135c375SStefano Zampini   PetscFunctionReturn(0);
2498135c375SStefano Zampini }
2508135c375SStefano Zampini 
251cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[])
2528135c375SStefano Zampini {
253cc0d3ed7SStefano Zampini   PetscInt       dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL;
2548135c375SStefano Zampini   PetscErrorCode ierr;
2558135c375SStefano Zampini 
2568135c375SStefano Zampini   PetscFunctionBegin;
2578135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
258cc0d3ed7SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
259cc0d3ed7SStefano Zampini   sdim = dim;
260cc0d3ed7SStefano Zampini   if (csec) {
261cc0d3ed7SStefano Zampini     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
262cc0d3ed7SStefano Zampini     ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr);
263cc0d3ed7SStefano Zampini     off  = off/sdim;
264cc0d3ed7SStefano Zampini     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
265cc0d3ed7SStefano Zampini   }
266cc0d3ed7SStefano Zampini   if (!dof) {
2678135c375SStefano Zampini     ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
2688135c375SStefano Zampini     for (i=0,q=0;i<numPoints*2;i+= 2)
2698135c375SStefano Zampini       if ((points[i] >= vStart) && (points[i] < vEnd))
270cc0d3ed7SStefano Zampini         vids[q++] = points[i]-vStart+off;
2718135c375SStefano Zampini     ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
272cc0d3ed7SStefano Zampini   } else {
273cc0d3ed7SStefano Zampini     ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr);
274cc0d3ed7SStefano Zampini     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
275cc0d3ed7SStefano Zampini     for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q;
276cc0d3ed7SStefano Zampini   }
2778135c375SStefano Zampini   *nv  = q;
2788135c375SStefano Zampini   PetscFunctionReturn(0);
2798135c375SStefano Zampini }
2808135c375SStefano Zampini 
28177eacf09SStefano Zampini /*
28277eacf09SStefano Zampini    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
28377eacf09SStefano Zampini    Higher order meshes are also supported
28477eacf09SStefano Zampini */
2858135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
2868135c375SStefano Zampini {
2878135c375SStefano Zampini   DMLabel              label;
2888135c375SStefano Zampini   PetscSection         coordSection,parentSection;
28977eacf09SStefano Zampini   Vec                  coordinates,hovec;
2908135c375SStefano Zampini   const PetscScalar    *array;
2918135c375SStefano Zampini   PetscInt             bf,p,sdim,dim,depth,novl;
292cc0d3ed7SStefano Zampini   PetscInt             cStart,cEnd,cEndInterior,vStart,vEnd,nvert;
2938135c375SStefano Zampini   PetscMPIInt          commsize;
2943e6c54aaSStefano Zampini   PetscBool            localized,isascii;
2953e6c54aaSStefano Zampini   PetscBool            enable_mfem,enable_boundary,enable_ncmesh;
2963e6c54aaSStefano Zampini   PetscBT              pown,vown;
2978135c375SStefano Zampini   PetscErrorCode       ierr;
2988135c375SStefano Zampini   PetscContainer       glvis_container;
29977eacf09SStefano Zampini   PetscBool            periodic, enabled = PETSC_TRUE;
30077eacf09SStefano Zampini   const char           *fmt;
301*7bf4dd16SStefano Zampini   char                 emark[64] = "",bmark[64] = "";
3028135c375SStefano Zampini 
3038135c375SStefano Zampini   PetscFunctionBegin;
3048135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3058135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3068135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
3078135c375SStefano Zampini   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
3088135c375SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
3098135c375SStefano Zampini   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
3108135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
3118135c375SStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
3128135c375SStefano Zampini   if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
3138135c375SStefano Zampini 
3148135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
3158135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
3168135c375SStefano Zampini   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
3178135c375SStefano Zampini   {
3188135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
3198135c375SStefano Zampini     ierr    = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
3208135c375SStefano Zampini     enabled = glvis_info->enabled;
32177eacf09SStefano Zampini     fmt     = glvis_info->fmt;
3228135c375SStefano Zampini   }
323cc0d3ed7SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3248135c375SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
3258135c375SStefano Zampini   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
3268135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
32777eacf09SStefano Zampini   ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr);
328cc0d3ed7SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
32977eacf09SStefano Zampini   if (periodic && !localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
330cc0d3ed7SStefano Zampini   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
33177eacf09SStefano Zampini   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
33277eacf09SStefano Zampini   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
33377eacf09SStefano Zampini 
33477eacf09SStefano Zampini   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
33577eacf09SStefano Zampini      DMPlex does not currently support HO meshes, so there's no API for this */
33677eacf09SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr);
3378135c375SStefano Zampini 
3388135c375SStefano Zampini   /*
3398135c375SStefano Zampini      a couple of sections of the mesh specification are disabled
34077eacf09SStefano Zampini        - boundary: unless we want to visualize boundary attributes or we have a periodic mesh
34177eacf09SStefano Zampini                    the boundary is not needed for proper mesh visualization
34277eacf09SStefano Zampini        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
3433e6c54aaSStefano Zampini                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
3448135c375SStefano Zampini   */
345e74666bcSStefano Zampini   enable_boundary = periodic;
3468135c375SStefano Zampini   enable_ncmesh   = PETSC_FALSE;
3473e6c54aaSStefano Zampini   enable_mfem     = PETSC_FALSE;
348*7bf4dd16SStefano Zampini   /* I'm tired of problems with negative values in the markers, disable them */
3498135c375SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
3503e6c54aaSStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr);
3513e6c54aaSStefano 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);
3523e6c54aaSStefano 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);
353*7bf4dd16SStefano Zampini   ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),NULL);CHKERRQ(ierr);
354*7bf4dd16SStefano Zampini   ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),NULL);CHKERRQ(ierr);
3558135c375SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
3563e6c54aaSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
3573e6c54aaSStefano Zampini   if (enable_ncmesh && commsize > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
3588135c375SStefano Zampini 
3598135c375SStefano Zampini   /* Identify possible cells in the overlap */
3608135c375SStefano Zampini   novl = 0;
3618135c375SStefano Zampini   pown = NULL;
3628135c375SStefano Zampini   if (commsize > 1) {
3633e6c54aaSStefano Zampini     IS             globalNum = NULL;
3643e6c54aaSStefano Zampini     const PetscInt *gNum;
3653e6c54aaSStefano Zampini     PetscBool      ovl  = PETSC_FALSE;
3663e6c54aaSStefano Zampini 
3678135c375SStefano Zampini     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
3688135c375SStefano Zampini     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
3698135c375SStefano Zampini     for (p=cStart; p<cEnd; p++) {
3708135c375SStefano Zampini       if (gNum[p-cStart] < 0) {
3718135c375SStefano Zampini         ovl = PETSC_TRUE;
3728135c375SStefano Zampini         novl++;
3738135c375SStefano Zampini       }
3748135c375SStefano Zampini     }
3758135c375SStefano Zampini     if (ovl) {
3768135c375SStefano Zampini       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
3778135c375SStefano Zampini          TODO: garbage collector? attach pown to dm?  */
3783e6c54aaSStefano Zampini       ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr);
3793e6c54aaSStefano Zampini       for (p=cStart; p<cEnd; p++) {
3803e6c54aaSStefano Zampini         if (gNum[p-cStart] < 0) continue;
3813e6c54aaSStefano Zampini         else {
3823e6c54aaSStefano Zampini           ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
3838135c375SStefano Zampini         }
3848135c375SStefano Zampini       }
3853e6c54aaSStefano Zampini     }
3863e6c54aaSStefano Zampini     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
3873e6c54aaSStefano Zampini   }
3888135c375SStefano Zampini 
3893e6c54aaSStefano Zampini   /* return if this process is disabled */
3908135c375SStefano Zampini   if (!enabled) {
3918135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
3928135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
3938135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
3948135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
3958135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3968135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
3978135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3988135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
3998135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
4008135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
4018135c375SStefano Zampini     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
4028135c375SStefano Zampini     PetscFunctionReturn(0);
4038135c375SStefano Zampini   }
4048135c375SStefano Zampini 
4053e6c54aaSStefano Zampini   if (enable_mfem) {
4063e6c54aaSStefano Zampini     if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
4073e6c54aaSStefano Zampini       PetscInt    vpc = 0;
4083e6c54aaSStefano Zampini       char        fec[64];
4093e6c54aaSStefano Zampini       int         vids[8] = {0,1,2,3,4,5,6,7};
4103e6c54aaSStefano Zampini       int         hexv[8] = {0,1,3,2,4,5,7,6},*dof;
4113e6c54aaSStefano Zampini       PetscScalar *array,*ptr;
4123e6c54aaSStefano Zampini 
4133e6c54aaSStefano Zampini       ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%dD_P1",dim);CHKERRQ(ierr);
4143e6c54aaSStefano Zampini       if (cEnd-cStart) {
4153e6c54aaSStefano Zampini         PetscInt fpc;
4163e6c54aaSStefano Zampini 
4173e6c54aaSStefano Zampini         ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
4183e6c54aaSStefano Zampini         switch(dim) {
4193e6c54aaSStefano Zampini           case 1:
4203e6c54aaSStefano Zampini             vpc = 2;
4213e6c54aaSStefano Zampini             dof = hexv;
4223e6c54aaSStefano Zampini             break;
4233e6c54aaSStefano Zampini           case 2:
4243e6c54aaSStefano Zampini             switch (fpc) {
4253e6c54aaSStefano Zampini               case 3:
4263e6c54aaSStefano Zampini                 vpc = 3;
4273e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4283e6c54aaSStefano Zampini                 break;
4293e6c54aaSStefano Zampini               case 4:
4303e6c54aaSStefano Zampini                 vpc = 4;
4313e6c54aaSStefano Zampini                 dof = hexv;
4323e6c54aaSStefano Zampini                 break;
4333e6c54aaSStefano Zampini               default:
4343e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4353e6c54aaSStefano Zampini                 break;
4363e6c54aaSStefano Zampini             }
4373e6c54aaSStefano Zampini             break;
4383e6c54aaSStefano Zampini           case 3:
4393e6c54aaSStefano Zampini             switch (fpc) {
4403e6c54aaSStefano Zampini               case 4:
4413e6c54aaSStefano Zampini                 vpc = 4;
4423e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4433e6c54aaSStefano Zampini                 break;
4443e6c54aaSStefano Zampini               case 6:
4453e6c54aaSStefano Zampini                 vpc = 8;
4463e6c54aaSStefano Zampini                 dof = hexv;
4473e6c54aaSStefano Zampini                 break;
4483e6c54aaSStefano Zampini               default:
4493e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4503e6c54aaSStefano Zampini                 break;
4513e6c54aaSStefano Zampini             }
4523e6c54aaSStefano Zampini             break;
4533e6c54aaSStefano Zampini           default:
4543e6c54aaSStefano Zampini             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
4553e6c54aaSStefano Zampini             break;
4563e6c54aaSStefano Zampini         }
4573e6c54aaSStefano Zampini         ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr);
4583e6c54aaSStefano Zampini       }
4593e6c54aaSStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr);
4603e6c54aaSStefano Zampini       ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr);
4613e6c54aaSStefano Zampini       ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr);
4623e6c54aaSStefano Zampini       ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr);
4633e6c54aaSStefano Zampini       ierr = VecGetArray(hovec,&array);CHKERRQ(ierr);
4643e6c54aaSStefano Zampini       ptr  = array;
4653e6c54aaSStefano Zampini       for (p=cStart;p<cEnd;p++) {
4663e6c54aaSStefano Zampini         PetscInt    csize,v,d;
4673e6c54aaSStefano Zampini         PetscScalar *vals = NULL;
4683e6c54aaSStefano Zampini 
4693e6c54aaSStefano Zampini         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
4703e6c54aaSStefano Zampini         ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
4713e6c54aaSStefano 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);
4723e6c54aaSStefano Zampini         for (v=0;v<vpc;v++) {
4733e6c54aaSStefano Zampini           for (d=0;d<sdim;d++) {
4743e6c54aaSStefano Zampini             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
4753e6c54aaSStefano Zampini           }
4763e6c54aaSStefano Zampini         }
4773e6c54aaSStefano Zampini         ptr += vpc*sdim;
4783e6c54aaSStefano Zampini         ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
4793e6c54aaSStefano Zampini       }
4803e6c54aaSStefano Zampini       ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr);
4813e6c54aaSStefano Zampini     }
4823e6c54aaSStefano Zampini   }
4833e6c54aaSStefano Zampini 
4843e6c54aaSStefano Zampini 
4858135c375SStefano Zampini   /* header */
4868135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
4878135c375SStefano Zampini 
4888135c375SStefano Zampini   /* topological dimension */
4898135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
4908135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
4918135c375SStefano Zampini 
4928135c375SStefano Zampini   /* elements */
493*7bf4dd16SStefano Zampini   ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr);
4948135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
4958135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
4968135c375SStefano Zampini   for (p=cStart;p<cEnd;p++) {
49711a4995dSStefano Zampini     int      vids[8];
49811a4995dSStefano Zampini     PetscInt i,nv = 0,cid = -1,mid = 1;
4998135c375SStefano Zampini 
5003e6c54aaSStefano Zampini     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
5018135c375SStefano Zampini     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
50277eacf09SStefano Zampini     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
50377eacf09SStefano Zampini     ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
5048135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
5058135c375SStefano Zampini     for (i=0;i<nv;i++) {
50611a4995dSStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
5078135c375SStefano Zampini     }
5088135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
5098135c375SStefano Zampini   }
5108135c375SStefano Zampini 
511cc0d3ed7SStefano Zampini   /* boundary */
512cc0d3ed7SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
513cc0d3ed7SStefano Zampini   if (!enable_boundary) {
514cc0d3ed7SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
515cc0d3ed7SStefano Zampini   } else {
51677eacf09SStefano Zampini     DMLabel  perLabel;
51777eacf09SStefano Zampini     PetscBT  bfaces;
51877eacf09SStefano Zampini     PetscInt fStart,fEnd,fEndInterior,*fcells;
51977eacf09SStefano Zampini     PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0;
520cc0d3ed7SStefano Zampini     PetscInt fv1[]     = {0,1},
521cc0d3ed7SStefano Zampini              fv2tri[]  = {0,1,
522cc0d3ed7SStefano Zampini                           1,2,
523cc0d3ed7SStefano Zampini                           2,0},
524cc0d3ed7SStefano Zampini              fv2quad[] = {0,1,
525cc0d3ed7SStefano Zampini                           1,2,
526cc0d3ed7SStefano Zampini                           2,3,
527cc0d3ed7SStefano Zampini                           3,0},
528cc0d3ed7SStefano Zampini              fv3tet[]  = {0,1,2,
529cc0d3ed7SStefano Zampini                           0,3,1,
530cc0d3ed7SStefano Zampini                           0,2,3,
531cc0d3ed7SStefano Zampini                           2,1,3},
532cc0d3ed7SStefano Zampini              fv3hex[]  = {0,1,2,3,
533cc0d3ed7SStefano Zampini                        4,5,6,7,
534cc0d3ed7SStefano Zampini                        0,3,5,4,
535cc0d3ed7SStefano Zampini                        2,1,7,6,
536cc0d3ed7SStefano Zampini                        3,2,6,5,
537cc0d3ed7SStefano Zampini                        0,4,7,1};
538cc0d3ed7SStefano Zampini 
5398135c375SStefano Zampini     /* determine orientation of boundary mesh */
5408135c375SStefano Zampini     if (cEnd-cStart) {
5418135c375SStefano Zampini       ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
5428135c375SStefano Zampini       switch(dim) {
5438135c375SStefano Zampini         case 1:
5448135c375SStefano Zampini           if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
5458135c375SStefano Zampini           faces = fv1;
5468135c375SStefano Zampini           vpf = 1;
54777eacf09SStefano Zampini           vpc = 2;
5488135c375SStefano Zampini           break;
5498135c375SStefano Zampini         case 2:
5508135c375SStefano Zampini           switch (fpc) {
5518135c375SStefano Zampini             case 3:
5528135c375SStefano Zampini               faces = fv2tri;
5538135c375SStefano Zampini               vpf   = 2;
55477eacf09SStefano Zampini               vpc   = 3;
5558135c375SStefano Zampini               break;
5568135c375SStefano Zampini             case 4:
5578135c375SStefano Zampini               faces = fv2quad;
5588135c375SStefano Zampini               vpf   = 2;
55977eacf09SStefano Zampini               vpc   = 4;
5608135c375SStefano Zampini               break;
5618135c375SStefano Zampini             default:
5628135c375SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
5638135c375SStefano Zampini               break;
5648135c375SStefano Zampini           }
5658135c375SStefano Zampini           break;
5668135c375SStefano Zampini         case 3:
5678135c375SStefano Zampini           switch (fpc) {
5688135c375SStefano Zampini             case 4:
5698135c375SStefano Zampini               faces = fv3tet;
5708135c375SStefano Zampini               vpf   = 3;
57177eacf09SStefano Zampini               vpc   = 4;
5728135c375SStefano Zampini               break;
5738135c375SStefano Zampini             case 6:
5748135c375SStefano Zampini               faces = fv3hex;
5758135c375SStefano Zampini               vpf   = 4;
57677eacf09SStefano Zampini               vpc   = 8;
5778135c375SStefano Zampini               break;
5788135c375SStefano Zampini             default:
5798135c375SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
5808135c375SStefano Zampini               break;
5818135c375SStefano Zampini           }
5828135c375SStefano Zampini           break;
5838135c375SStefano Zampini         default:
5848135c375SStefano Zampini           SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
5858135c375SStefano Zampini           break;
5868135c375SStefano Zampini       }
5878135c375SStefano Zampini     }
588cc0d3ed7SStefano Zampini     ierr = DMPlexGetHybridBounds(dm,NULL,&fEndInterior,NULL,NULL);CHKERRQ(ierr);
5898135c375SStefano Zampini     ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
5908135c375SStefano Zampini     fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
59177eacf09SStefano Zampini     ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr);
59277eacf09SStefano Zampini     ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr);
59377eacf09SStefano Zampini     ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr);
59477eacf09SStefano Zampini     ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
59577eacf09SStefano Zampini     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
59677eacf09SStefano Zampini       ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr);
59777eacf09SStefano Zampini       ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
59877eacf09SStefano Zampini       ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr);
59977eacf09SStefano Zampini       for (p=cStart;p<cEnd;p++) {
60077eacf09SStefano Zampini         PetscInt dof;
60177eacf09SStefano Zampini         ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
60277eacf09SStefano Zampini         if (dof) {
60377eacf09SStefano Zampini           PetscInt    v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
60477eacf09SStefano Zampini           PetscScalar *vals = NULL;
60577eacf09SStefano Zampini 
60677eacf09SStefano Zampini           if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
607*7bf4dd16SStefano Zampini           if (dof/sdim != vpc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible number of cell dofs %D, vertices %D and space dim %D",dof/sdim,vpc,sdim);
60877eacf09SStefano Zampini           ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
60977eacf09SStefano Zampini           ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
61077eacf09SStefano Zampini           for (v=0;v<cellClosureSize;v++)
61177eacf09SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
61277eacf09SStefano Zampini               vidxs = cellClosure + 2*v;
61377eacf09SStefano Zampini               break;
61477eacf09SStefano Zampini             }
61577eacf09SStefano Zampini           if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
61677eacf09SStefano Zampini           for (v=0;v<vpc;v++) {
61777eacf09SStefano Zampini             PetscInt s;
61877eacf09SStefano Zampini             for (s=0;s<sdim;s++) {
61977eacf09SStefano Zampini               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+vpc*sdim])>PETSC_MACHINE_EPSILON) {
62077eacf09SStefano Zampini                 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr);
62177eacf09SStefano Zampini               }
62277eacf09SStefano Zampini             }
62377eacf09SStefano Zampini           }
62477eacf09SStefano Zampini           ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
62577eacf09SStefano Zampini           ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
62677eacf09SStefano Zampini         }
62777eacf09SStefano Zampini       }
62877eacf09SStefano Zampini       if (dim > 1) {
62977eacf09SStefano Zampini         PetscInt eEnd,eStart,eEndInterior;
63077eacf09SStefano Zampini         ierr = DMPlexGetHybridBounds(dm,NULL,NULL,&eEndInterior,NULL);CHKERRQ(ierr);
63177eacf09SStefano Zampini         ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr);
63277eacf09SStefano Zampini         eEnd = eEndInterior < 0 ? eEnd : eEndInterior;
63377eacf09SStefano Zampini         for (p=eStart;p<eEnd;p++) {
63477eacf09SStefano Zampini           const PetscInt *cone;
63577eacf09SStefano Zampini           PetscInt       coneSize,i;
63677eacf09SStefano Zampini           PetscBool      ispe = PETSC_TRUE;
63777eacf09SStefano Zampini 
63877eacf09SStefano Zampini           ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
63977eacf09SStefano Zampini           ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
64077eacf09SStefano Zampini           for (i=0;i<coneSize;i++) {
64177eacf09SStefano Zampini             PetscInt v;
64277eacf09SStefano Zampini 
64377eacf09SStefano Zampini             ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
64477eacf09SStefano Zampini             ispe = (PetscBool)(ispe && (v==2));
64577eacf09SStefano Zampini           }
64677eacf09SStefano Zampini           if (ispe && coneSize) {
64777eacf09SStefano Zampini             ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
64877eacf09SStefano Zampini           }
64977eacf09SStefano Zampini         }
65077eacf09SStefano Zampini         if (dim > 2) {
65177eacf09SStefano Zampini           for (p=fStart;p<fEnd;p++) {
65277eacf09SStefano Zampini             const PetscInt *cone;
65377eacf09SStefano Zampini             PetscInt       coneSize,i;
65477eacf09SStefano Zampini             PetscBool      ispe = PETSC_TRUE;
65577eacf09SStefano Zampini 
65677eacf09SStefano Zampini             ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
65777eacf09SStefano Zampini             ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
65877eacf09SStefano Zampini             for (i=0;i<coneSize;i++) {
65977eacf09SStefano Zampini               PetscInt v;
66077eacf09SStefano Zampini 
66177eacf09SStefano Zampini               ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
66277eacf09SStefano Zampini               ispe = (PetscBool)(ispe && (v==2));
66377eacf09SStefano Zampini             }
66477eacf09SStefano Zampini             if (ispe && coneSize) {
66577eacf09SStefano Zampini               ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
66677eacf09SStefano Zampini             }
66777eacf09SStefano Zampini           }
66877eacf09SStefano Zampini         }
66977eacf09SStefano Zampini       }
67077eacf09SStefano Zampini     }
67177eacf09SStefano Zampini     for (p=fStart;p<fEnd;p++) {
67277eacf09SStefano Zampini       const PetscInt *support;
6738135c375SStefano Zampini       PetscInt       supportSize;
67477eacf09SStefano Zampini       PetscBool      isbf = PETSC_FALSE;
6758135c375SStefano Zampini 
6768135c375SStefano Zampini       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
6773e6c54aaSStefano Zampini       if (pown) {
6788135c375SStefano Zampini         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
67977eacf09SStefano Zampini         PetscInt  i;
68077eacf09SStefano Zampini 
68177eacf09SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
68277eacf09SStefano Zampini         for (i=0;i<supportSize;i++) {
68377eacf09SStefano Zampini           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
68477eacf09SStefano Zampini           else has_ghost = PETSC_TRUE;
68577eacf09SStefano Zampini         }
68677eacf09SStefano Zampini         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
68777eacf09SStefano Zampini       } else {
68877eacf09SStefano Zampini         isbf = (PetscBool)(supportSize == 1);
68977eacf09SStefano Zampini       }
69077eacf09SStefano Zampini       if (!isbf && perLabel) {
69177eacf09SStefano Zampini         const PetscInt *cone;
69277eacf09SStefano Zampini         PetscInt       coneSize,i;
69377eacf09SStefano Zampini 
69477eacf09SStefano Zampini         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
69577eacf09SStefano Zampini         ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
69677eacf09SStefano Zampini         isbf = PETSC_TRUE;
69777eacf09SStefano Zampini         for (i=0;i<coneSize;i++) {
69877eacf09SStefano Zampini           PetscInt v,d;
69977eacf09SStefano Zampini 
70077eacf09SStefano Zampini           ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
70177eacf09SStefano Zampini           ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr);
70277eacf09SStefano Zampini           isbf = (PetscBool)(isbf && v != d);
70377eacf09SStefano Zampini         }
70477eacf09SStefano Zampini       }
70577eacf09SStefano Zampini       if (isbf) {
70677eacf09SStefano Zampini         ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr);
70777eacf09SStefano Zampini       }
70877eacf09SStefano Zampini     }
70977eacf09SStefano Zampini     /* count boundary faces */
71077eacf09SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
71177eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
71277eacf09SStefano Zampini         const PetscInt *support;
71377eacf09SStefano Zampini         PetscInt       supportSize,c;
7148135c375SStefano Zampini 
7158135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
7168135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
7173e6c54aaSStefano Zampini         if (pown) {
71877eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
71977eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
72077eacf09SStefano Zampini               bf++;
7218135c375SStefano Zampini             }
72277eacf09SStefano Zampini           }
72377eacf09SStefano Zampini         } else bf += supportSize;
7248135c375SStefano Zampini       }
7258135c375SStefano Zampini     }
726*7bf4dd16SStefano Zampini     ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr);
7278135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
7288135c375SStefano Zampini     for (p=fStart;p<fEnd;p++) {
72977eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
7308135c375SStefano Zampini         const PetscInt *support;
73177eacf09SStefano Zampini         PetscInt       supportSize,c,nc = 0;
7328135c375SStefano Zampini 
7338135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
7348135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
7353e6c54aaSStefano Zampini         if (pown) {
73677eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
73777eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
73877eacf09SStefano Zampini               fcells[nc++] = support[c];
7398135c375SStefano Zampini             }
74077eacf09SStefano Zampini           }
74177eacf09SStefano Zampini         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
74277eacf09SStefano Zampini         for (c=0;c<nc;c++) {
74377eacf09SStefano Zampini           const    PetscInt *cone;
74477eacf09SStefano Zampini           int      vids[8];
74577eacf09SStefano Zampini           PetscInt i,cell,cl,nv,cid = -1,mid = -1;
7468135c375SStefano Zampini 
74777eacf09SStefano Zampini           cell = fcells[c];
74877eacf09SStefano Zampini           ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr);
74977eacf09SStefano Zampini           for (cl=0;cl<fpc;cl++)
75077eacf09SStefano Zampini             if (cone[cl] == p)
7518135c375SStefano Zampini               break;
7528135c375SStefano Zampini 
75377eacf09SStefano Zampini           /* face material id and type */
7548135c375SStefano Zampini           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
7558135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
75677eacf09SStefano Zampini           /* vertex ids */
75777eacf09SStefano Zampini           ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
7588135c375SStefano Zampini           for (i=0;i<vpf;i++) {
75977eacf09SStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer," %D",vids[faces[cl*vpf+i]]);CHKERRQ(ierr);
7608135c375SStefano Zampini           }
7618135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
76277eacf09SStefano Zampini         }
76377eacf09SStefano Zampini         bf = bf-nc;
7648135c375SStefano Zampini       }
7658135c375SStefano Zampini     }
76677eacf09SStefano Zampini     if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
76777eacf09SStefano Zampini     ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr);
76877eacf09SStefano Zampini     ierr = PetscFree(fcells);CHKERRQ(ierr);
7698135c375SStefano Zampini   }
7708135c375SStefano Zampini 
7718135c375SStefano Zampini   /* mark owned vertices */
7723e6c54aaSStefano Zampini   vown = NULL;
7733e6c54aaSStefano Zampini   if (pown) {
7743e6c54aaSStefano Zampini     ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
7758135c375SStefano Zampini     for (p=cStart;p<cEnd;p++) {
7768135c375SStefano Zampini       PetscInt i,closureSize,*closure = NULL;
7778135c375SStefano Zampini 
7783e6c54aaSStefano Zampini       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
7798135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
7808135c375SStefano Zampini       for (i=0;i<closureSize;i++) {
7818135c375SStefano Zampini         const PetscInt pp = closure[2*i];
7828135c375SStefano Zampini 
7838135c375SStefano Zampini         if (pp >= vStart && pp < vEnd) {
7843e6c54aaSStefano Zampini           ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr);
7858135c375SStefano Zampini         }
7868135c375SStefano Zampini       }
7878135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
7888135c375SStefano Zampini     }
7898135c375SStefano Zampini   }
7908135c375SStefano Zampini 
7918135c375SStefano Zampini   /* vertex_parents (Non-conforming meshes) */
7928135c375SStefano Zampini   parentSection  = NULL;
7938135c375SStefano Zampini   if (enable_ncmesh) {
7948135c375SStefano Zampini     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
7958135c375SStefano Zampini   }
7968135c375SStefano Zampini   if (parentSection) {
7978135c375SStefano Zampini     PetscInt vp,gvp;
7988135c375SStefano Zampini 
7998135c375SStefano Zampini     for (vp=0,p=vStart;p<vEnd;p++) {
8008135c375SStefano Zampini       DMLabel  dlabel;
8018135c375SStefano Zampini       PetscInt parent,depth;
8028135c375SStefano Zampini 
8033e6c54aaSStefano Zampini       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
8048135c375SStefano Zampini       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
8058135c375SStefano Zampini       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
8068135c375SStefano Zampini       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
8078135c375SStefano Zampini       if (parent != p) vp++;
8088135c375SStefano Zampini     }
8098135c375SStefano Zampini     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
8108135c375SStefano Zampini     if (gvp) {
8118135c375SStefano Zampini       PetscInt  maxsupp;
8128135c375SStefano Zampini       PetscBool *skip = NULL;
8138135c375SStefano Zampini 
8148135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
8158135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
8168135c375SStefano Zampini       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
8178135c375SStefano Zampini       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
8188135c375SStefano Zampini       for (p=vStart;p<vEnd;p++) {
8198135c375SStefano Zampini         DMLabel  dlabel;
8208135c375SStefano Zampini         PetscInt parent;
8218135c375SStefano Zampini 
8223e6c54aaSStefano Zampini         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
8238135c375SStefano Zampini         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
8248135c375SStefano Zampini         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
8258135c375SStefano Zampini         if (parent != p) {
826d42aff11SStefano Zampini           int            vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
82711a4995dSStefano Zampini           PetscInt       i,nv,size,n,numChildren,depth = -1;
8288135c375SStefano Zampini           const PetscInt *children;
8298135c375SStefano Zampini           ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr);
8308135c375SStefano Zampini           switch (size) {
8318135c375SStefano Zampini             case 2: /* edge */
8328135c375SStefano Zampini               nv   = 0;
833cc0d3ed7SStefano Zampini               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
83477eacf09SStefano Zampini               ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
8358135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
8368135c375SStefano Zampini               for (i=0;i<nv;i++) {
83711a4995dSStefano Zampini                 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr);
8388135c375SStefano Zampini               }
8398135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
8408135c375SStefano Zampini               vp--;
8418135c375SStefano Zampini               break;
8428135c375SStefano Zampini             case 4: /* face */
8438135c375SStefano Zampini               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
8448135c375SStefano Zampini               for (n=0;n<numChildren;n++) {
8458135c375SStefano Zampini                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
8468135c375SStefano Zampini                 if (!depth) {
8478135c375SStefano Zampini                   const PetscInt *hvsupp,*hesupp,*cone;
8488135c375SStefano Zampini                   PetscInt       hvsuppSize,hesuppSize,coneSize;
849451a39c7SStefano Zampini                   PetscInt       hv = children[n],he = -1,f;
8508135c375SStefano Zampini 
851451a39c7SStefano Zampini                   ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
8528135c375SStefano Zampini                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
8538135c375SStefano Zampini                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
8548135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
8558135c375SStefano Zampini                     PetscInt ep;
8568135c375SStefano Zampini                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
8578135c375SStefano Zampini                     if (ep != hvsupp[i]) {
8588135c375SStefano Zampini                       he = hvsupp[i];
8598135c375SStefano Zampini                     } else {
8608135c375SStefano Zampini                       skip[i] = PETSC_TRUE;
8618135c375SStefano Zampini                     }
8628135c375SStefano Zampini                   }
863451a39c7SStefano Zampini                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
8648135c375SStefano Zampini                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
86511a4995dSStefano Zampini                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
8668135c375SStefano Zampini                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
8678135c375SStefano Zampini                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
8688135c375SStefano Zampini                   for (f=0;f<hesuppSize;f++) {
8698135c375SStefano Zampini                     PetscInt j;
8708135c375SStefano Zampini 
8718135c375SStefano Zampini                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
8728135c375SStefano Zampini                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
8738135c375SStefano Zampini                     for (j=0;j<coneSize;j++) {
8748135c375SStefano Zampini                       PetscInt k;
8758135c375SStefano Zampini                       for (k=0;k<hvsuppSize;k++) {
8768135c375SStefano Zampini                         if (hvsupp[k] == cone[j]) {
8778135c375SStefano Zampini                           skip[k] = PETSC_TRUE;
8788135c375SStefano Zampini                           break;
8798135c375SStefano Zampini                         }
8808135c375SStefano Zampini                       }
8818135c375SStefano Zampini                     }
8828135c375SStefano Zampini                   }
8838135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
8848135c375SStefano Zampini                     if (!skip[i]) {
8858135c375SStefano Zampini                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
88611a4995dSStefano Zampini                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
8878135c375SStefano Zampini                     }
8888135c375SStefano Zampini                   }
8898135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
8908135c375SStefano Zampini                   for (i=0;i<2;i++) {
89111a4995dSStefano Zampini                     ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);CHKERRQ(ierr);
8928135c375SStefano Zampini                   }
8938135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
8948135c375SStefano Zampini                   vp--;
8958135c375SStefano Zampini                 }
8968135c375SStefano Zampini               }
8978135c375SStefano Zampini               break;
8988135c375SStefano Zampini             default:
8998135c375SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size);
9008135c375SStefano Zampini           }
9018135c375SStefano Zampini         }
9028135c375SStefano Zampini       }
9038135c375SStefano Zampini       ierr = PetscFree(skip);CHKERRQ(ierr);
9048135c375SStefano Zampini     }
9058135c375SStefano Zampini     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
9068135c375SStefano Zampini   }
9078135c375SStefano Zampini   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
9083e6c54aaSStefano Zampini   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
9098135c375SStefano Zampini 
9108135c375SStefano Zampini   /* vertices */
91177eacf09SStefano Zampini   if (hovec) { /* higher-order meshes */
91277eacf09SStefano Zampini     const char *fec;
91377eacf09SStefano Zampini     PetscInt   i,n;
91477eacf09SStefano Zampini 
91577eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
91677eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
91777eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
91877eacf09SStefano Zampini     ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr);
91977eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
9203e7cab22SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr);
92177eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
92277eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
92377eacf09SStefano Zampini     ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr);
92477eacf09SStefano Zampini     ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr);
92577eacf09SStefano Zampini     if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
92677eacf09SStefano Zampini     for (i=0;i<n/sdim;i++) {
92777eacf09SStefano Zampini       PetscInt s;
92877eacf09SStefano Zampini       for (s=0;s<sdim;s++) {
92977eacf09SStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr);
93077eacf09SStefano Zampini       }
93177eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
93277eacf09SStefano Zampini     }
93377eacf09SStefano Zampini     ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr);
93477eacf09SStefano Zampini   } else {
935cc0d3ed7SStefano Zampini     ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr);
9368135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
937cc0d3ed7SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr);
9388135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
9398135c375SStefano Zampini     ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
940cc0d3ed7SStefano Zampini     for (p=0;p<nvert/sdim;p++) {
941cc0d3ed7SStefano Zampini       PetscInt s;
942cc0d3ed7SStefano Zampini       for (s=0;s<sdim;s++) {
94377eacf09SStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr);
9448135c375SStefano Zampini       }
945cc0d3ed7SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9468135c375SStefano Zampini     }
9478135c375SStefano Zampini     ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
94877eacf09SStefano Zampini   }
9498135c375SStefano Zampini   PetscFunctionReturn(0);
9508135c375SStefano Zampini }
9518135c375SStefano Zampini 
9528135c375SStefano Zampini /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
9538135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
9548135c375SStefano Zampini {
9558135c375SStefano Zampini   PetscErrorCode ierr;
9568135c375SStefano Zampini   PetscBool      isglvis,isascii;
9578135c375SStefano Zampini 
9588135c375SStefano Zampini   PetscFunctionBegin;
9598135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9608135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
9618135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
9628135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
9638135c375SStefano Zampini   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
9648135c375SStefano Zampini   if (isglvis) {
9658135c375SStefano Zampini     PetscViewer          view;
9668135c375SStefano Zampini     PetscViewerGLVisType type;
9678135c375SStefano Zampini 
9688135c375SStefano Zampini     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
9698135c375SStefano Zampini     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
9708135c375SStefano Zampini     if (view) { /* in the socket case, it may happen that the connection failed */
9718135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
9728135c375SStefano Zampini         PetscMPIInt size,rank;
9738135c375SStefano Zampini         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
9748135c375SStefano Zampini         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
9758135c375SStefano Zampini         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
9768135c375SStefano Zampini       }
9778135c375SStefano Zampini       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
9788135c375SStefano Zampini       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
9798135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
9808135c375SStefano Zampini         PetscInt    dim;
9818135c375SStefano Zampini         const char* name;
9828135c375SStefano Zampini 
9838135c375SStefano Zampini         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
9848135c375SStefano Zampini         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
9858135c375SStefano Zampini         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
9868135c375SStefano Zampini         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
9878135c375SStefano Zampini       }
9888135c375SStefano Zampini     }
9898135c375SStefano Zampini   } else {
9908135c375SStefano Zampini     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
9918135c375SStefano Zampini   }
9928135c375SStefano Zampini   PetscFunctionReturn(0);
9938135c375SStefano Zampini }
994