xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 8183e3f6d727e23efcbafaec026bca9e6938d73d)
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 {
111bfadf1d8SStefano 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) {
136bfadf1d8SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr);
1378135c375SStefano Zampini           } else {
1388135c375SStefano Zampini             H1   = PETSC_FALSE;
139bfadf1d8SStefano 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);
177bfadf1d8SStefano 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];
180bfadf1d8SStefano 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 
234044a5661SStefano Zampini MFEM_cid mfem_table_cid_unint[4][9] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF},
235044a5661SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF},
236044a5661SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF},
237044a5661SStefano Zampini                                         {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_UNDEF,MFEM_UNDEF,MFEM_CUBE } };
238044a5661SStefano Zampini 
239cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid)
2408135c375SStefano Zampini {
2418135c375SStefano Zampini   DMLabel        dlabel;
242044a5661SStefano Zampini   PetscInt       depth,csize,pdepth,dim;
2438135c375SStefano Zampini   PetscErrorCode ierr;
2448135c375SStefano Zampini 
2458135c375SStefano Zampini   PetscFunctionBegin;
2468135c375SStefano Zampini   ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
247044a5661SStefano Zampini   ierr = DMLabelGetValue(dlabel,p,&pdepth);CHKERRQ(ierr);
2488135c375SStefano Zampini   ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr);
249044a5661SStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
250044a5661SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
2518135c375SStefano Zampini   if (label) {
2528135c375SStefano Zampini     ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr);
25377eacf09SStefano Zampini   } else *mid = 1;
254044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
255044a5661SStefano Zampini #if defined PETSC_USE_DEBUG
256044a5661SStefano Zampini     if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim);
257044a5661SStefano Zampini     if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p);
258044a5661SStefano 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);
259044a5661SStefano Zampini #endif
260044a5661SStefano Zampini     *cid = mfem_table_cid_unint[dim][csize];
261044a5661SStefano Zampini   } else {
262044a5661SStefano Zampini #if defined PETSC_USE_DEBUG
263044a5661SStefano Zampini     if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p);
264044a5661SStefano Zampini     if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p);
265044a5661SStefano Zampini #endif
266044a5661SStefano Zampini     *cid = mfem_table_cid[pdepth][csize];
267044a5661SStefano Zampini   }
2688135c375SStefano Zampini   PetscFunctionReturn(0);
2698135c375SStefano Zampini }
2708135c375SStefano Zampini 
271cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int 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) {
281cc0d3ed7SStefano Zampini     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
282cc0d3ed7SStefano Zampini     ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr);
283cc0d3ed7SStefano Zampini     off  = off/sdim;
284cc0d3ed7SStefano Zampini     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
285cc0d3ed7SStefano Zampini   }
286cc0d3ed7SStefano Zampini   if (!dof) {
2878135c375SStefano Zampini     ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
2888135c375SStefano Zampini     for (i=0,q=0;i<numPoints*2;i+= 2)
2898135c375SStefano Zampini       if ((points[i] >= vStart) && (points[i] < vEnd))
290*8183e3f6SStefano Zampini         vids[q++] = (int)(points[i]-vStart+off);
2918135c375SStefano Zampini     ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
292cc0d3ed7SStefano Zampini   } else {
293cc0d3ed7SStefano Zampini     ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr);
294cc0d3ed7SStefano Zampini     ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr);
295*8183e3f6SStefano Zampini     for (q=0;q<dof/sdim;q++) vids[q] = (int)(off/sdim + q);
296cc0d3ed7SStefano Zampini   }
2978135c375SStefano Zampini   *nv = q;
2988135c375SStefano Zampini   PetscFunctionReturn(0);
2998135c375SStefano Zampini }
3008135c375SStefano Zampini 
30177eacf09SStefano Zampini /*
30277eacf09SStefano Zampini    ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR
30377eacf09SStefano Zampini    Higher order meshes are also supported
30477eacf09SStefano Zampini */
3058135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
3068135c375SStefano Zampini {
3078135c375SStefano Zampini   DMLabel              label;
3088135c375SStefano Zampini   PetscSection         coordSection,parentSection;
30977eacf09SStefano Zampini   Vec                  coordinates,hovec;
3108135c375SStefano Zampini   const PetscScalar    *array;
3118135c375SStefano Zampini   PetscInt             bf,p,sdim,dim,depth,novl;
312cc0d3ed7SStefano Zampini   PetscInt             cStart,cEnd,cEndInterior,vStart,vEnd,nvert;
3133924b612SStefano Zampini   PetscMPIInt          size;
3143e6c54aaSStefano Zampini   PetscBool            localized,isascii;
3153e6c54aaSStefano Zampini   PetscBool            enable_mfem,enable_boundary,enable_ncmesh;
3163e6c54aaSStefano Zampini   PetscBT              pown,vown;
3178135c375SStefano Zampini   PetscErrorCode       ierr;
3188135c375SStefano Zampini   PetscContainer       glvis_container;
319044a5661SStefano Zampini   PetscBool            cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE;
32077eacf09SStefano Zampini   const char           *fmt;
3217bf4dd16SStefano Zampini   char                 emark[64] = "",bmark[64] = "";
3228135c375SStefano Zampini 
3238135c375SStefano Zampini   PetscFunctionBegin;
3248135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3258135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3268135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
3278135c375SStefano Zampini   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
3283924b612SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
3293924b612SStefano Zampini   if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
3308135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
3318135c375SStefano Zampini 
3328135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
3338135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
3348135c375SStefano Zampini   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
3358135c375SStefano Zampini   {
3368135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
3378135c375SStefano Zampini     ierr    = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
3388135c375SStefano Zampini     enabled = glvis_info->enabled;
33977eacf09SStefano Zampini     fmt     = glvis_info->fmt;
3408135c375SStefano Zampini   }
341cc0d3ed7SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3428135c375SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
3438135c375SStefano Zampini   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
3448135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
34577eacf09SStefano Zampini   ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr);
346cc0d3ed7SStefano Zampini   ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr);
34777eacf09SStefano Zampini   if (periodic && !localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized");
348cc0d3ed7SStefano Zampini   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
34977eacf09SStefano Zampini   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
35077eacf09SStefano Zampini   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
35143eeeb2dSStefano Zampini   if (!coordinates) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector");
35277eacf09SStefano Zampini 
35377eacf09SStefano Zampini   /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh
35477eacf09SStefano Zampini      DMPlex does not currently support HO meshes, so there's no API for this */
35577eacf09SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr);
3568135c375SStefano Zampini 
3578135c375SStefano Zampini   /*
3588135c375SStefano Zampini      a couple of sections of the mesh specification are disabled
35977eacf09SStefano Zampini        - boundary: unless we want to visualize boundary attributes or we have a periodic mesh
36077eacf09SStefano Zampini                    the boundary is not needed for proper mesh visualization
36177eacf09SStefano Zampini        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package
3623e6c54aaSStefano Zampini                          and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel)
3638135c375SStefano Zampini   */
364e74666bcSStefano Zampini   enable_boundary = periodic;
3658135c375SStefano Zampini   enable_ncmesh   = PETSC_FALSE;
3663e6c54aaSStefano Zampini   enable_mfem     = PETSC_FALSE;
3677bf4dd16SStefano Zampini   /* I'm tired of problems with negative values in the markers, disable them */
3688135c375SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
3693e6c54aaSStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr);
3703e6c54aaSStefano 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);
3713e6c54aaSStefano 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);
3727bf4dd16SStefano Zampini   ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),NULL);CHKERRQ(ierr);
3737bf4dd16SStefano Zampini   ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),NULL);CHKERRQ(ierr);
3748135c375SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
3753924b612SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
3763924b612SStefano Zampini   if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel");
3777e1aca4eSStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
3787e1aca4eSStefano Zampini   if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
3797e1aca4eSStefano Zampini                                                              "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0");
3807e1aca4eSStefano Zampini   if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. "
3817e1aca4eSStefano Zampini                                                            "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0");
382044a5661SStefano Zampini   if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */
383044a5661SStefano Zampini     if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth);
384044a5661SStefano Zampini     cellvertex = PETSC_TRUE;
385044a5661SStefano Zampini   }
3868135c375SStefano Zampini 
3878135c375SStefano Zampini   /* Identify possible cells in the overlap */
3888135c375SStefano Zampini   novl = 0;
3898135c375SStefano Zampini   pown = NULL;
3903924b612SStefano Zampini   if (size > 1) {
3913e6c54aaSStefano Zampini     IS             globalNum = NULL;
3923e6c54aaSStefano Zampini     const PetscInt *gNum;
3933e6c54aaSStefano Zampini     PetscBool      ovl  = PETSC_FALSE;
3943e6c54aaSStefano Zampini 
3958135c375SStefano Zampini     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
3968135c375SStefano Zampini     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
3978135c375SStefano Zampini     for (p=cStart; p<cEnd; p++) {
3988135c375SStefano Zampini       if (gNum[p-cStart] < 0) {
3998135c375SStefano Zampini         ovl = PETSC_TRUE;
4008135c375SStefano Zampini         novl++;
4018135c375SStefano Zampini       }
4028135c375SStefano Zampini     }
4038135c375SStefano Zampini     if (ovl) {
4048135c375SStefano Zampini       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
4058135c375SStefano Zampini          TODO: garbage collector? attach pown to dm?  */
4063e6c54aaSStefano Zampini       ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr);
4073e6c54aaSStefano Zampini       for (p=cStart; p<cEnd; p++) {
4083e6c54aaSStefano Zampini         if (gNum[p-cStart] < 0) continue;
4093e6c54aaSStefano Zampini         else {
4103e6c54aaSStefano Zampini           ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
4118135c375SStefano Zampini         }
4128135c375SStefano Zampini       }
4133e6c54aaSStefano Zampini     }
4143e6c54aaSStefano Zampini     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
4153e6c54aaSStefano Zampini   }
4168135c375SStefano Zampini 
4173e6c54aaSStefano Zampini   /* return if this process is disabled */
4188135c375SStefano Zampini   if (!enabled) {
4198135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
4208135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
4218135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
4228135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
4238135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
4248135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
4258135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
4268135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
4278135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
4288135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
4298135c375SStefano Zampini     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
4308135c375SStefano Zampini     PetscFunctionReturn(0);
4318135c375SStefano Zampini   }
4328135c375SStefano Zampini 
4333e6c54aaSStefano Zampini   if (enable_mfem) {
4343e6c54aaSStefano Zampini     if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */
4353e6c54aaSStefano Zampini       PetscInt    vpc = 0;
4363e6c54aaSStefano Zampini       char        fec[64];
4373e6c54aaSStefano Zampini       int         vids[8] = {0,1,2,3,4,5,6,7};
438044a5661SStefano Zampini       int         hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3};
439044a5661SStefano Zampini       int         quadv[8] = {0,1,3,2}, triv[3] = {0,1,2};
440044a5661SStefano Zampini       int         *dof = NULL;
4413e6c54aaSStefano Zampini       PetscScalar *array,*ptr;
4423e6c54aaSStefano Zampini 
443bfadf1d8SStefano Zampini       ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr);
4443e6c54aaSStefano Zampini       if (cEnd-cStart) {
4453e6c54aaSStefano Zampini         PetscInt fpc;
4463e6c54aaSStefano Zampini 
4473e6c54aaSStefano Zampini         ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
4483e6c54aaSStefano Zampini         switch(dim) {
4493e6c54aaSStefano Zampini           case 1:
4503e6c54aaSStefano Zampini             vpc = 2;
4513e6c54aaSStefano Zampini             dof = hexv;
4523e6c54aaSStefano Zampini             break;
4533e6c54aaSStefano Zampini           case 2:
4543e6c54aaSStefano Zampini             switch (fpc) {
4553e6c54aaSStefano Zampini               case 3:
4563e6c54aaSStefano Zampini                 vpc = 3;
457044a5661SStefano Zampini                 dof = triv;
4583e6c54aaSStefano Zampini                 break;
4593e6c54aaSStefano Zampini               case 4:
4603e6c54aaSStefano Zampini                 vpc = 4;
461044a5661SStefano Zampini                 dof = quadv;
4623e6c54aaSStefano Zampini                 break;
4633e6c54aaSStefano Zampini               default:
4643e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4653e6c54aaSStefano Zampini                 break;
4663e6c54aaSStefano Zampini             }
4673e6c54aaSStefano Zampini             break;
4683e6c54aaSStefano Zampini           case 3:
4693e6c54aaSStefano Zampini             switch (fpc) {
470044a5661SStefano Zampini               case 4: /* TODO: still need to understand L2 ordering for tets */
4713e6c54aaSStefano Zampini                 vpc = 4;
472044a5661SStefano Zampini                 dof = tetv;
473044a5661SStefano Zampini                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case");
4743e6c54aaSStefano Zampini                 break;
4753e6c54aaSStefano Zampini               case 6:
476044a5661SStefano Zampini                 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc);
477044a5661SStefano Zampini                 vpc = 8;
478044a5661SStefano Zampini                 dof = hexv;
479044a5661SStefano Zampini                 break;
480044a5661SStefano Zampini               case 8:
481044a5661SStefano Zampini                 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4823e6c54aaSStefano Zampini                 vpc = 8;
4833e6c54aaSStefano Zampini                 dof = hexv;
4843e6c54aaSStefano Zampini                 break;
4853e6c54aaSStefano Zampini               default:
4863e6c54aaSStefano Zampini                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
4873e6c54aaSStefano Zampini                 break;
4883e6c54aaSStefano Zampini             }
4893e6c54aaSStefano Zampini             break;
4903e6c54aaSStefano Zampini           default:
4913e6c54aaSStefano Zampini             SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
4923e6c54aaSStefano Zampini             break;
4933e6c54aaSStefano Zampini         }
4943e6c54aaSStefano Zampini         ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr);
4953e6c54aaSStefano Zampini       }
496044a5661SStefano Zampini       if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs");
4973e6c54aaSStefano Zampini       ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr);
4983e6c54aaSStefano Zampini       ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr);
4993e6c54aaSStefano Zampini       ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr);
5003e6c54aaSStefano Zampini       ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr);
5013e6c54aaSStefano Zampini       ierr = VecGetArray(hovec,&array);CHKERRQ(ierr);
5023e6c54aaSStefano Zampini       ptr  = array;
5033e6c54aaSStefano Zampini       for (p=cStart;p<cEnd;p++) {
5043e6c54aaSStefano Zampini         PetscInt    csize,v,d;
5053e6c54aaSStefano Zampini         PetscScalar *vals = NULL;
5063e6c54aaSStefano Zampini 
5073e6c54aaSStefano Zampini         if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
5083e6c54aaSStefano Zampini         ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
5093e6c54aaSStefano 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);
5103e6c54aaSStefano Zampini         for (v=0;v<vpc;v++) {
5113e6c54aaSStefano Zampini           for (d=0;d<sdim;d++) {
5123e6c54aaSStefano Zampini             ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d];
5133e6c54aaSStefano Zampini           }
5143e6c54aaSStefano Zampini         }
5153e6c54aaSStefano Zampini         ptr += vpc*sdim;
5163e6c54aaSStefano Zampini         ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
5173e6c54aaSStefano Zampini       }
5183e6c54aaSStefano Zampini       ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr);
5193e6c54aaSStefano Zampini     }
5203e6c54aaSStefano Zampini   }
5213e6c54aaSStefano Zampini 
5223e6c54aaSStefano Zampini 
5238135c375SStefano Zampini   /* header */
5248135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
5258135c375SStefano Zampini 
5268135c375SStefano Zampini   /* topological dimension */
5278135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
5288135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
5298135c375SStefano Zampini 
5308135c375SStefano Zampini   /* elements */
5317bf4dd16SStefano Zampini   ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr);
5328135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
5338135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
5348135c375SStefano Zampini   for (p=cStart;p<cEnd;p++) {
53511a4995dSStefano Zampini     int      vids[8];
53611a4995dSStefano Zampini     PetscInt i,nv = 0,cid = -1,mid = 1;
5378135c375SStefano Zampini 
5383e6c54aaSStefano Zampini     if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue;
5398135c375SStefano Zampini     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
54077eacf09SStefano Zampini     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
54177eacf09SStefano Zampini     ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
5428135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
5438135c375SStefano Zampini     for (i=0;i<nv;i++) {
544bfadf1d8SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr);
5458135c375SStefano Zampini     }
5468135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
5478135c375SStefano Zampini   }
5488135c375SStefano Zampini 
549cc0d3ed7SStefano Zampini   /* boundary */
550cc0d3ed7SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
551cc0d3ed7SStefano Zampini   if (!enable_boundary) {
552cc0d3ed7SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
553cc0d3ed7SStefano Zampini   } else {
55477eacf09SStefano Zampini     DMLabel  perLabel;
55577eacf09SStefano Zampini     PetscBT  bfaces;
55677eacf09SStefano Zampini     PetscInt fStart,fEnd,fEndInterior,*fcells;
55777eacf09SStefano Zampini     PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0;
558cc0d3ed7SStefano Zampini     PetscInt fv1[]     = {0,1},
559cc0d3ed7SStefano Zampini              fv2tri[]  = {0,1,
560cc0d3ed7SStefano Zampini                           1,2,
561cc0d3ed7SStefano Zampini                           2,0},
562cc0d3ed7SStefano Zampini              fv2quad[] = {0,1,
563cc0d3ed7SStefano Zampini                           1,2,
564cc0d3ed7SStefano Zampini                           2,3,
565cc0d3ed7SStefano Zampini                           3,0},
566cc0d3ed7SStefano Zampini              fv3tet[]  = {0,1,2,
567cc0d3ed7SStefano Zampini                           0,3,1,
568cc0d3ed7SStefano Zampini                           0,2,3,
569cc0d3ed7SStefano Zampini                           2,1,3},
570cc0d3ed7SStefano Zampini              fv3hex[]  = {0,1,2,3,
571cc0d3ed7SStefano Zampini                        4,5,6,7,
572cc0d3ed7SStefano Zampini                        0,3,5,4,
573cc0d3ed7SStefano Zampini                        2,1,7,6,
574cc0d3ed7SStefano Zampini                        3,2,6,5,
575cc0d3ed7SStefano Zampini                        0,4,7,1};
576cc0d3ed7SStefano Zampini 
5778135c375SStefano Zampini     /* determine orientation of boundary mesh */
5788135c375SStefano Zampini     if (cEnd-cStart) {
5798135c375SStefano Zampini       ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
5808135c375SStefano Zampini       switch(dim) {
5818135c375SStefano Zampini         case 1:
5828135c375SStefano Zampini           if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
5838135c375SStefano Zampini           faces = fv1;
5848135c375SStefano Zampini           vpf = 1;
58577eacf09SStefano Zampini           vpc = 2;
5868135c375SStefano Zampini           break;
5878135c375SStefano Zampini         case 2:
5888135c375SStefano Zampini           switch (fpc) {
5898135c375SStefano Zampini             case 3:
5908135c375SStefano Zampini               faces = fv2tri;
5918135c375SStefano Zampini               vpf   = 2;
59277eacf09SStefano Zampini               vpc   = 3;
5938135c375SStefano Zampini               break;
5948135c375SStefano Zampini             case 4:
5958135c375SStefano Zampini               faces = fv2quad;
5968135c375SStefano Zampini               vpf   = 2;
59777eacf09SStefano Zampini               vpc   = 4;
5988135c375SStefano Zampini               break;
5998135c375SStefano Zampini             default:
6008135c375SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
6018135c375SStefano Zampini               break;
6028135c375SStefano Zampini           }
6038135c375SStefano Zampini           break;
6048135c375SStefano Zampini         case 3:
6058135c375SStefano Zampini           switch (fpc) {
6068135c375SStefano Zampini             case 4:
6078135c375SStefano Zampini               faces = fv3tet;
6088135c375SStefano Zampini               vpf   = 3;
60977eacf09SStefano Zampini               vpc   = 4;
6108135c375SStefano Zampini               break;
6118135c375SStefano Zampini             case 6:
6128135c375SStefano Zampini               faces = fv3hex;
6138135c375SStefano Zampini               vpf   = 4;
61477eacf09SStefano Zampini               vpc   = 8;
6158135c375SStefano Zampini               break;
6168135c375SStefano Zampini             default:
6178135c375SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
6188135c375SStefano Zampini               break;
6198135c375SStefano Zampini           }
6208135c375SStefano Zampini           break;
6218135c375SStefano Zampini         default:
6228135c375SStefano Zampini           SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
6238135c375SStefano Zampini           break;
6248135c375SStefano Zampini       }
6258135c375SStefano Zampini     }
626cc0d3ed7SStefano Zampini     ierr = DMPlexGetHybridBounds(dm,NULL,&fEndInterior,NULL,NULL);CHKERRQ(ierr);
6278135c375SStefano Zampini     ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
6288135c375SStefano Zampini     fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
62977eacf09SStefano Zampini     ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr);
63077eacf09SStefano Zampini     ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr);
63177eacf09SStefano Zampini     ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr);
63277eacf09SStefano Zampini     ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
63377eacf09SStefano Zampini     if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */
63477eacf09SStefano Zampini       ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr);
63577eacf09SStefano Zampini       ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr);
63677eacf09SStefano Zampini       ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr);
63777eacf09SStefano Zampini       for (p=cStart;p<cEnd;p++) {
63877eacf09SStefano Zampini         PetscInt dof;
63977eacf09SStefano Zampini         ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
64077eacf09SStefano Zampini         if (dof) {
64177eacf09SStefano Zampini           PetscInt    v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL;
64277eacf09SStefano Zampini           PetscScalar *vals = NULL;
64377eacf09SStefano Zampini 
64477eacf09SStefano Zampini           if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim);
6457bf4dd16SStefano 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);
64677eacf09SStefano Zampini           ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
64777eacf09SStefano Zampini           ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
64877eacf09SStefano Zampini           for (v=0;v<cellClosureSize;v++)
64977eacf09SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) {
65077eacf09SStefano Zampini               vidxs = cellClosure + 2*v;
65177eacf09SStefano Zampini               break;
65277eacf09SStefano Zampini             }
65377eacf09SStefano Zampini           if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices");
65477eacf09SStefano Zampini           for (v=0;v<vpc;v++) {
65577eacf09SStefano Zampini             PetscInt s;
656044a5661SStefano Zampini 
65777eacf09SStefano Zampini             for (s=0;s<sdim;s++) {
65877eacf09SStefano Zampini               if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+vpc*sdim])>PETSC_MACHINE_EPSILON) {
65977eacf09SStefano Zampini                 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr);
66077eacf09SStefano Zampini               }
66177eacf09SStefano Zampini             }
66277eacf09SStefano Zampini           }
66377eacf09SStefano Zampini           ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
66477eacf09SStefano Zampini           ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr);
66577eacf09SStefano Zampini         }
66677eacf09SStefano Zampini       }
66777eacf09SStefano Zampini       if (dim > 1) {
66877eacf09SStefano Zampini         PetscInt eEnd,eStart,eEndInterior;
669044a5661SStefano Zampini 
67077eacf09SStefano Zampini         ierr = DMPlexGetHybridBounds(dm,NULL,NULL,&eEndInterior,NULL);CHKERRQ(ierr);
67177eacf09SStefano Zampini         ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr);
67277eacf09SStefano Zampini         eEnd = eEndInterior < 0 ? eEnd : eEndInterior;
67377eacf09SStefano Zampini         for (p=eStart;p<eEnd;p++) {
67477eacf09SStefano Zampini           const PetscInt *cone;
67577eacf09SStefano Zampini           PetscInt       coneSize,i;
67677eacf09SStefano Zampini           PetscBool      ispe = PETSC_TRUE;
67777eacf09SStefano Zampini 
67877eacf09SStefano Zampini           ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
67977eacf09SStefano Zampini           ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
68077eacf09SStefano Zampini           for (i=0;i<coneSize;i++) {
68177eacf09SStefano Zampini             PetscInt v;
68277eacf09SStefano Zampini 
68377eacf09SStefano Zampini             ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
68477eacf09SStefano Zampini             ispe = (PetscBool)(ispe && (v==2));
68577eacf09SStefano Zampini           }
68677eacf09SStefano Zampini           if (ispe && coneSize) {
68777eacf09SStefano Zampini             ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
68877eacf09SStefano Zampini           }
68977eacf09SStefano Zampini         }
69077eacf09SStefano Zampini         if (dim > 2) {
69177eacf09SStefano Zampini           for (p=fStart;p<fEnd;p++) {
69277eacf09SStefano Zampini             const PetscInt *cone;
69377eacf09SStefano Zampini             PetscInt       coneSize,i;
69477eacf09SStefano Zampini             PetscBool      ispe = PETSC_TRUE;
69577eacf09SStefano Zampini 
69677eacf09SStefano Zampini             ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
69777eacf09SStefano Zampini             ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
69877eacf09SStefano Zampini             for (i=0;i<coneSize;i++) {
69977eacf09SStefano Zampini               PetscInt v;
70077eacf09SStefano Zampini 
70177eacf09SStefano Zampini               ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
70277eacf09SStefano Zampini               ispe = (PetscBool)(ispe && (v==2));
70377eacf09SStefano Zampini             }
70477eacf09SStefano Zampini             if (ispe && coneSize) {
70577eacf09SStefano Zampini               ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr);
70677eacf09SStefano Zampini             }
70777eacf09SStefano Zampini           }
70877eacf09SStefano Zampini         }
70977eacf09SStefano Zampini       }
71077eacf09SStefano Zampini     }
71177eacf09SStefano Zampini     for (p=fStart;p<fEnd;p++) {
71277eacf09SStefano Zampini       const PetscInt *support;
7138135c375SStefano Zampini       PetscInt       supportSize;
71477eacf09SStefano Zampini       PetscBool      isbf = PETSC_FALSE;
7158135c375SStefano Zampini 
7168135c375SStefano Zampini       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
7173e6c54aaSStefano Zampini       if (pown) {
7188135c375SStefano Zampini         PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
71977eacf09SStefano Zampini         PetscInt  i;
72077eacf09SStefano Zampini 
72177eacf09SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
72277eacf09SStefano Zampini         for (i=0;i<supportSize;i++) {
72377eacf09SStefano Zampini           if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE;
72477eacf09SStefano Zampini           else has_ghost = PETSC_TRUE;
72577eacf09SStefano Zampini         }
72677eacf09SStefano Zampini         isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost));
72777eacf09SStefano Zampini       } else {
72877eacf09SStefano Zampini         isbf = (PetscBool)(supportSize == 1);
72977eacf09SStefano Zampini       }
73077eacf09SStefano Zampini       if (!isbf && perLabel) {
73177eacf09SStefano Zampini         const PetscInt *cone;
73277eacf09SStefano Zampini         PetscInt       coneSize,i;
73377eacf09SStefano Zampini 
73477eacf09SStefano Zampini         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
73577eacf09SStefano Zampini         ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr);
73677eacf09SStefano Zampini         isbf = PETSC_TRUE;
73777eacf09SStefano Zampini         for (i=0;i<coneSize;i++) {
73877eacf09SStefano Zampini           PetscInt v,d;
73977eacf09SStefano Zampini 
74077eacf09SStefano Zampini           ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr);
74177eacf09SStefano Zampini           ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr);
74277eacf09SStefano Zampini           isbf = (PetscBool)(isbf && v != d);
74377eacf09SStefano Zampini         }
74477eacf09SStefano Zampini       }
74577eacf09SStefano Zampini       if (isbf) {
74677eacf09SStefano Zampini         ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr);
74777eacf09SStefano Zampini       }
74877eacf09SStefano Zampini     }
74977eacf09SStefano Zampini     /* count boundary faces */
75077eacf09SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
75177eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
75277eacf09SStefano Zampini         const PetscInt *support;
75377eacf09SStefano Zampini         PetscInt       supportSize,c;
7548135c375SStefano Zampini 
7558135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
7568135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
7573e6c54aaSStefano Zampini         if (pown) {
75877eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
75977eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
76077eacf09SStefano Zampini               bf++;
7618135c375SStefano Zampini             }
76277eacf09SStefano Zampini           }
76377eacf09SStefano Zampini         } else bf += supportSize;
7648135c375SStefano Zampini       }
7658135c375SStefano Zampini     }
7667bf4dd16SStefano Zampini     ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr);
7678135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
7688135c375SStefano Zampini     for (p=fStart;p<fEnd;p++) {
76977eacf09SStefano Zampini       if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) {
7708135c375SStefano Zampini         const PetscInt *support;
77177eacf09SStefano Zampini         PetscInt       supportSize,c,nc = 0;
7728135c375SStefano Zampini 
7738135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
7748135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
7753e6c54aaSStefano Zampini         if (pown) {
77677eacf09SStefano Zampini           for (c=0;c<supportSize;c++) {
77777eacf09SStefano Zampini             if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) {
77877eacf09SStefano Zampini               fcells[nc++] = support[c];
7798135c375SStefano Zampini             }
78077eacf09SStefano Zampini           }
78177eacf09SStefano Zampini         } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c];
78277eacf09SStefano Zampini         for (c=0;c<nc;c++) {
78377eacf09SStefano Zampini           const    PetscInt *cone;
78477eacf09SStefano Zampini           int      vids[8];
78577eacf09SStefano Zampini           PetscInt i,cell,cl,nv,cid = -1,mid = -1;
7868135c375SStefano Zampini 
78777eacf09SStefano Zampini           cell = fcells[c];
78877eacf09SStefano Zampini           ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr);
78977eacf09SStefano Zampini           for (cl=0;cl<fpc;cl++)
79077eacf09SStefano Zampini             if (cone[cl] == p)
7918135c375SStefano Zampini               break;
7928135c375SStefano Zampini 
79377eacf09SStefano Zampini           /* face material id and type */
7948135c375SStefano Zampini           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
7958135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
79677eacf09SStefano Zampini           /* vertex ids */
79777eacf09SStefano Zampini           ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
7988135c375SStefano Zampini           for (i=0;i<vpf;i++) {
799*8183e3f6SStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[faces[cl*vpf+i]]);CHKERRQ(ierr);
8008135c375SStefano Zampini           }
8018135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
80277eacf09SStefano Zampini         }
80377eacf09SStefano Zampini         bf = bf-nc;
8048135c375SStefano Zampini       }
8058135c375SStefano Zampini     }
80677eacf09SStefano Zampini     if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf);
80777eacf09SStefano Zampini     ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr);
80877eacf09SStefano Zampini     ierr = PetscFree(fcells);CHKERRQ(ierr);
8098135c375SStefano Zampini   }
8108135c375SStefano Zampini 
8118135c375SStefano Zampini   /* mark owned vertices */
8123e6c54aaSStefano Zampini   vown = NULL;
8133e6c54aaSStefano Zampini   if (pown) {
8143e6c54aaSStefano Zampini     ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
8158135c375SStefano Zampini     for (p=cStart;p<cEnd;p++) {
8168135c375SStefano Zampini       PetscInt i,closureSize,*closure = NULL;
8178135c375SStefano Zampini 
8183e6c54aaSStefano Zampini       if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue;
8198135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
8208135c375SStefano Zampini       for (i=0;i<closureSize;i++) {
8218135c375SStefano Zampini         const PetscInt pp = closure[2*i];
8228135c375SStefano Zampini 
8238135c375SStefano Zampini         if (pp >= vStart && pp < vEnd) {
8243e6c54aaSStefano Zampini           ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr);
8258135c375SStefano Zampini         }
8268135c375SStefano Zampini       }
8278135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
8288135c375SStefano Zampini     }
8298135c375SStefano Zampini   }
8308135c375SStefano Zampini 
8318135c375SStefano Zampini   /* vertex_parents (Non-conforming meshes) */
8328135c375SStefano Zampini   parentSection  = NULL;
8338135c375SStefano Zampini   if (enable_ncmesh) {
8348135c375SStefano Zampini     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
8358135c375SStefano Zampini   }
8368135c375SStefano Zampini   if (parentSection) {
8378135c375SStefano Zampini     PetscInt vp,gvp;
8388135c375SStefano Zampini 
8398135c375SStefano Zampini     for (vp=0,p=vStart;p<vEnd;p++) {
8408135c375SStefano Zampini       DMLabel  dlabel;
8418135c375SStefano Zampini       PetscInt parent,depth;
8428135c375SStefano Zampini 
8433e6c54aaSStefano Zampini       if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
8448135c375SStefano Zampini       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
8458135c375SStefano Zampini       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
8468135c375SStefano Zampini       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
8478135c375SStefano Zampini       if (parent != p) vp++;
8488135c375SStefano Zampini     }
8498135c375SStefano Zampini     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
8508135c375SStefano Zampini     if (gvp) {
8518135c375SStefano Zampini       PetscInt  maxsupp;
8528135c375SStefano Zampini       PetscBool *skip = NULL;
8538135c375SStefano Zampini 
8548135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
8558135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
8568135c375SStefano Zampini       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
8578135c375SStefano Zampini       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
8588135c375SStefano Zampini       for (p=vStart;p<vEnd;p++) {
8598135c375SStefano Zampini         DMLabel  dlabel;
8608135c375SStefano Zampini         PetscInt parent;
8618135c375SStefano Zampini 
8623e6c54aaSStefano Zampini         if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue;
8638135c375SStefano Zampini         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
8648135c375SStefano Zampini         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
8658135c375SStefano Zampini         if (parent != p) {
866d42aff11SStefano Zampini           int            vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */
8673924b612SStefano Zampini           PetscInt       i,nv,ssize,n,numChildren,depth = -1;
8688135c375SStefano Zampini           const PetscInt *children;
8693924b612SStefano Zampini 
8703924b612SStefano Zampini           ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr);
8713924b612SStefano Zampini           switch (ssize) {
8728135c375SStefano Zampini             case 2: /* edge */
8738135c375SStefano Zampini               nv   = 0;
874cc0d3ed7SStefano Zampini               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr);
87577eacf09SStefano Zampini               ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr);
8768135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
8778135c375SStefano Zampini               for (i=0;i<nv;i++) {
878bfadf1d8SStefano Zampini                 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr);
8798135c375SStefano Zampini               }
8808135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
8818135c375SStefano Zampini               vp--;
8828135c375SStefano Zampini               break;
8838135c375SStefano Zampini             case 4: /* face */
8848135c375SStefano Zampini               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
8858135c375SStefano Zampini               for (n=0;n<numChildren;n++) {
8868135c375SStefano Zampini                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
8878135c375SStefano Zampini                 if (!depth) {
8888135c375SStefano Zampini                   const PetscInt *hvsupp,*hesupp,*cone;
8898135c375SStefano Zampini                   PetscInt       hvsuppSize,hesuppSize,coneSize;
890451a39c7SStefano Zampini                   PetscInt       hv = children[n],he = -1,f;
8918135c375SStefano Zampini 
892451a39c7SStefano Zampini                   ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
8938135c375SStefano Zampini                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
8948135c375SStefano Zampini                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
8958135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
8968135c375SStefano Zampini                     PetscInt ep;
8978135c375SStefano Zampini                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
8988135c375SStefano Zampini                     if (ep != hvsupp[i]) {
8998135c375SStefano Zampini                       he = hvsupp[i];
9008135c375SStefano Zampini                     } else {
9018135c375SStefano Zampini                       skip[i] = PETSC_TRUE;
9028135c375SStefano Zampini                     }
9038135c375SStefano Zampini                   }
904451a39c7SStefano Zampini                   if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
9058135c375SStefano Zampini                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
90611a4995dSStefano Zampini                   vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
9078135c375SStefano Zampini                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
9088135c375SStefano Zampini                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
9098135c375SStefano Zampini                   for (f=0;f<hesuppSize;f++) {
9108135c375SStefano Zampini                     PetscInt j;
9118135c375SStefano Zampini 
9128135c375SStefano Zampini                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
9138135c375SStefano Zampini                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
9148135c375SStefano Zampini                     for (j=0;j<coneSize;j++) {
9158135c375SStefano Zampini                       PetscInt k;
9168135c375SStefano Zampini                       for (k=0;k<hvsuppSize;k++) {
9178135c375SStefano Zampini                         if (hvsupp[k] == cone[j]) {
9188135c375SStefano Zampini                           skip[k] = PETSC_TRUE;
9198135c375SStefano Zampini                           break;
9208135c375SStefano Zampini                         }
9218135c375SStefano Zampini                       }
9228135c375SStefano Zampini                     }
9238135c375SStefano Zampini                   }
9248135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
9258135c375SStefano Zampini                     if (!skip[i]) {
9268135c375SStefano Zampini                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
92711a4995dSStefano Zampini                       vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]);
9288135c375SStefano Zampini                     }
9298135c375SStefano Zampini                   }
9308135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
9318135c375SStefano Zampini                   for (i=0;i<2;i++) {
932bfadf1d8SStefano Zampini                     ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)(vids[i]-vStart));CHKERRQ(ierr);
9338135c375SStefano Zampini                   }
9348135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9358135c375SStefano Zampini                   vp--;
9368135c375SStefano Zampini                 }
9378135c375SStefano Zampini               }
9388135c375SStefano Zampini               break;
9398135c375SStefano Zampini             default:
9403924b612SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize);
9418135c375SStefano Zampini           }
9428135c375SStefano Zampini         }
9438135c375SStefano Zampini       }
9448135c375SStefano Zampini       ierr = PetscFree(skip);CHKERRQ(ierr);
9458135c375SStefano Zampini     }
9468135c375SStefano Zampini     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
9478135c375SStefano Zampini   }
9488135c375SStefano Zampini   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
9493e6c54aaSStefano Zampini   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
9508135c375SStefano Zampini 
9518135c375SStefano Zampini   /* vertices */
95277eacf09SStefano Zampini   if (hovec) { /* higher-order meshes */
95377eacf09SStefano Zampini     const char *fec;
95477eacf09SStefano Zampini     PetscInt   i,n;
95577eacf09SStefano Zampini 
95677eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
95777eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
95877eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
95977eacf09SStefano Zampini     ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr);
96077eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
9613e7cab22SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr);
96277eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
96377eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
96477eacf09SStefano Zampini     ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr);
96577eacf09SStefano Zampini     ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr);
96677eacf09SStefano Zampini     if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim);
96777eacf09SStefano Zampini     for (i=0;i<n/sdim;i++) {
96877eacf09SStefano Zampini       PetscInt s;
96977eacf09SStefano Zampini       for (s=0;s<sdim;s++) {
97077eacf09SStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr);
97177eacf09SStefano Zampini       }
97277eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
97377eacf09SStefano Zampini     }
97477eacf09SStefano Zampini     ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr);
97577eacf09SStefano Zampini   } else {
976cc0d3ed7SStefano Zampini     ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr);
9778135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
978cc0d3ed7SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr);
9798135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
9808135c375SStefano Zampini     ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
981cc0d3ed7SStefano Zampini     for (p=0;p<nvert/sdim;p++) {
982cc0d3ed7SStefano Zampini       PetscInt s;
983cc0d3ed7SStefano Zampini       for (s=0;s<sdim;s++) {
98477eacf09SStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr);
9858135c375SStefano Zampini       }
986cc0d3ed7SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
9878135c375SStefano Zampini     }
9888135c375SStefano Zampini     ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
98977eacf09SStefano Zampini   }
9908135c375SStefano Zampini   PetscFunctionReturn(0);
9918135c375SStefano Zampini }
9928135c375SStefano Zampini 
9938135c375SStefano Zampini /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
9948135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
9958135c375SStefano Zampini {
9968135c375SStefano Zampini   PetscErrorCode ierr;
9978135c375SStefano Zampini   PetscBool      isglvis,isascii;
9988135c375SStefano Zampini 
9998135c375SStefano Zampini   PetscFunctionBegin;
10008135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10018135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
10028135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
10038135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
10048135c375SStefano Zampini   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
10058135c375SStefano Zampini   if (isglvis) {
10068135c375SStefano Zampini     PetscViewer          view;
10078135c375SStefano Zampini     PetscViewerGLVisType type;
10088135c375SStefano Zampini 
10098135c375SStefano Zampini     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
10108135c375SStefano Zampini     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
10118135c375SStefano Zampini     if (view) { /* in the socket case, it may happen that the connection failed */
10128135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
10138135c375SStefano Zampini         PetscMPIInt size,rank;
10143924b612SStefano Zampini 
10158135c375SStefano Zampini         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
10168135c375SStefano Zampini         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
10178135c375SStefano Zampini         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
10188135c375SStefano Zampini       }
10198135c375SStefano Zampini       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
10208135c375SStefano Zampini       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
10218135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
10228135c375SStefano Zampini         PetscInt    dim;
10238135c375SStefano Zampini         const char* name;
10248135c375SStefano Zampini 
10258135c375SStefano Zampini         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
10268135c375SStefano Zampini         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
10278135c375SStefano Zampini         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
10288135c375SStefano Zampini         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
10298135c375SStefano Zampini       }
10308135c375SStefano Zampini     }
10318135c375SStefano Zampini   } else {
10328135c375SStefano Zampini     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
10338135c375SStefano Zampini   }
10348135c375SStefano Zampini   PetscFunctionReturn(0);
10358135c375SStefano Zampini }
1036