18135c375SStefano Zampini #include <petsc/private/glvisviewerimpl.h> 28135c375SStefano Zampini #include <petsc/private/petscimpl.h> 38135c375SStefano Zampini #include <petsc/private/dmpleximpl.h> 48135c375SStefano Zampini #include <petscbt.h> 58135c375SStefano Zampini #include <petscdmplex.h> 68135c375SStefano Zampini #include <petscsf.h> 78135c375SStefano Zampini #include <petscds.h> 88135c375SStefano Zampini 98135c375SStefano Zampini typedef struct { 108135c375SStefano Zampini PetscInt nf; 118135c375SStefano Zampini VecScatter *scctx; 128135c375SStefano Zampini } GLVisViewerCtx; 138135c375SStefano Zampini 148135c375SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 158135c375SStefano Zampini { 168135c375SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 178135c375SStefano Zampini PetscInt i; 188135c375SStefano Zampini PetscErrorCode ierr; 198135c375SStefano Zampini 208135c375SStefano Zampini PetscFunctionBegin; 218135c375SStefano Zampini for (i=0;i<ctx->nf;i++) { 228135c375SStefano Zampini ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr); 238135c375SStefano Zampini } 248135c375SStefano Zampini ierr = PetscFree(ctx->scctx);CHKERRQ(ierr); 2551f51421SSatish Balay ierr = PetscFree(vctx);CHKERRQ(ierr); 268135c375SStefano Zampini PetscFunctionReturn(0); 278135c375SStefano Zampini } 288135c375SStefano Zampini 298135c375SStefano Zampini static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 308135c375SStefano Zampini { 318135c375SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 328135c375SStefano Zampini PetscInt f; 338135c375SStefano Zampini PetscErrorCode ierr; 348135c375SStefano Zampini 358135c375SStefano Zampini PetscFunctionBegin; 368135c375SStefano Zampini for (f=0;f<nf;f++) { 378135c375SStefano Zampini ierr = VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 388135c375SStefano Zampini ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 398135c375SStefano Zampini } 408135c375SStefano Zampini PetscFunctionReturn(0); 418135c375SStefano Zampini } 428135c375SStefano Zampini 438135c375SStefano Zampini /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */ 448135c375SStefano Zampini PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer) 458135c375SStefano Zampini { 468135c375SStefano Zampini DM dm = (DM)odm; 474cac2994SStefano Zampini Vec xlocal,xfield,*Ufield; 488135c375SStefano Zampini PetscDS ds; 498135c375SStefano Zampini IS globalNum,isfield; 508135c375SStefano Zampini PetscBT vown; 518135c375SStefano Zampini char **fieldname = NULL,**fec_type = NULL; 528135c375SStefano Zampini const PetscInt *gNum; 53bb77a09fSStefano Zampini PetscInt *nlocal,*bs,*idxs,*dims; 548135c375SStefano Zampini PetscInt f,maxfields,nfields,c,totc,totdofs,Nv,cum,i; 55b135d7daSStefano Zampini PetscInt dim,cStart,cEnd,vStart,vEnd; 568135c375SStefano Zampini GLVisViewerCtx *ctx; 578135c375SStefano Zampini PetscSection s; 588135c375SStefano Zampini PetscErrorCode ierr; 598135c375SStefano Zampini 608135c375SStefano Zampini PetscFunctionBegin; 618135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 628135c375SStefano Zampini ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 638135c375SStefano Zampini ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 64b135d7daSStefano Zampini ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 65b135d7daSStefano Zampini if (!globalNum) { 66b135d7daSStefano Zampini ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 67b135d7daSStefano Zampini ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 68b135d7daSStefano Zampini ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 69b135d7daSStefano Zampini } 708135c375SStefano Zampini ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 718135c375SStefano Zampini ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 728135c375SStefano Zampini for (c = cStart, totc = 0; c < cEnd; c++) { 738135c375SStefano Zampini if (gNum[c-cStart] >= 0) { 748135c375SStefano Zampini PetscInt i,numPoints,*points = NULL; 758135c375SStefano Zampini 768135c375SStefano Zampini totc++; 778135c375SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 788135c375SStefano Zampini for (i=0;i<numPoints*2;i+= 2) { 798135c375SStefano Zampini if ((points[i] >= vStart) && (points[i] < vEnd)) { 808135c375SStefano Zampini ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr); 818135c375SStefano Zampini } 828135c375SStefano Zampini } 838135c375SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 848135c375SStefano Zampini } 858135c375SStefano Zampini } 8677eacf09SStefano Zampini for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++; 878135c375SStefano Zampini 888135c375SStefano Zampini ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr); 898135c375SStefano Zampini ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr); 90e87a4003SBarry Smith ierr = DMGetSection(dm,&s);CHKERRQ(ierr); 918135c375SStefano Zampini ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr); 928135c375SStefano Zampini for (f=0,maxfields=0;f<nfields;f++) { 938135c375SStefano Zampini PetscInt bs; 948135c375SStefano Zampini 958135c375SStefano Zampini ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr); 968135c375SStefano Zampini maxfields += bs; 978135c375SStefano Zampini } 984cac2994SStefano Zampini ierr = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);CHKERRQ(ierr); 998135c375SStefano Zampini ierr = PetscNew(&ctx);CHKERRQ(ierr); 1008135c375SStefano Zampini ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr); 1018135c375SStefano Zampini ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 1028135c375SStefano Zampini if (ds) { 1038135c375SStefano Zampini for (f=0;f<nfields;f++) { 1048135c375SStefano Zampini const char* fname; 1058135c375SStefano Zampini char name[256]; 1068135c375SStefano Zampini PetscObject disc; 1078135c375SStefano Zampini size_t len; 1088135c375SStefano Zampini 1098135c375SStefano Zampini ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr); 1108135c375SStefano Zampini ierr = PetscStrlen(fname,&len);CHKERRQ(ierr); 1118135c375SStefano Zampini if (len) { 1128135c375SStefano Zampini ierr = PetscStrcpy(name,fname);CHKERRQ(ierr); 1138135c375SStefano Zampini } else { 114bfadf1d8SStefano Zampini ierr = PetscSNPrintf(name,256,"Field%D",f);CHKERRQ(ierr); 1158135c375SStefano Zampini } 1168135c375SStefano Zampini ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 1178135c375SStefano Zampini if (disc) { 1188135c375SStefano Zampini PetscClassId id; 1198135c375SStefano Zampini PetscInt Nc; 1208135c375SStefano Zampini char fec[64]; 1218135c375SStefano Zampini 1228135c375SStefano Zampini ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 1238135c375SStefano Zampini if (id == PETSCFE_CLASSID) { 1248135c375SStefano Zampini PetscFE fem = (PetscFE)disc; 1258135c375SStefano Zampini PetscDualSpace sp; 1268135c375SStefano Zampini PetscDualSpaceType spname; 1278135c375SStefano Zampini PetscInt order; 1288135c375SStefano Zampini PetscBool islag,continuous,H1 = PETSC_TRUE; 1298135c375SStefano Zampini 1308135c375SStefano Zampini ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr); 1318135c375SStefano Zampini ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr); 1328135c375SStefano Zampini ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr); 1338135c375SStefano Zampini ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr); 1348135c375SStefano Zampini if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space"); 1358135c375SStefano Zampini ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 1368135c375SStefano Zampini if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported"); 1378135c375SStefano Zampini ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 1388135c375SStefano Zampini if (continuous && order > 0) { 139bfadf1d8SStefano Zampini ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr); 1408135c375SStefano Zampini } else { 1418135c375SStefano Zampini H1 = PETSC_FALSE; 142bfadf1d8SStefano Zampini ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order);CHKERRQ(ierr); 1438135c375SStefano Zampini } 1448135c375SStefano Zampini ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr); 1458135c375SStefano Zampini bs[ctx->nf] = Nc; 146bb77a09fSStefano Zampini dims[ctx->nf] = dim; 1478135c375SStefano Zampini if (H1) { 1488135c375SStefano Zampini nlocal[ctx->nf] = Nc * Nv; 1498135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 1508135c375SStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr); 1518135c375SStefano Zampini for (i=0,cum=0;i<vEnd-vStart;i++) { 1528135c375SStefano Zampini PetscInt j,off; 1538135c375SStefano Zampini 1548135c375SStefano Zampini if (PetscUnlikely(!PetscBTLookup(vown,i))) continue; 1558135c375SStefano Zampini ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr); 1568135c375SStefano Zampini for (j=0;j<Nc;j++) idxs[cum++] = off + j; 1578135c375SStefano Zampini } 1588135c375SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 1598135c375SStefano Zampini } else { 1608135c375SStefano Zampini nlocal[ctx->nf] = Nc * totc; 1618135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 1628135c375SStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr); 1638135c375SStefano Zampini for (i=0,cum=0;i<cEnd-cStart;i++) { 1648135c375SStefano Zampini PetscInt j,off; 1658135c375SStefano Zampini 1668135c375SStefano Zampini if (PetscUnlikely(gNum[i] < 0)) continue; 1678135c375SStefano Zampini ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 1688135c375SStefano Zampini for (j=0;j<Nc;j++) idxs[cum++] = off + j; 1698135c375SStefano Zampini } 1708135c375SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 1718135c375SStefano Zampini } 17235928de7SBarry Smith ierr = VecScatterCreateWithData(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 1738135c375SStefano Zampini ierr = VecDestroy(&xfield);CHKERRQ(ierr); 1748135c375SStefano Zampini ierr = ISDestroy(&isfield);CHKERRQ(ierr); 1758135c375SStefano Zampini ctx->nf++; 1768135c375SStefano Zampini } else if (id == PETSCFV_CLASSID) { 1778135c375SStefano Zampini PetscInt c; 1788135c375SStefano Zampini 1798135c375SStefano Zampini ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr); 180bfadf1d8SStefano Zampini ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim);CHKERRQ(ierr); 1818135c375SStefano Zampini for (c = 0; c < Nc; c++) { 1828135c375SStefano Zampini char comp[256]; 183bfadf1d8SStefano Zampini ierr = PetscSNPrintf(comp,256,"%s-Comp%D",name,c);CHKERRQ(ierr); 1848135c375SStefano Zampini ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr); 1858135c375SStefano Zampini bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ 1868135c375SStefano Zampini nlocal[ctx->nf] = totc; 187bb77a09fSStefano Zampini dims[ctx->nf] = dim; 1888135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 1898135c375SStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr); 1908135c375SStefano Zampini for (i=0,cum=0;i<cEnd-cStart;i++) { 1918135c375SStefano Zampini PetscInt off; 1928135c375SStefano Zampini 1938135c375SStefano Zampini if (PetscUnlikely(gNum[i])<0) continue; 1948135c375SStefano Zampini ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 1958135c375SStefano Zampini idxs[cum++] = off + c; 1968135c375SStefano Zampini } 1978135c375SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 19835928de7SBarry Smith ierr = VecScatterCreateWithData(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 1998135c375SStefano Zampini ierr = VecDestroy(&xfield);CHKERRQ(ierr); 2008135c375SStefano Zampini ierr = ISDestroy(&isfield);CHKERRQ(ierr); 2018135c375SStefano Zampini ctx->nf++; 2028135c375SStefano Zampini } 2038135c375SStefano Zampini } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f); 2048135c375SStefano Zampini } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f); 2058135c375SStefano Zampini } 2068135c375SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM"); 2078135c375SStefano Zampini ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 2088135c375SStefano Zampini ierr = VecDestroy(&xlocal);CHKERRQ(ierr); 2098135c375SStefano Zampini ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 2108135c375SStefano Zampini 2114cac2994SStefano Zampini /* create work vectors */ 2124cac2994SStefano Zampini for (f=0;f<ctx->nf;f++) { 2134cac2994SStefano Zampini ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr); 2144cac2994SStefano Zampini ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr); 2154cac2994SStefano Zampini ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr); 2164cac2994SStefano Zampini ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr); 2174cac2994SStefano Zampini } 2184cac2994SStefano Zampini 2198135c375SStefano Zampini /* customize the viewer */ 2204cac2994SStefano Zampini ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 2218135c375SStefano Zampini for (f=0;f<ctx->nf;f++) { 2228135c375SStefano Zampini ierr = PetscFree(fieldname[f]);CHKERRQ(ierr); 2238135c375SStefano Zampini ierr = PetscFree(fec_type[f]);CHKERRQ(ierr); 2244cac2994SStefano Zampini ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr); 2258135c375SStefano Zampini } 2264cac2994SStefano Zampini ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr); 2278135c375SStefano Zampini PetscFunctionReturn(0); 2288135c375SStefano Zampini } 2298135c375SStefano Zampini 230b135d7daSStefano Zampini typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid; 2318135c375SStefano Zampini 2328135c375SStefano Zampini MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, 2338135c375SStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, 2348135c375SStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_UNDEF}, 235b135d7daSStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_PRISM,MFEM_CUBE } }; 2368135c375SStefano Zampini 237b135d7daSStefano Zampini MFEM_cid mfem_table_cid_unint[4][9] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 238b135d7daSStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 239b135d7daSStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 240b135d7daSStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_CUBE } }; 241044a5661SStefano Zampini 242*f86f7544SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid) 2438135c375SStefano Zampini { 2448135c375SStefano Zampini DMLabel dlabel; 245044a5661SStefano Zampini PetscInt depth,csize,pdepth,dim; 2468135c375SStefano Zampini PetscErrorCode ierr; 2478135c375SStefano Zampini 2488135c375SStefano Zampini PetscFunctionBegin; 2498135c375SStefano Zampini ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 250044a5661SStefano Zampini ierr = DMLabelGetValue(dlabel,p,&pdepth);CHKERRQ(ierr); 2518135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr); 252044a5661SStefano Zampini ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 253044a5661SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 2548135c375SStefano Zampini if (label) { 2558135c375SStefano Zampini ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr); 256*f86f7544SStefano Zampini *mid = *mid - minl + 1; /* MFEM does not like negative markers */ 25777eacf09SStefano Zampini } else *mid = 1; 258044a5661SStefano Zampini if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 259044a5661SStefano Zampini #if defined PETSC_USE_DEBUG 260044a5661SStefano Zampini if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim); 261044a5661SStefano Zampini if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p); 262044a5661SStefano 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); 263044a5661SStefano Zampini #endif 264044a5661SStefano Zampini *cid = mfem_table_cid_unint[dim][csize]; 265044a5661SStefano Zampini } else { 266044a5661SStefano Zampini #if defined PETSC_USE_DEBUG 267044a5661SStefano Zampini if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p); 268044a5661SStefano Zampini if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p); 269044a5661SStefano Zampini #endif 270044a5661SStefano Zampini *cid = mfem_table_cid[pdepth][csize]; 271044a5661SStefano Zampini } 2728135c375SStefano Zampini PetscFunctionReturn(0); 2738135c375SStefano Zampini } 2748135c375SStefano Zampini 275cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[]) 2768135c375SStefano Zampini { 277cc0d3ed7SStefano Zampini PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; 2788135c375SStefano Zampini PetscErrorCode ierr; 2798135c375SStefano Zampini 2808135c375SStefano Zampini PetscFunctionBegin; 2818135c375SStefano Zampini ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 282cc0d3ed7SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 283cc0d3ed7SStefano Zampini sdim = dim; 284cc0d3ed7SStefano Zampini if (csec) { 28584f354e3SLisandro Dalcin PetscInt sStart,sEnd; 28684f354e3SLisandro Dalcin 287cc0d3ed7SStefano Zampini ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 28884f354e3SLisandro Dalcin ierr = PetscSectionGetChart(csec,&sStart,&sEnd);CHKERRQ(ierr); 289cc0d3ed7SStefano Zampini ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); 290cc0d3ed7SStefano Zampini off = off/sdim; 29184f354e3SLisandro Dalcin if (p >= sStart && p < sEnd) { 292cc0d3ed7SStefano Zampini ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 293cc0d3ed7SStefano Zampini } 29484f354e3SLisandro Dalcin } 295cc0d3ed7SStefano Zampini if (!dof) { 2968135c375SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 2978135c375SStefano Zampini for (i=0,q=0;i<numPoints*2;i+= 2) 2988135c375SStefano Zampini if ((points[i] >= vStart) && (points[i] < vEnd)) 2998183e3f6SStefano Zampini vids[q++] = (int)(points[i]-vStart+off); 3008135c375SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 301cc0d3ed7SStefano Zampini } else { 302cc0d3ed7SStefano Zampini ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); 303cc0d3ed7SStefano Zampini ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 3048183e3f6SStefano Zampini for (q=0;q<dof/sdim;q++) vids[q] = (int)(off/sdim + q); 305cc0d3ed7SStefano Zampini } 3068135c375SStefano Zampini *nv = q; 3078135c375SStefano Zampini PetscFunctionReturn(0); 3088135c375SStefano Zampini } 3098135c375SStefano Zampini 310b135d7daSStefano Zampini static PetscErrorCode DMPlexGlvisInvertHybrid(PetscInt cid,int vids[]) 311b135d7daSStefano Zampini { 312b135d7daSStefano Zampini int tmp; 313b135d7daSStefano Zampini 314b135d7daSStefano Zampini PetscFunctionBegin; 315b135d7daSStefano Zampini if (cid == MFEM_SQUARE) { /* PETSc stores hybrid quads not as counter-clockwise quad */ 316b135d7daSStefano Zampini tmp = vids[2]; 317b135d7daSStefano Zampini vids[2] = vids[3]; 318b135d7daSStefano Zampini vids[3] = tmp; 319b135d7daSStefano Zampini } else if (cid == MFEM_PRISM) { /* MFEM uses a different orientation for the base and top triangles of the wedge */ 320b135d7daSStefano Zampini tmp = vids[1]; 321b135d7daSStefano Zampini vids[1] = vids[2]; 322b135d7daSStefano Zampini vids[2] = tmp; 323b135d7daSStefano Zampini tmp = vids[4]; 324b135d7daSStefano Zampini vids[4] = vids[5]; 325b135d7daSStefano Zampini vids[5] = tmp; 326b135d7daSStefano Zampini } 327b135d7daSStefano Zampini PetscFunctionReturn(0); 328b135d7daSStefano Zampini } 329b135d7daSStefano Zampini 33077eacf09SStefano Zampini /* 33177eacf09SStefano Zampini ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 33277eacf09SStefano Zampini Higher order meshes are also supported 33377eacf09SStefano Zampini */ 3348135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 3358135c375SStefano Zampini { 3368135c375SStefano Zampini DMLabel label; 3378135c375SStefano Zampini PetscSection coordSection,parentSection; 33877eacf09SStefano Zampini Vec coordinates,hovec; 3398135c375SStefano Zampini const PetscScalar *array; 340*f86f7544SStefano Zampini PetscInt bf,p,sdim,dim,depth,novl,minl; 341cc0d3ed7SStefano Zampini PetscInt cStart,cEnd,cEndInterior,vStart,vEnd,nvert; 3423924b612SStefano Zampini PetscMPIInt size; 3433e6c54aaSStefano Zampini PetscBool localized,isascii; 3443e6c54aaSStefano Zampini PetscBool enable_mfem,enable_boundary,enable_ncmesh; 3453e6c54aaSStefano Zampini PetscBT pown,vown; 3468135c375SStefano Zampini PetscErrorCode ierr; 3478135c375SStefano Zampini PetscContainer glvis_container; 348044a5661SStefano Zampini PetscBool cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE; 349*f86f7544SStefano Zampini PetscBool enable_emark,enable_bmark; 35077eacf09SStefano Zampini const char *fmt; 3517bf4dd16SStefano Zampini char emark[64] = "",bmark[64] = ""; 3528135c375SStefano Zampini 3538135c375SStefano Zampini PetscFunctionBegin; 3548135c375SStefano Zampini PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3558135c375SStefano Zampini PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3568135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 3578135c375SStefano Zampini if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 3583924b612SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr); 3593924b612SStefano Zampini if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 3608135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 3618135c375SStefano Zampini 3628135c375SStefano Zampini /* get container: determines if a process visualizes is portion of the data or not */ 3638135c375SStefano Zampini ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 3648135c375SStefano Zampini if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 3658135c375SStefano Zampini { 3668135c375SStefano Zampini PetscViewerGLVisInfo glvis_info; 3678135c375SStefano Zampini ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 3688135c375SStefano Zampini enabled = glvis_info->enabled; 36977eacf09SStefano Zampini fmt = glvis_info->fmt; 3708135c375SStefano Zampini } 37121414b21SStefano Zampini 37221414b21SStefano Zampini /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh 37321414b21SStefano Zampini DMPlex does not currently support HO meshes, so there's no API for this */ 37421414b21SStefano Zampini ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); 37521414b21SStefano Zampini 376cc0d3ed7SStefano Zampini ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 3778135c375SStefano Zampini ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 378b135d7daSStefano Zampini cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior; 3798135c375SStefano Zampini ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 38077eacf09SStefano Zampini ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); 381cc0d3ed7SStefano Zampini ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 38221414b21SStefano Zampini if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); 383cc0d3ed7SStefano Zampini ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 38477eacf09SStefano Zampini ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 38577eacf09SStefano Zampini ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 38621414b21SStefano Zampini if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector"); 3878135c375SStefano Zampini 3888135c375SStefano Zampini /* 3898135c375SStefano Zampini a couple of sections of the mesh specification are disabled 3903e96840aSStefano Zampini - boundary: the boundary is not needed for proper mesh visualization unless we want to visualize boundary attributes or we have high-order coordinates in 3D (topologically) 39177eacf09SStefano Zampini - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 3923e6c54aaSStefano Zampini and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) 3938135c375SStefano Zampini */ 3943e96840aSStefano Zampini enable_boundary = PETSC_FALSE; 3958135c375SStefano Zampini enable_ncmesh = PETSC_FALSE; 3963e6c54aaSStefano Zampini enable_mfem = PETSC_FALSE; 397*f86f7544SStefano Zampini enable_emark = PETSC_FALSE; 398*f86f7544SStefano Zampini enable_bmark = PETSC_FALSE; 3997bf4dd16SStefano Zampini /* I'm tired of problems with negative values in the markers, disable them */ 4008135c375SStefano Zampini ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr); 4013e6c54aaSStefano Zampini ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 4023e6c54aaSStefano 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); 4033e6c54aaSStefano 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); 404*f86f7544SStefano Zampini ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);CHKERRQ(ierr); 405*f86f7544SStefano Zampini ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);CHKERRQ(ierr); 4068135c375SStefano Zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 407*f86f7544SStefano Zampini if (enable_bmark) enable_boundary = PETSC_TRUE; 408*f86f7544SStefano Zampini 4093924b612SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 4103924b612SStefano Zampini if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); 4117e1aca4eSStefano Zampini ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 4127e1aca4eSStefano Zampini if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 4137e1aca4eSStefano Zampini "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); 4147e1aca4eSStefano Zampini if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 4157e1aca4eSStefano Zampini "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0"); 416044a5661SStefano Zampini if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 417044a5661SStefano Zampini if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth); 418044a5661SStefano Zampini cellvertex = PETSC_TRUE; 419044a5661SStefano Zampini } 4208135c375SStefano Zampini 4218135c375SStefano Zampini /* Identify possible cells in the overlap */ 4228135c375SStefano Zampini novl = 0; 4238135c375SStefano Zampini pown = NULL; 4243924b612SStefano Zampini if (size > 1) { 4253e6c54aaSStefano Zampini IS globalNum = NULL; 4263e6c54aaSStefano Zampini const PetscInt *gNum; 4273e6c54aaSStefano Zampini PetscBool ovl = PETSC_FALSE; 4283e6c54aaSStefano Zampini 429b135d7daSStefano Zampini ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 430b135d7daSStefano Zampini if (!globalNum) { 431b135d7daSStefano Zampini ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 432b135d7daSStefano Zampini ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 433b135d7daSStefano Zampini ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 434b135d7daSStefano Zampini } 4358135c375SStefano Zampini ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 4368135c375SStefano Zampini for (p=cStart; p<cEnd; p++) { 4378135c375SStefano Zampini if (gNum[p-cStart] < 0) { 4388135c375SStefano Zampini ovl = PETSC_TRUE; 4398135c375SStefano Zampini novl++; 4408135c375SStefano Zampini } 4418135c375SStefano Zampini } 4428135c375SStefano Zampini if (ovl) { 4438135c375SStefano Zampini /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 4448135c375SStefano Zampini TODO: garbage collector? attach pown to dm? */ 4453e6c54aaSStefano Zampini ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr); 4463e6c54aaSStefano Zampini for (p=cStart; p<cEnd; p++) { 4473e6c54aaSStefano Zampini if (gNum[p-cStart] < 0) continue; 4483e6c54aaSStefano Zampini else { 4493e6c54aaSStefano Zampini ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr); 4508135c375SStefano Zampini } 4518135c375SStefano Zampini } 4523e6c54aaSStefano Zampini } 4533e6c54aaSStefano Zampini ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 4543e6c54aaSStefano Zampini } 4558135c375SStefano Zampini 4563e6c54aaSStefano Zampini /* return if this process is disabled */ 4578135c375SStefano Zampini if (!enabled) { 4588135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 4598135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 4608135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 4618135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 4628135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 4638135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 4648135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 4658135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 4668135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 4678135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 4688135c375SStefano Zampini ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 4698135c375SStefano Zampini PetscFunctionReturn(0); 4708135c375SStefano Zampini } 4718135c375SStefano Zampini 4723e6c54aaSStefano Zampini if (enable_mfem) { 4733e6c54aaSStefano Zampini if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */ 4743e6c54aaSStefano Zampini PetscInt vpc = 0; 4753e6c54aaSStefano Zampini char fec[64]; 4763e6c54aaSStefano Zampini int vids[8] = {0,1,2,3,4,5,6,7}; 477044a5661SStefano Zampini int hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3}; 478044a5661SStefano Zampini int quadv[8] = {0,1,3,2}, triv[3] = {0,1,2}; 479044a5661SStefano Zampini int *dof = NULL; 4803e6c54aaSStefano Zampini PetscScalar *array,*ptr; 4813e6c54aaSStefano Zampini 482bfadf1d8SStefano Zampini ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr); 483b135d7daSStefano Zampini if (cEndInterior < cEnd) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Support for hybrid meshed not currently implemented"); 4843e6c54aaSStefano Zampini if (cEnd-cStart) { 4853e6c54aaSStefano Zampini PetscInt fpc; 4863e6c54aaSStefano Zampini 4873e6c54aaSStefano Zampini ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 4883e6c54aaSStefano Zampini switch(dim) { 4893e6c54aaSStefano Zampini case 1: 4903e6c54aaSStefano Zampini vpc = 2; 4913e6c54aaSStefano Zampini dof = hexv; 4923e6c54aaSStefano Zampini break; 4933e6c54aaSStefano Zampini case 2: 4943e6c54aaSStefano Zampini switch (fpc) { 4953e6c54aaSStefano Zampini case 3: 4963e6c54aaSStefano Zampini vpc = 3; 497044a5661SStefano Zampini dof = triv; 4983e6c54aaSStefano Zampini break; 4993e6c54aaSStefano Zampini case 4: 5003e6c54aaSStefano Zampini vpc = 4; 501044a5661SStefano Zampini dof = quadv; 5023e6c54aaSStefano Zampini break; 5033e6c54aaSStefano Zampini default: 5043e6c54aaSStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 5053e6c54aaSStefano Zampini break; 5063e6c54aaSStefano Zampini } 5073e6c54aaSStefano Zampini break; 5083e6c54aaSStefano Zampini case 3: 5093e6c54aaSStefano Zampini switch (fpc) { 510044a5661SStefano Zampini case 4: /* TODO: still need to understand L2 ordering for tets */ 5113e6c54aaSStefano Zampini vpc = 4; 512044a5661SStefano Zampini dof = tetv; 513044a5661SStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case"); 5143e6c54aaSStefano Zampini break; 5153e6c54aaSStefano Zampini case 6: 516044a5661SStefano Zampini if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc); 517044a5661SStefano Zampini vpc = 8; 518044a5661SStefano Zampini dof = hexv; 519044a5661SStefano Zampini break; 520044a5661SStefano Zampini case 8: 521044a5661SStefano Zampini if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 5223e6c54aaSStefano Zampini vpc = 8; 5233e6c54aaSStefano Zampini dof = hexv; 5243e6c54aaSStefano Zampini break; 5253e6c54aaSStefano Zampini default: 5263e6c54aaSStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 5273e6c54aaSStefano Zampini break; 5283e6c54aaSStefano Zampini } 5293e6c54aaSStefano Zampini break; 5303e6c54aaSStefano Zampini default: 5313e6c54aaSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 5323e6c54aaSStefano Zampini break; 5333e6c54aaSStefano Zampini } 5343e6c54aaSStefano Zampini ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr); 5353e6c54aaSStefano Zampini } 536044a5661SStefano Zampini if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs"); 5373e6c54aaSStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr); 5383e6c54aaSStefano Zampini ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr); 5393e6c54aaSStefano Zampini ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr); 5403e6c54aaSStefano Zampini ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr); 5413e6c54aaSStefano Zampini ierr = VecGetArray(hovec,&array);CHKERRQ(ierr); 5423e6c54aaSStefano Zampini ptr = array; 5433e6c54aaSStefano Zampini for (p=cStart;p<cEnd;p++) { 5443e6c54aaSStefano Zampini PetscInt csize,v,d; 5453e6c54aaSStefano Zampini PetscScalar *vals = NULL; 5463e6c54aaSStefano Zampini 5473e6c54aaSStefano Zampini if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 5483e6c54aaSStefano Zampini ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 5493e6c54aaSStefano 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); 5503e6c54aaSStefano Zampini for (v=0;v<vpc;v++) { 5513e6c54aaSStefano Zampini for (d=0;d<sdim;d++) { 5523e6c54aaSStefano Zampini ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d]; 5533e6c54aaSStefano Zampini } 5543e6c54aaSStefano Zampini } 5553e6c54aaSStefano Zampini ptr += vpc*sdim; 5563e6c54aaSStefano Zampini ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 5573e6c54aaSStefano Zampini } 5583e6c54aaSStefano Zampini ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr); 5593e6c54aaSStefano Zampini } 5603e6c54aaSStefano Zampini } 5613e96840aSStefano Zampini /* if we have high-order coordinates in 3D, we need to specify the boundary */ 5623e96840aSStefano Zampini if (hovec && dim == 3) enable_boundary = PETSC_TRUE; 5633e6c54aaSStefano Zampini 5648135c375SStefano Zampini /* header */ 5658135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 5668135c375SStefano Zampini 5678135c375SStefano Zampini /* topological dimension */ 5688135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 5698135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 5708135c375SStefano Zampini 5718135c375SStefano Zampini /* elements */ 572*f86f7544SStefano Zampini minl = 1; 573*f86f7544SStefano Zampini label = NULL; 574*f86f7544SStefano Zampini if (enable_emark) { 575*f86f7544SStefano Zampini PetscInt lminl = PETSC_MAX_INT; 576*f86f7544SStefano Zampini 5777bf4dd16SStefano Zampini ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr); 578*f86f7544SStefano Zampini if (label) { 579*f86f7544SStefano Zampini IS vals; 580*f86f7544SStefano Zampini PetscInt ldef; 581*f86f7544SStefano Zampini 582*f86f7544SStefano Zampini ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr); 583*f86f7544SStefano Zampini ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr); 584*f86f7544SStefano Zampini ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr); 585*f86f7544SStefano Zampini ierr = ISDestroy(&vals);CHKERRQ(ierr); 586*f86f7544SStefano Zampini lminl = PetscMin(ldef,lminl); 587*f86f7544SStefano Zampini } 588*f86f7544SStefano Zampini ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 589*f86f7544SStefano Zampini if (minl == PETSC_MAX_INT) minl = 1; 590*f86f7544SStefano Zampini } 5918135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 5928135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr); 5938135c375SStefano Zampini for (p=cStart;p<cEnd;p++) { 59411a4995dSStefano Zampini int vids[8]; 59511a4995dSStefano Zampini PetscInt i,nv = 0,cid = -1,mid = 1; 5968135c375SStefano Zampini 5973e6c54aaSStefano Zampini if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 598*f86f7544SStefano Zampini ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr); 59977eacf09SStefano Zampini ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 60077eacf09SStefano Zampini ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 601b135d7daSStefano Zampini if (p >= cEndInterior) { 602b135d7daSStefano Zampini ierr = DMPlexGlvisInvertHybrid(cid,vids);CHKERRQ(ierr); 603b135d7daSStefano Zampini } 6048135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 6058135c375SStefano Zampini for (i=0;i<nv;i++) { 606bfadf1d8SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr); 6078135c375SStefano Zampini } 6088135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6098135c375SStefano Zampini } 6108135c375SStefano Zampini 611cc0d3ed7SStefano Zampini /* boundary */ 612cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 613cc0d3ed7SStefano Zampini if (!enable_boundary) { 614cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 615cc0d3ed7SStefano Zampini } else { 61677eacf09SStefano Zampini DMLabel perLabel; 61777eacf09SStefano Zampini PetscBT bfaces; 618b135d7daSStefano Zampini PetscInt fStart,fEnd,*fcells; 61977eacf09SStefano Zampini PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0; 620b135d7daSStefano Zampini PetscInt *facesH = NULL,fpcH = 0,vpfH = 0, vpcH = 0; 621cc0d3ed7SStefano Zampini PetscInt fv1[] = {0,1}, 622cc0d3ed7SStefano Zampini fv2tri[] = {0,1, 623cc0d3ed7SStefano Zampini 1,2, 624cc0d3ed7SStefano Zampini 2,0}, 625cc0d3ed7SStefano Zampini fv2quad[] = {0,1, 626cc0d3ed7SStefano Zampini 1,2, 627cc0d3ed7SStefano Zampini 2,3, 628cc0d3ed7SStefano Zampini 3,0}, 629b135d7daSStefano Zampini fv2quadH[] = {0,1, 630b135d7daSStefano Zampini 2,3, 631b135d7daSStefano Zampini 0,2, 632b135d7daSStefano Zampini 1,3}, 633cc0d3ed7SStefano Zampini fv3tet[] = {0,1,2, 634cc0d3ed7SStefano Zampini 0,3,1, 635cc0d3ed7SStefano Zampini 0,2,3, 636cc0d3ed7SStefano Zampini 2,1,3}, 637b135d7daSStefano Zampini fv3wedge[] = {0,1,2,-1, 638b135d7daSStefano Zampini 3,4,5,-1, 639b135d7daSStefano Zampini 0,1,3,4, 640b135d7daSStefano Zampini 1,2,4,5, 641b135d7daSStefano Zampini 2,0,5,3}, 642cc0d3ed7SStefano Zampini fv3hex[] = {0,1,2,3, 643cc0d3ed7SStefano Zampini 4,5,6,7, 644cc0d3ed7SStefano Zampini 0,3,5,4, 645cc0d3ed7SStefano Zampini 2,1,7,6, 646cc0d3ed7SStefano Zampini 3,2,6,5, 647cc0d3ed7SStefano Zampini 0,4,7,1}; 648cc0d3ed7SStefano Zampini 6498135c375SStefano Zampini /* determine orientation of boundary mesh */ 6508135c375SStefano Zampini if (cEnd-cStart) { 651b135d7daSStefano Zampini if (cEndInterior < cEnd) { 652b135d7daSStefano Zampini ierr = DMPlexGetConeSize(dm,cEndInterior,&fpcH);CHKERRQ(ierr); 653b135d7daSStefano Zampini } 654b135d7daSStefano Zampini if (cEndInterior > cStart) { 6558135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 656b135d7daSStefano Zampini } 6578135c375SStefano Zampini switch(dim) { 6588135c375SStefano Zampini case 1: 6598135c375SStefano Zampini if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc); 6608135c375SStefano Zampini faces = fv1; 6618135c375SStefano Zampini vpf = 1; 66277eacf09SStefano Zampini vpc = 2; 6638135c375SStefano Zampini break; 6648135c375SStefano Zampini case 2: 6658135c375SStefano Zampini switch (fpc) { 666b135d7daSStefano Zampini case 0: 6678135c375SStefano Zampini case 3: 6688135c375SStefano Zampini faces = fv2tri; 6698135c375SStefano Zampini vpf = 2; 67077eacf09SStefano Zampini vpc = 3; 671b135d7daSStefano Zampini if (fpcH == 4) { 672b135d7daSStefano Zampini facesH = fv2quadH; 673b135d7daSStefano Zampini vpfH = 2; 674b135d7daSStefano Zampini vpcH = 4; 675b135d7daSStefano Zampini } else if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 6768135c375SStefano Zampini break; 6778135c375SStefano Zampini case 4: 6788135c375SStefano Zampini faces = fv2quad; 6798135c375SStefano Zampini vpf = 2; 68077eacf09SStefano Zampini vpc = 4; 681b135d7daSStefano Zampini if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 6828135c375SStefano Zampini break; 6838135c375SStefano Zampini default: 6848135c375SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 6858135c375SStefano Zampini break; 6868135c375SStefano Zampini } 6878135c375SStefano Zampini break; 6888135c375SStefano Zampini case 3: 6898135c375SStefano Zampini switch (fpc) { 690b135d7daSStefano Zampini case 0: 6918135c375SStefano Zampini case 4: 6928135c375SStefano Zampini faces = fv3tet; 6938135c375SStefano Zampini vpf = 3; 69477eacf09SStefano Zampini vpc = 4; 695b135d7daSStefano Zampini if (fpcH == 5) { 696b135d7daSStefano Zampini facesH = fv3wedge; 697b135d7daSStefano Zampini vpfH = -4; 698b135d7daSStefano Zampini vpcH = 6; 699b135d7daSStefano Zampini } else if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 7008135c375SStefano Zampini break; 7018135c375SStefano Zampini case 6: 7028135c375SStefano Zampini faces = fv3hex; 7038135c375SStefano Zampini vpf = 4; 70477eacf09SStefano Zampini vpc = 8; 705b135d7daSStefano Zampini if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 7068135c375SStefano Zampini break; 7078135c375SStefano Zampini default: 7088135c375SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 7098135c375SStefano Zampini break; 7108135c375SStefano Zampini } 7118135c375SStefano Zampini break; 7128135c375SStefano Zampini default: 7138135c375SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 7148135c375SStefano Zampini break; 7158135c375SStefano Zampini } 7168135c375SStefano Zampini } 7178135c375SStefano Zampini ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 71877eacf09SStefano Zampini ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr); 71977eacf09SStefano Zampini ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr); 72077eacf09SStefano Zampini ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr); 72177eacf09SStefano Zampini ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 72277eacf09SStefano Zampini if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 72377eacf09SStefano Zampini ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr); 72477eacf09SStefano Zampini ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 72577eacf09SStefano Zampini ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr); 72677eacf09SStefano Zampini for (p=cStart;p<cEnd;p++) { 727b135d7daSStefano Zampini PetscInt dof, uvpc; 728b135d7daSStefano Zampini 72977eacf09SStefano Zampini ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr); 73077eacf09SStefano Zampini if (dof) { 73177eacf09SStefano Zampini PetscInt v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL; 73277eacf09SStefano Zampini PetscScalar *vals = NULL; 733b135d7daSStefano Zampini if (p < cEndInterior) uvpc = vpc; 734b135d7daSStefano Zampini else uvpc = vpcH; 73577eacf09SStefano Zampini if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim); 736b135d7daSStefano Zampini if (dof/sdim != uvpc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible number of cell dofs %D, vertices %D and space dim %D",dof/sdim,uvpc,sdim); 73777eacf09SStefano Zampini ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 73877eacf09SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 73977eacf09SStefano Zampini for (v=0;v<cellClosureSize;v++) 74077eacf09SStefano Zampini if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) { 74177eacf09SStefano Zampini vidxs = cellClosure + 2*v; 74277eacf09SStefano Zampini break; 74377eacf09SStefano Zampini } 74477eacf09SStefano Zampini if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); 745b135d7daSStefano Zampini for (v=0;v<uvpc;v++) { 74677eacf09SStefano Zampini PetscInt s; 747044a5661SStefano Zampini 74877eacf09SStefano Zampini for (s=0;s<sdim;s++) { 749b135d7daSStefano Zampini if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) { 75077eacf09SStefano Zampini ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); 75177eacf09SStefano Zampini } 75277eacf09SStefano Zampini } 75377eacf09SStefano Zampini } 75477eacf09SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 75577eacf09SStefano Zampini ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 75677eacf09SStefano Zampini } 75777eacf09SStefano Zampini } 75877eacf09SStefano Zampini if (dim > 1) { 759b135d7daSStefano Zampini PetscInt eEnd,eStart; 760044a5661SStefano Zampini 76177eacf09SStefano Zampini ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 76277eacf09SStefano Zampini for (p=eStart;p<eEnd;p++) { 76377eacf09SStefano Zampini const PetscInt *cone; 76477eacf09SStefano Zampini PetscInt coneSize,i; 76577eacf09SStefano Zampini PetscBool ispe = PETSC_TRUE; 76677eacf09SStefano Zampini 76777eacf09SStefano Zampini ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 76877eacf09SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 76977eacf09SStefano Zampini for (i=0;i<coneSize;i++) { 77077eacf09SStefano Zampini PetscInt v; 77177eacf09SStefano Zampini 77277eacf09SStefano Zampini ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 77377eacf09SStefano Zampini ispe = (PetscBool)(ispe && (v==2)); 77477eacf09SStefano Zampini } 77577eacf09SStefano Zampini if (ispe && coneSize) { 7763e96840aSStefano Zampini PetscInt ch, numChildren; 7773e96840aSStefano Zampini const PetscInt *children; 7783e96840aSStefano Zampini 77977eacf09SStefano Zampini ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 7803e96840aSStefano Zampini ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 7813e96840aSStefano Zampini for (ch = 0; ch < numChildren; ch++) { 7823e96840aSStefano Zampini ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 7833e96840aSStefano Zampini } 78477eacf09SStefano Zampini } 78577eacf09SStefano Zampini } 78677eacf09SStefano Zampini if (dim > 2) { 78777eacf09SStefano Zampini for (p=fStart;p<fEnd;p++) { 78877eacf09SStefano Zampini const PetscInt *cone; 78977eacf09SStefano Zampini PetscInt coneSize,i; 79077eacf09SStefano Zampini PetscBool ispe = PETSC_TRUE; 79177eacf09SStefano Zampini 79277eacf09SStefano Zampini ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 79377eacf09SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 79477eacf09SStefano Zampini for (i=0;i<coneSize;i++) { 79577eacf09SStefano Zampini PetscInt v; 79677eacf09SStefano Zampini 79777eacf09SStefano Zampini ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 79877eacf09SStefano Zampini ispe = (PetscBool)(ispe && (v==2)); 79977eacf09SStefano Zampini } 80077eacf09SStefano Zampini if (ispe && coneSize) { 8013e96840aSStefano Zampini PetscInt ch, numChildren; 8023e96840aSStefano Zampini const PetscInt *children; 8033e96840aSStefano Zampini 80477eacf09SStefano Zampini ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 8053e96840aSStefano Zampini ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 8063e96840aSStefano Zampini for (ch = 0; ch < numChildren; ch++) { 8073e96840aSStefano Zampini ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 8083e96840aSStefano Zampini } 80977eacf09SStefano Zampini } 81077eacf09SStefano Zampini } 81177eacf09SStefano Zampini } 81277eacf09SStefano Zampini } 81377eacf09SStefano Zampini } 81477eacf09SStefano Zampini for (p=fStart;p<fEnd;p++) { 81577eacf09SStefano Zampini const PetscInt *support; 8168135c375SStefano Zampini PetscInt supportSize; 81777eacf09SStefano Zampini PetscBool isbf = PETSC_FALSE; 8188135c375SStefano Zampini 8198135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 8203e6c54aaSStefano Zampini if (pown) { 8218135c375SStefano Zampini PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 82277eacf09SStefano Zampini PetscInt i; 82377eacf09SStefano Zampini 82477eacf09SStefano Zampini ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 82577eacf09SStefano Zampini for (i=0;i<supportSize;i++) { 82677eacf09SStefano Zampini if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE; 82777eacf09SStefano Zampini else has_ghost = PETSC_TRUE; 82877eacf09SStefano Zampini } 82977eacf09SStefano Zampini isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 83077eacf09SStefano Zampini } else { 83177eacf09SStefano Zampini isbf = (PetscBool)(supportSize == 1); 83277eacf09SStefano Zampini } 83377eacf09SStefano Zampini if (!isbf && perLabel) { 83477eacf09SStefano Zampini const PetscInt *cone; 83577eacf09SStefano Zampini PetscInt coneSize,i; 83677eacf09SStefano Zampini 83777eacf09SStefano Zampini ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 83877eacf09SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 83977eacf09SStefano Zampini isbf = PETSC_TRUE; 84077eacf09SStefano Zampini for (i=0;i<coneSize;i++) { 84177eacf09SStefano Zampini PetscInt v,d; 84277eacf09SStefano Zampini 84377eacf09SStefano Zampini ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 84477eacf09SStefano Zampini ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr); 84577eacf09SStefano Zampini isbf = (PetscBool)(isbf && v != d); 84677eacf09SStefano Zampini } 84777eacf09SStefano Zampini } 84877eacf09SStefano Zampini if (isbf) { 84977eacf09SStefano Zampini ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr); 85077eacf09SStefano Zampini } 85177eacf09SStefano Zampini } 85277eacf09SStefano Zampini /* count boundary faces */ 85377eacf09SStefano Zampini for (p=fStart,bf=0;p<fEnd;p++) { 85477eacf09SStefano Zampini if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 85577eacf09SStefano Zampini const PetscInt *support; 85677eacf09SStefano Zampini PetscInt supportSize,c; 8578135c375SStefano Zampini 8588135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 8598135c375SStefano Zampini ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 86077eacf09SStefano Zampini for (c=0;c<supportSize;c++) { 8613e96840aSStefano Zampini const PetscInt *cone; 862b135d7daSStefano Zampini PetscInt cell,cl,coneSize; 8633e96840aSStefano Zampini 8643e96840aSStefano Zampini cell = support[c]; 8653e96840aSStefano Zampini if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue; 8663e96840aSStefano Zampini ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 867b135d7daSStefano Zampini ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 868b135d7daSStefano Zampini for (cl=0;cl<coneSize;cl++) { 8693e96840aSStefano Zampini if (cone[cl] == p) { 8703e96840aSStefano Zampini bf += 1; 8713e96840aSStefano Zampini break; 8728135c375SStefano Zampini } 87377eacf09SStefano Zampini } 8743e96840aSStefano Zampini } 8758135c375SStefano Zampini } 8768135c375SStefano Zampini } 877*f86f7544SStefano Zampini minl = 1; 878*f86f7544SStefano Zampini label = NULL; 879*f86f7544SStefano Zampini if (enable_bmark) { 880*f86f7544SStefano Zampini PetscInt lminl = PETSC_MAX_INT; 881*f86f7544SStefano Zampini 8827bf4dd16SStefano Zampini ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr); 883*f86f7544SStefano Zampini if (label) { 884*f86f7544SStefano Zampini IS vals; 885*f86f7544SStefano Zampini PetscInt ldef; 886*f86f7544SStefano Zampini 887*f86f7544SStefano Zampini ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr); 888*f86f7544SStefano Zampini ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr); 889*f86f7544SStefano Zampini ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr); 890*f86f7544SStefano Zampini ierr = ISDestroy(&vals);CHKERRQ(ierr); 891*f86f7544SStefano Zampini lminl = PetscMin(ldef,lminl); 892*f86f7544SStefano Zampini } 893*f86f7544SStefano Zampini ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 894*f86f7544SStefano Zampini if (minl == PETSC_MAX_INT) minl = 1; 895*f86f7544SStefano Zampini } 8968135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr); 8978135c375SStefano Zampini for (p=fStart;p<fEnd;p++) { 89877eacf09SStefano Zampini if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 8998135c375SStefano Zampini const PetscInt *support; 90077eacf09SStefano Zampini PetscInt supportSize,c,nc = 0; 9018135c375SStefano Zampini 9028135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 9038135c375SStefano Zampini ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 9043e6c54aaSStefano Zampini if (pown) { 90577eacf09SStefano Zampini for (c=0;c<supportSize;c++) { 90677eacf09SStefano Zampini if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 90777eacf09SStefano Zampini fcells[nc++] = support[c]; 9088135c375SStefano Zampini } 90977eacf09SStefano Zampini } 91077eacf09SStefano Zampini } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c]; 91177eacf09SStefano Zampini for (c=0;c<nc;c++) { 91277eacf09SStefano Zampini const PetscInt *cone; 91377eacf09SStefano Zampini int vids[8]; 914b135d7daSStefano Zampini PetscInt i,coneSize,cell,cl,nv,cid = -1,mid = -1; 9158135c375SStefano Zampini 91677eacf09SStefano Zampini cell = fcells[c]; 91777eacf09SStefano Zampini ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 918b135d7daSStefano Zampini ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 919b135d7daSStefano Zampini for (cl=0;cl<coneSize;cl++) 92077eacf09SStefano Zampini if (cone[cl] == p) 9218135c375SStefano Zampini break; 922b135d7daSStefano Zampini if (cl == coneSize) continue; 9238135c375SStefano Zampini 92477eacf09SStefano Zampini /* face material id and type */ 925*f86f7544SStefano Zampini ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr); 9268135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 92777eacf09SStefano Zampini /* vertex ids */ 92877eacf09SStefano Zampini ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 929b135d7daSStefano Zampini if (cell >= cEndInterior) { 930b135d7daSStefano Zampini PetscInt nv = vpfH, inc = vpfH; 931b135d7daSStefano Zampini if (vpfH < 0) { /* Wedge */ 932b135d7daSStefano Zampini if (cl == 0 || cl == 1) nv = 3; 933b135d7daSStefano Zampini else nv = 4; 934b135d7daSStefano Zampini inc = -vpfH; 935b135d7daSStefano Zampini } 936b135d7daSStefano Zampini for (i=0;i<nv;i++) { 937b135d7daSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %d",vids[facesH[cl*inc+i]]);CHKERRQ(ierr); 938b135d7daSStefano Zampini } 939b135d7daSStefano Zampini } else { 9408135c375SStefano Zampini for (i=0;i<vpf;i++) { 941b56f74b9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d",vids[faces[cl*vpf+i]]);CHKERRQ(ierr); 9428135c375SStefano Zampini } 943b135d7daSStefano Zampini } 9448135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9453e96840aSStefano Zampini bf -= 1; 94677eacf09SStefano Zampini } 9478135c375SStefano Zampini } 9488135c375SStefano Zampini } 94977eacf09SStefano Zampini if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf); 95077eacf09SStefano Zampini ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr); 95177eacf09SStefano Zampini ierr = PetscFree(fcells);CHKERRQ(ierr); 9528135c375SStefano Zampini } 9538135c375SStefano Zampini 9548135c375SStefano Zampini /* mark owned vertices */ 9553e6c54aaSStefano Zampini vown = NULL; 9563e6c54aaSStefano Zampini if (pown) { 9573e6c54aaSStefano Zampini ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 9588135c375SStefano Zampini for (p=cStart;p<cEnd;p++) { 9598135c375SStefano Zampini PetscInt i,closureSize,*closure = NULL; 9608135c375SStefano Zampini 9613e6c54aaSStefano Zampini if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue; 9628135c375SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 9638135c375SStefano Zampini for (i=0;i<closureSize;i++) { 9648135c375SStefano Zampini const PetscInt pp = closure[2*i]; 9658135c375SStefano Zampini 9668135c375SStefano Zampini if (pp >= vStart && pp < vEnd) { 9673e6c54aaSStefano Zampini ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr); 9688135c375SStefano Zampini } 9698135c375SStefano Zampini } 9708135c375SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 9718135c375SStefano Zampini } 9728135c375SStefano Zampini } 9738135c375SStefano Zampini 9748135c375SStefano Zampini /* vertex_parents (Non-conforming meshes) */ 9758135c375SStefano Zampini parentSection = NULL; 9768135c375SStefano Zampini if (enable_ncmesh) { 9778135c375SStefano Zampini ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 9788135c375SStefano Zampini } 9798135c375SStefano Zampini if (parentSection) { 9808135c375SStefano Zampini PetscInt vp,gvp; 9818135c375SStefano Zampini 9828135c375SStefano Zampini for (vp=0,p=vStart;p<vEnd;p++) { 9838135c375SStefano Zampini DMLabel dlabel; 9848135c375SStefano Zampini PetscInt parent,depth; 9858135c375SStefano Zampini 9863e6c54aaSStefano Zampini if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 9878135c375SStefano Zampini ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 9888135c375SStefano Zampini ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 9898135c375SStefano Zampini ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 9908135c375SStefano Zampini if (parent != p) vp++; 9918135c375SStefano Zampini } 9928135c375SStefano Zampini ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 9938135c375SStefano Zampini if (gvp) { 9948135c375SStefano Zampini PetscInt maxsupp; 9958135c375SStefano Zampini PetscBool *skip = NULL; 9968135c375SStefano Zampini 9978135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr); 9988135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr); 9998135c375SStefano Zampini ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr); 10008135c375SStefano Zampini ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr); 10018135c375SStefano Zampini for (p=vStart;p<vEnd;p++) { 10028135c375SStefano Zampini DMLabel dlabel; 10038135c375SStefano Zampini PetscInt parent; 10048135c375SStefano Zampini 10053e6c54aaSStefano Zampini if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 10068135c375SStefano Zampini ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 10078135c375SStefano Zampini ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 10088135c375SStefano Zampini if (parent != p) { 1009d42aff11SStefano Zampini int vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */ 10103924b612SStefano Zampini PetscInt i,nv,ssize,n,numChildren,depth = -1; 10118135c375SStefano Zampini const PetscInt *children; 10123924b612SStefano Zampini 10133924b612SStefano Zampini ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr); 10143924b612SStefano Zampini switch (ssize) { 10158135c375SStefano Zampini case 2: /* edge */ 10168135c375SStefano Zampini nv = 0; 1017cc0d3ed7SStefano Zampini ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 101877eacf09SStefano Zampini ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 10198135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr); 10208135c375SStefano Zampini for (i=0;i<nv;i++) { 1021bfadf1d8SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr); 10228135c375SStefano Zampini } 10238135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 10248135c375SStefano Zampini vp--; 10258135c375SStefano Zampini break; 10268135c375SStefano Zampini case 4: /* face */ 10278135c375SStefano Zampini ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 10288135c375SStefano Zampini for (n=0;n<numChildren;n++) { 10298135c375SStefano Zampini ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr); 10308135c375SStefano Zampini if (!depth) { 10318135c375SStefano Zampini const PetscInt *hvsupp,*hesupp,*cone; 10328135c375SStefano Zampini PetscInt hvsuppSize,hesuppSize,coneSize; 1033451a39c7SStefano Zampini PetscInt hv = children[n],he = -1,f; 10348135c375SStefano Zampini 1035451a39c7SStefano Zampini ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr); 10368135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr); 10378135c375SStefano Zampini ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr); 10388135c375SStefano Zampini for (i=0;i<hvsuppSize;i++) { 10398135c375SStefano Zampini PetscInt ep; 10408135c375SStefano Zampini ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr); 10418135c375SStefano Zampini if (ep != hvsupp[i]) { 10428135c375SStefano Zampini he = hvsupp[i]; 10438135c375SStefano Zampini } else { 10448135c375SStefano Zampini skip[i] = PETSC_TRUE; 10458135c375SStefano Zampini } 10468135c375SStefano Zampini } 1047451a39c7SStefano Zampini if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize); 10488135c375SStefano Zampini ierr = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr); 104911a4995dSStefano Zampini vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 10508135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr); 10518135c375SStefano Zampini ierr = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr); 10528135c375SStefano Zampini for (f=0;f<hesuppSize;f++) { 10538135c375SStefano Zampini PetscInt j; 10548135c375SStefano Zampini 10558135c375SStefano Zampini ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr); 10568135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr); 10578135c375SStefano Zampini for (j=0;j<coneSize;j++) { 10588135c375SStefano Zampini PetscInt k; 10598135c375SStefano Zampini for (k=0;k<hvsuppSize;k++) { 10608135c375SStefano Zampini if (hvsupp[k] == cone[j]) { 10618135c375SStefano Zampini skip[k] = PETSC_TRUE; 10628135c375SStefano Zampini break; 10638135c375SStefano Zampini } 10648135c375SStefano Zampini } 10658135c375SStefano Zampini } 10668135c375SStefano Zampini } 10678135c375SStefano Zampini for (i=0;i<hvsuppSize;i++) { 10688135c375SStefano Zampini if (!skip[i]) { 10698135c375SStefano Zampini ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr); 107011a4995dSStefano Zampini vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 10718135c375SStefano Zampini } 10728135c375SStefano Zampini } 10738135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr); 10748135c375SStefano Zampini for (i=0;i<2;i++) { 1075bfadf1d8SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)(vids[i]-vStart));CHKERRQ(ierr); 10768135c375SStefano Zampini } 10778135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 10788135c375SStefano Zampini vp--; 10798135c375SStefano Zampini } 10808135c375SStefano Zampini } 10818135c375SStefano Zampini break; 10828135c375SStefano Zampini default: 10833924b612SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize); 10848135c375SStefano Zampini } 10858135c375SStefano Zampini } 10868135c375SStefano Zampini } 10878135c375SStefano Zampini ierr = PetscFree(skip);CHKERRQ(ierr); 10888135c375SStefano Zampini } 10898135c375SStefano Zampini if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp); 10908135c375SStefano Zampini } 10918135c375SStefano Zampini ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 10923e6c54aaSStefano Zampini ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 10938135c375SStefano Zampini 10948135c375SStefano Zampini /* vertices */ 109577eacf09SStefano Zampini if (hovec) { /* higher-order meshes */ 109677eacf09SStefano Zampini const char *fec; 10970286d493SLisandro Dalcin PetscInt i,n,s; 109877eacf09SStefano Zampini 109977eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 110077eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr); 110177eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 110277eacf09SStefano Zampini ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr); 110377eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 11043e7cab22SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr); 110577eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 110677eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 110777eacf09SStefano Zampini ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr); 110877eacf09SStefano Zampini ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr); 110977eacf09SStefano Zampini if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim); 111077eacf09SStefano Zampini for (i=0;i<n/sdim;i++) { 111177eacf09SStefano Zampini for (s=0;s<sdim;s++) { 111277eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr); 111377eacf09SStefano Zampini } 111477eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 111577eacf09SStefano Zampini } 111677eacf09SStefano Zampini ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr); 111777eacf09SStefano Zampini } else { 1118cc0d3ed7SStefano Zampini ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr); 11198135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1120cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr); 11218135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 11228135c375SStefano Zampini ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr); 1123cc0d3ed7SStefano Zampini for (p=0;p<nvert/sdim;p++) { 1124cc0d3ed7SStefano Zampini PetscInt s; 1125cc0d3ed7SStefano Zampini for (s=0;s<sdim;s++) { 11263e96840aSStefano Zampini PetscReal v = PetscRealPart(array[p*sdim+s]); 11273e96840aSStefano Zampini 11283e96840aSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);CHKERRQ(ierr); 11298135c375SStefano Zampini } 1130cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 11318135c375SStefano Zampini } 11328135c375SStefano Zampini ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr); 113377eacf09SStefano Zampini } 11348135c375SStefano Zampini PetscFunctionReturn(0); 11358135c375SStefano Zampini } 11368135c375SStefano Zampini 11370286d493SLisandro Dalcin PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 11388135c375SStefano Zampini { 11398135c375SStefano Zampini PetscErrorCode ierr; 11408135c375SStefano Zampini PetscFunctionBegin; 11410286d493SLisandro Dalcin ierr = DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);CHKERRQ(ierr); 11428135c375SStefano Zampini PetscFunctionReturn(0); 11438135c375SStefano Zampini } 1144