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 14d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 15d71ae5a4SJacob Faibussowitsch { 168135c375SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 178135c375SStefano Zampini PetscInt i; 188135c375SStefano Zampini 198135c375SStefano Zampini PetscFunctionBegin; 2048a46eb9SPierre Jolivet for (i = 0; i < ctx->nf; i++) PetscCall(VecScatterDestroy(&ctx->scctx[i])); 219566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx->scctx)); 229566063dSJacob Faibussowitsch PetscCall(PetscFree(vctx)); 233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 248135c375SStefano Zampini } 258135c375SStefano Zampini 26d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 27d71ae5a4SJacob Faibussowitsch { 288135c375SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 298135c375SStefano Zampini PetscInt f; 308135c375SStefano Zampini 318135c375SStefano Zampini PetscFunctionBegin; 328135c375SStefano Zampini for (f = 0; f < nf; f++) { 339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx->scctx[f], (Vec)oX, (Vec)oXfield[f], INSERT_VALUES, SCATTER_FORWARD)); 349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx->scctx[f], (Vec)oX, (Vec)oXfield[f], INSERT_VALUES, SCATTER_FORWARD)); 358135c375SStefano Zampini } 363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 378135c375SStefano Zampini } 388135c375SStefano Zampini 398135c375SStefano Zampini /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */ 40d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer) 41d71ae5a4SJacob Faibussowitsch { 428135c375SStefano Zampini DM dm = (DM)odm; 434cac2994SStefano Zampini Vec xlocal, xfield, *Ufield; 448135c375SStefano Zampini PetscDS ds; 458135c375SStefano Zampini IS globalNum, isfield; 468135c375SStefano Zampini PetscBT vown; 478135c375SStefano Zampini char **fieldname = NULL, **fec_type = NULL; 488135c375SStefano Zampini const PetscInt *gNum; 49bb77a09fSStefano Zampini PetscInt *nlocal, *bs, *idxs, *dims; 508135c375SStefano Zampini PetscInt f, maxfields, nfields, c, totc, totdofs, Nv, cum, i; 51b135d7daSStefano Zampini PetscInt dim, cStart, cEnd, vStart, vEnd; 528135c375SStefano Zampini GLVisViewerCtx *ctx; 538135c375SStefano Zampini PetscSection s; 548135c375SStefano Zampini 558135c375SStefano Zampini PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 579566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 589566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "_glvis_plex_gnum", (PetscObject *)&globalNum)); 60b135d7daSStefano Zampini if (!globalNum) { 61*e2b8d0fcSMatthew G. Knepley PetscCall(DMPlexCreateCellNumbering(dm, PETSC_TRUE, &globalNum)); 629566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "_glvis_plex_gnum", (PetscObject)globalNum)); 639566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)globalNum)); 64b135d7daSStefano Zampini } 659566063dSJacob Faibussowitsch PetscCall(ISGetIndices(globalNum, &gNum)); 669566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(vEnd - vStart, &vown)); 678135c375SStefano Zampini for (c = cStart, totc = 0; c < cEnd; c++) { 688135c375SStefano Zampini if (gNum[c - cStart] >= 0) { 698135c375SStefano Zampini PetscInt i, numPoints, *points = NULL; 708135c375SStefano Zampini 718135c375SStefano Zampini totc++; 729566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &numPoints, &points)); 738135c375SStefano Zampini for (i = 0; i < numPoints * 2; i += 2) { 7448a46eb9SPierre Jolivet if ((points[i] >= vStart) && (points[i] < vEnd)) PetscCall(PetscBTSet(vown, points[i] - vStart)); 758135c375SStefano Zampini } 769566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &numPoints, &points)); 778135c375SStefano Zampini } 788135c375SStefano Zampini } 799371c9d4SSatish Balay for (f = 0, Nv = 0; f < vEnd - vStart; f++) 809371c9d4SSatish Balay if (PetscLikely(PetscBTLookup(vown, f))) Nv++; 818135c375SStefano Zampini 829566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &xlocal)); 839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xlocal, &totdofs)); 849566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &nfields)); 868135c375SStefano Zampini for (f = 0, maxfields = 0; f < nfields; f++) { 878135c375SStefano Zampini PetscInt bs; 888135c375SStefano Zampini 899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(s, f, &bs)); 908135c375SStefano Zampini maxfields += bs; 918135c375SStefano Zampini } 929566063dSJacob Faibussowitsch PetscCall(PetscCalloc7(maxfields, &fieldname, maxfields, &nlocal, maxfields, &bs, maxfields, &dims, maxfields, &fec_type, totdofs, &idxs, maxfields, &Ufield)); 939566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 949566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(maxfields, &ctx->scctx)); 959566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 968135c375SStefano Zampini if (ds) { 978135c375SStefano Zampini for (f = 0; f < nfields; f++) { 988135c375SStefano Zampini const char *fname; 998135c375SStefano Zampini char name[256]; 1008135c375SStefano Zampini PetscObject disc; 1018135c375SStefano Zampini size_t len; 1028135c375SStefano Zampini 1039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldName(s, f, &fname)); 1049566063dSJacob Faibussowitsch PetscCall(PetscStrlen(fname, &len)); 1058135c375SStefano Zampini if (len) { 106c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(name, fname, sizeof(name))); 1078135c375SStefano Zampini } else { 10863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, 256, "Field%" PetscInt_FMT, f)); 1098135c375SStefano Zampini } 1109566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, &disc)); 1118135c375SStefano Zampini if (disc) { 1128135c375SStefano Zampini PetscClassId id; 1138135c375SStefano Zampini PetscInt Nc; 1148135c375SStefano Zampini char fec[64]; 1158135c375SStefano Zampini 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(disc, &id)); 1178135c375SStefano Zampini if (id == PETSCFE_CLASSID) { 1188135c375SStefano Zampini PetscFE fem = (PetscFE)disc; 1198135c375SStefano Zampini PetscDualSpace sp; 1208135c375SStefano Zampini PetscDualSpaceType spname; 1218135c375SStefano Zampini PetscInt order; 1228135c375SStefano Zampini PetscBool islag, continuous, H1 = PETSC_TRUE; 1238135c375SStefano Zampini 1249566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fem, &Nc)); 1259566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fem, &sp)); 1269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetType(sp, &spname)); 1279566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(spname, PETSCDUALSPACELAGRANGE, &islag)); 12828b400f6SJacob Faibussowitsch PetscCheck(islag, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported dual space"); 1299566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &continuous)); 1309566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetOrder(sp, &order)); 13128d58a37SPierre Jolivet if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */ 13263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec, 64, "FiniteElementCollection: H1_%" PetscInt_FMT "D_P1", dim)); 1338135c375SStefano Zampini } else { 13463a3b9bcSJacob Faibussowitsch PetscCheck(continuous || !order, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Discontinuous space visualization currently unsupported for order %" PetscInt_FMT, order); 1358135c375SStefano Zampini H1 = PETSC_FALSE; 13663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec, 64, "FiniteElementCollection: L2_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, order)); 1378135c375SStefano Zampini } 1389566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &fieldname[ctx->nf])); 1398135c375SStefano Zampini bs[ctx->nf] = Nc; 140bb77a09fSStefano Zampini dims[ctx->nf] = dim; 1418135c375SStefano Zampini if (H1) { 1428135c375SStefano Zampini nlocal[ctx->nf] = Nc * Nv; 1439566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec, &fec_type[ctx->nf])); 1449566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nv * Nc, &xfield)); 1458135c375SStefano Zampini for (i = 0, cum = 0; i < vEnd - vStart; i++) { 1468135c375SStefano Zampini PetscInt j, off; 1478135c375SStefano Zampini 1488135c375SStefano Zampini if (PetscUnlikely(!PetscBTLookup(vown, i))) continue; 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(s, i + vStart, f, &off)); 1508135c375SStefano Zampini for (j = 0; j < Nc; j++) idxs[cum++] = off + j; 1518135c375SStefano Zampini } 1529566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal), Nv * Nc, idxs, PETSC_USE_POINTER, &isfield)); 1538135c375SStefano Zampini } else { 1548135c375SStefano Zampini nlocal[ctx->nf] = Nc * totc; 1559566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec, &fec_type[ctx->nf])); 1569566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nc * totc, &xfield)); 1578135c375SStefano Zampini for (i = 0, cum = 0; i < cEnd - cStart; i++) { 1588135c375SStefano Zampini PetscInt j, off; 1598135c375SStefano Zampini 1608135c375SStefano Zampini if (PetscUnlikely(gNum[i] < 0)) continue; 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(s, i + cStart, f, &off)); 1628135c375SStefano Zampini for (j = 0; j < Nc; j++) idxs[cum++] = off + j; 1638135c375SStefano Zampini } 1649566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal), totc * Nc, idxs, PETSC_USE_POINTER, &isfield)); 1658135c375SStefano Zampini } 1669566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(xlocal, isfield, xfield, NULL, &ctx->scctx[ctx->nf])); 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xfield)); 1689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfield)); 1698135c375SStefano Zampini ctx->nf++; 1708135c375SStefano Zampini } else if (id == PETSCFV_CLASSID) { 1718135c375SStefano Zampini PetscInt c; 1728135c375SStefano Zampini 1739566063dSJacob Faibussowitsch PetscCall(PetscFVGetNumComponents((PetscFV)disc, &Nc)); 17463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec, 64, "FiniteElementCollection: L2_%" PetscInt_FMT "D_P0", dim)); 1758135c375SStefano Zampini for (c = 0; c < Nc; c++) { 1768135c375SStefano Zampini char comp[256]; 17763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(comp, 256, "%s-Comp%" PetscInt_FMT, name, c)); 1789566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(comp, &fieldname[ctx->nf])); 1798135c375SStefano Zampini bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ 1808135c375SStefano Zampini nlocal[ctx->nf] = totc; 181bb77a09fSStefano Zampini dims[ctx->nf] = dim; 1829566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec, &fec_type[ctx->nf])); 1839566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, totc, &xfield)); 1848135c375SStefano Zampini for (i = 0, cum = 0; i < cEnd - cStart; i++) { 1858135c375SStefano Zampini PetscInt off; 1868135c375SStefano Zampini 1878135c375SStefano Zampini if (PetscUnlikely(gNum[i]) < 0) continue; 1889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(s, i + cStart, f, &off)); 1898135c375SStefano Zampini idxs[cum++] = off + c; 1908135c375SStefano Zampini } 1919566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)xlocal), totc, idxs, PETSC_USE_POINTER, &isfield)); 1929566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(xlocal, isfield, xfield, NULL, &ctx->scctx[ctx->nf])); 1939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xfield)); 1949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfield)); 1958135c375SStefano Zampini ctx->nf++; 1968135c375SStefano Zampini } 19763a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f); 19863a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Missing discretization for field %" PetscInt_FMT, f); 1998135c375SStefano Zampini } 2008135c375SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Needs a DS attached to the DM"); 2019566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&vown)); 2029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xlocal)); 2039566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(globalNum, &gNum)); 2048135c375SStefano Zampini 2054cac2994SStefano Zampini /* create work vectors */ 2064cac2994SStefano Zampini for (f = 0; f < ctx->nf; f++) { 2079566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm), nlocal[f], PETSC_DECIDE, &Ufield[f])); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Ufield[f], fieldname[f])); 2099566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(Ufield[f], bs[f])); 2109566063dSJacob Faibussowitsch PetscCall(VecSetDM(Ufield[f], dm)); 2114cac2994SStefano Zampini } 2124cac2994SStefano Zampini 2138135c375SStefano Zampini /* customize the viewer */ 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, ctx->nf, (const char **)fec_type, dims, DMPlexSampleGLVisFields_Private, (PetscObject *)Ufield, ctx, DestroyGLVisViewerCtx_Private)); 2158135c375SStefano Zampini for (f = 0; f < ctx->nf; f++) { 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldname[f])); 2179566063dSJacob Faibussowitsch PetscCall(PetscFree(fec_type[f])); 2189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ufield[f])); 2198135c375SStefano Zampini } 2209566063dSJacob Faibussowitsch PetscCall(PetscFree7(fieldname, nlocal, bs, dims, fec_type, idxs, Ufield)); 2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2228135c375SStefano Zampini } 2238135c375SStefano Zampini 2249371c9d4SSatish Balay typedef enum { 2259371c9d4SSatish Balay MFEM_POINT = 0, 2269371c9d4SSatish Balay MFEM_SEGMENT, 2279371c9d4SSatish Balay MFEM_TRIANGLE, 2289371c9d4SSatish Balay MFEM_SQUARE, 2299371c9d4SSatish Balay MFEM_TETRAHEDRON, 2309371c9d4SSatish Balay MFEM_CUBE, 2319371c9d4SSatish Balay MFEM_PRISM, 2329371c9d4SSatish Balay MFEM_UNDEF 2339371c9d4SSatish Balay } MFEM_cid; 2348135c375SStefano Zampini 2359371c9d4SSatish Balay MFEM_cid mfem_table_cid[4][7] = { 2369371c9d4SSatish Balay {MFEM_POINT, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF}, 2378135c375SStefano Zampini {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF}, 2388135c375SStefano Zampini {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_TRIANGLE, MFEM_SQUARE, MFEM_UNDEF, MFEM_UNDEF}, 2399371c9d4SSatish Balay {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF, MFEM_TETRAHEDRON, MFEM_PRISM, MFEM_CUBE } 2409371c9d4SSatish Balay }; 2418135c375SStefano Zampini 2429371c9d4SSatish Balay MFEM_cid mfem_table_cid_unint[4][9] = { 2439371c9d4SSatish Balay {MFEM_POINT, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_UNDEF}, 244b135d7daSStefano Zampini {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF, MFEM_UNDEF, MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_UNDEF}, 245b135d7daSStefano Zampini {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_TRIANGLE, MFEM_SQUARE, MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_UNDEF}, 2469371c9d4SSatish Balay {MFEM_POINT, MFEM_UNDEF, MFEM_SEGMENT, MFEM_UNDEF, MFEM_TETRAHEDRON, MFEM_UNDEF, MFEM_PRISM, MFEM_UNDEF, MFEM_CUBE } 2479371c9d4SSatish Balay }; 248044a5661SStefano Zampini 249d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid) 250d71ae5a4SJacob Faibussowitsch { 2518135c375SStefano Zampini DMLabel dlabel; 252044a5661SStefano Zampini PetscInt depth, csize, pdepth, dim; 2538135c375SStefano Zampini 2548135c375SStefano Zampini PetscFunctionBegin; 2559566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &dlabel)); 2569566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(dlabel, p, &pdepth)); 2579566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &csize)); 2589566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2608135c375SStefano Zampini if (label) { 2619566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, p, mid)); 262f86f7544SStefano Zampini *mid = *mid - minl + 1; /* MFEM does not like negative markers */ 26377eacf09SStefano Zampini } else *mid = 1; 264044a5661SStefano Zampini if (depth >= 0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 2651dca8a05SBarry Smith PetscCheck(dim >= 0 && dim <= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %" PetscInt_FMT, dim); 26663a3b9bcSJacob Faibussowitsch PetscCheck(csize <= 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "Found cone size %" PetscInt_FMT " for point %" PetscInt_FMT, csize, p); 26763a3b9bcSJacob Faibussowitsch PetscCheck(depth == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Found depth %" PetscInt_FMT " for point %" PetscInt_FMT ". You should interpolate the mesh first", depth, p); 268044a5661SStefano Zampini *cid = mfem_table_cid_unint[dim][csize]; 269044a5661SStefano Zampini } else { 27063a3b9bcSJacob Faibussowitsch PetscCheck(csize <= 6, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cone size %" PetscInt_FMT " for point %" PetscInt_FMT, csize, p); 2711dca8a05SBarry Smith PetscCheck(pdepth >= 0 && pdepth <= 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Depth %" PetscInt_FMT " for point %" PetscInt_FMT, csize, p); 272044a5661SStefano Zampini *cid = mfem_table_cid[pdepth][csize]; 273044a5661SStefano Zampini } 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2758135c375SStefano Zampini } 2768135c375SStefano Zampini 277d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[]) 278d71ae5a4SJacob Faibussowitsch { 279cc0d3ed7SStefano Zampini PetscInt dim, sdim, dof = 0, off = 0, i, q, vStart, vEnd, numPoints, *points = NULL; 2808135c375SStefano Zampini 2818135c375SStefano Zampini PetscFunctionBegin; 2829566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 2839566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 284cc0d3ed7SStefano Zampini sdim = dim; 285cc0d3ed7SStefano Zampini if (csec) { 28684f354e3SLisandro Dalcin PetscInt sStart, sEnd; 28784f354e3SLisandro Dalcin 2889566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &sdim)); 2899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(csec, &sStart, &sEnd)); 2909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(csec, vStart, &off)); 291cc0d3ed7SStefano Zampini off = off / sdim; 29248a46eb9SPierre Jolivet if (p >= sStart && p < sEnd) PetscCall(PetscSectionGetDof(csec, p, &dof)); 29384f354e3SLisandro Dalcin } 294cc0d3ed7SStefano Zampini if (!dof) { 2959566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &numPoints, &points)); 2968135c375SStefano Zampini for (i = 0, q = 0; i < numPoints * 2; i += 2) 2979371c9d4SSatish Balay if ((points[i] >= vStart) && (points[i] < vEnd)) vids[q++] = points[i] - vStart + off; 2989566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &numPoints, &points)); 299cc0d3ed7SStefano Zampini } else { 3009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(csec, p, &off)); 3019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(csec, p, &dof)); 30296ca5757SLisandro Dalcin for (q = 0; q < dof / sdim; q++) vids[q] = off / sdim + q; 303cc0d3ed7SStefano Zampini } 3048135c375SStefano Zampini *nv = q; 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3068135c375SStefano Zampini } 3078135c375SStefano Zampini 308d71ae5a4SJacob Faibussowitsch static PetscErrorCode GLVisCreateFE(PetscFE femIn, char name[32], PetscFE *fem, IS *perm) 309d71ae5a4SJacob Faibussowitsch { 310066ea43fSLisandro Dalcin DM K; 311066ea43fSLisandro Dalcin PetscSpace P; 312066ea43fSLisandro Dalcin PetscDualSpace Q; 313066ea43fSLisandro Dalcin PetscQuadrature q, fq; 314066ea43fSLisandro Dalcin PetscInt dim, deg, dof; 315066ea43fSLisandro Dalcin DMPolytopeType ptype; 316066ea43fSLisandro Dalcin PetscBool isSimplex, isTensor; 317066ea43fSLisandro Dalcin PetscBool continuity = PETSC_FALSE; 318066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI; 319066ea43fSLisandro Dalcin PetscBool endpoint = PETSC_TRUE; 320066ea43fSLisandro Dalcin MPI_Comm comm; 321066ea43fSLisandro Dalcin 322066ea43fSLisandro Dalcin PetscFunctionBegin; 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)femIn, &comm)); 3249566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(femIn, &P)); 3259566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(femIn, &Q)); 3269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(Q, &K)); 3279566063dSJacob Faibussowitsch PetscCall(DMGetDimension(K, &dim)); 3289566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, °, NULL)); 3299566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(P, &dof)); 3309566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(K, 0, &ptype)); 331066ea43fSLisandro Dalcin switch (ptype) { 332066ea43fSLisandro Dalcin case DM_POLYTOPE_QUADRILATERAL: 333d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 334d71ae5a4SJacob Faibussowitsch isSimplex = PETSC_FALSE; 335d71ae5a4SJacob Faibussowitsch break; 336d71ae5a4SJacob Faibussowitsch default: 337d71ae5a4SJacob Faibussowitsch isSimplex = PETSC_TRUE; 338d71ae5a4SJacob Faibussowitsch break; 339066ea43fSLisandro Dalcin } 340066ea43fSLisandro Dalcin isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 3410c2bc6bfSStefano Zampini if (isSimplex) deg = PetscMin(deg, 3); /* Permutation not coded for degree higher than 3 */ 342066ea43fSLisandro Dalcin /* Create space */ 3439566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(comm, &P)); 3449566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 3459566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(P, isTensor)); 3469566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(P, dof)); 3479566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(P, dim)); 3489566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(P, deg, PETSC_DETERMINE)); 3499566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(P)); 350066ea43fSLisandro Dalcin /* Create dual space */ 3519566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(comm, &Q)); 3529566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 3539566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(Q, isTensor)); 3549566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetContinuity(Q, continuity)); 3559566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0)); 3569566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(Q, dof)); 3579566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(Q, deg)); 3589566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); 3599566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(Q, K)); 3609566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 3619566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(Q)); 362066ea43fSLisandro Dalcin /* Create quadrature */ 363066ea43fSLisandro Dalcin if (isSimplex) { 3649566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 1, -1, +1, &q)); 3659566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, deg + 1, -1, +1, &fq)); 366066ea43fSLisandro Dalcin } else { 3679566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, deg + 1, -1, +1, &q)); 3689566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, deg + 1, -1, +1, &fq)); 369066ea43fSLisandro Dalcin } 370066ea43fSLisandro Dalcin /* Create finite element */ 3719566063dSJacob Faibussowitsch PetscCall(PetscFECreate(comm, fem)); 37263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, 32, "L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, deg)); 3739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*fem, name)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 3759566063dSJacob Faibussowitsch PetscCall(PetscFESetNumComponents(*fem, dof)); 3769566063dSJacob Faibussowitsch PetscCall(PetscFESetBasisSpace(*fem, P)); 3779566063dSJacob Faibussowitsch PetscCall(PetscFESetDualSpace(*fem, Q)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(*fem, q)); 3799566063dSJacob Faibussowitsch PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 3809566063dSJacob Faibussowitsch PetscCall(PetscFESetUp(*fem)); 3810c2bc6bfSStefano Zampini 3820c2bc6bfSStefano Zampini /* Both MFEM and PETSc are lexicographic, but PLEX stores the swapped cone */ 3830c2bc6bfSStefano Zampini *perm = NULL; 3840c2bc6bfSStefano Zampini if (isSimplex && dim == 3) { 3850c2bc6bfSStefano Zampini PetscInt celldofs, *pidx; 3860c2bc6bfSStefano Zampini 3870c2bc6bfSStefano Zampini PetscCall(PetscDualSpaceGetDimension(Q, &celldofs)); 3880c2bc6bfSStefano Zampini celldofs /= dof; 3890c2bc6bfSStefano Zampini PetscCall(PetscMalloc1(celldofs, &pidx)); 3900c2bc6bfSStefano Zampini switch (celldofs) { 3910c2bc6bfSStefano Zampini case 4: 3920c2bc6bfSStefano Zampini pidx[0] = 2; 3930c2bc6bfSStefano Zampini pidx[1] = 0; 3940c2bc6bfSStefano Zampini pidx[2] = 1; 3950c2bc6bfSStefano Zampini pidx[3] = 3; 3960c2bc6bfSStefano Zampini break; 3970c2bc6bfSStefano Zampini case 10: 3980c2bc6bfSStefano Zampini pidx[0] = 5; 3990c2bc6bfSStefano Zampini pidx[1] = 3; 4000c2bc6bfSStefano Zampini pidx[2] = 0; 4010c2bc6bfSStefano Zampini pidx[3] = 4; 4020c2bc6bfSStefano Zampini pidx[4] = 1; 4030c2bc6bfSStefano Zampini pidx[5] = 2; 4040c2bc6bfSStefano Zampini pidx[6] = 8; 4050c2bc6bfSStefano Zampini pidx[7] = 6; 4060c2bc6bfSStefano Zampini pidx[8] = 7; 4070c2bc6bfSStefano Zampini pidx[9] = 9; 4080c2bc6bfSStefano Zampini break; 4090c2bc6bfSStefano Zampini case 20: 4100c2bc6bfSStefano Zampini pidx[0] = 9; 4110c2bc6bfSStefano Zampini pidx[1] = 7; 4120c2bc6bfSStefano Zampini pidx[2] = 4; 4130c2bc6bfSStefano Zampini pidx[3] = 0; 4140c2bc6bfSStefano Zampini pidx[4] = 8; 4150c2bc6bfSStefano Zampini pidx[5] = 5; 4160c2bc6bfSStefano Zampini pidx[6] = 1; 4170c2bc6bfSStefano Zampini pidx[7] = 6; 4180c2bc6bfSStefano Zampini pidx[8] = 2; 4190c2bc6bfSStefano Zampini pidx[9] = 3; 4200c2bc6bfSStefano Zampini pidx[10] = 15; 4210c2bc6bfSStefano Zampini pidx[11] = 13; 4220c2bc6bfSStefano Zampini pidx[12] = 10; 4230c2bc6bfSStefano Zampini pidx[13] = 14; 4240c2bc6bfSStefano Zampini pidx[14] = 11; 4250c2bc6bfSStefano Zampini pidx[15] = 12; 4260c2bc6bfSStefano Zampini pidx[16] = 18; 4270c2bc6bfSStefano Zampini pidx[17] = 16; 4280c2bc6bfSStefano Zampini pidx[18] = 17; 4290c2bc6bfSStefano Zampini pidx[19] = 19; 4300c2bc6bfSStefano Zampini break; 431d71ae5a4SJacob Faibussowitsch default: 432d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_SUP, "Unhandled degree,dof pair %" PetscInt_FMT ",%" PetscInt_FMT, deg, celldofs); 433d71ae5a4SJacob Faibussowitsch break; 4340c2bc6bfSStefano Zampini } 4350c2bc6bfSStefano Zampini PetscCall(ISCreateBlock(PETSC_COMM_SELF, dof, celldofs, pidx, PETSC_OWN_POINTER, perm)); 4360c2bc6bfSStefano Zampini } 4370c2bc6bfSStefano Zampini 438066ea43fSLisandro Dalcin /* Cleanup */ 4399566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&P)); 4409566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&Q)); 4419566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 4429566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&fq)); 4433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 444066ea43fSLisandro Dalcin } 445066ea43fSLisandro Dalcin 44677eacf09SStefano Zampini /* 44777eacf09SStefano Zampini ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 44877eacf09SStefano Zampini Higher order meshes are also supported 44977eacf09SStefano Zampini */ 450d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 451d71ae5a4SJacob Faibussowitsch { 4528135c375SStefano Zampini DMLabel label; 4536858538eSMatthew G. Knepley PetscSection coordSection, coordSectionCell, parentSection, hoSection = NULL; 4546858538eSMatthew G. Knepley Vec coordinates, coordinatesCell, hovec; 4558135c375SStefano Zampini const PetscScalar *array; 456f86f7544SStefano Zampini PetscInt bf, p, sdim, dim, depth, novl, minl; 457412e9a14SMatthew G. Knepley PetscInt cStart, cEnd, vStart, vEnd, nvert; 4583924b612SStefano Zampini PetscMPIInt size; 4593e6c54aaSStefano Zampini PetscBool localized, isascii; 46028d58a37SPierre Jolivet PetscBool enable_mfem, enable_boundary, enable_ncmesh, view_ovl = PETSC_FALSE; 4613e6c54aaSStefano Zampini PetscBT pown, vown; 4628135c375SStefano Zampini PetscContainer glvis_container; 4636858538eSMatthew G. Knepley PetscBool cellvertex = PETSC_FALSE, enabled = PETSC_TRUE; 464f86f7544SStefano Zampini PetscBool enable_emark, enable_bmark; 46577eacf09SStefano Zampini const char *fmt; 4667bf4dd16SStefano Zampini char emark[64] = "", bmark[64] = ""; 4678135c375SStefano Zampini 4688135c375SStefano Zampini PetscFunctionBegin; 4698135c375SStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4708135c375SStefano Zampini PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 4719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 47228b400f6SJacob Faibussowitsch PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer must be of type VIEWERASCII"); 4739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size)); 47408401ef6SPierre Jolivet PetscCheck(size <= 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Use single sequential viewers for parallel visualization"); 4759566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4760c2bc6bfSStefano Zampini PetscCall(DMPlexGetDepth(dm, &depth)); 4778135c375SStefano Zampini 4788135c375SStefano Zampini /* get container: determines if a process visualizes is portion of the data or not */ 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container)); 48028b400f6SJacob Faibussowitsch PetscCheck(glvis_container, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis container"); 4818135c375SStefano Zampini { 4828135c375SStefano Zampini PetscViewerGLVisInfo glvis_info; 4839566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info)); 4848135c375SStefano Zampini enabled = glvis_info->enabled; 48577eacf09SStefano Zampini fmt = glvis_info->fmt; 4868135c375SStefano Zampini } 48721414b21SStefano Zampini 4880c2bc6bfSStefano Zampini /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh */ 4899566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "_glvis_mesh_coords", (PetscObject *)&hovec)); 4909566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)hovec)); 491066ea43fSLisandro Dalcin if (!hovec) { 492066ea43fSLisandro Dalcin DM cdm; 493066ea43fSLisandro Dalcin PetscFE disc; 494066ea43fSLisandro Dalcin PetscClassId classid; 495066ea43fSLisandro Dalcin 4969566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 4979566063dSJacob Faibussowitsch PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&disc)); 4989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId((PetscObject)disc, &classid)); 499066ea43fSLisandro Dalcin if (classid == PETSCFE_CLASSID) { 500066ea43fSLisandro Dalcin DM hocdm; 501066ea43fSLisandro Dalcin PetscFE hodisc; 502066ea43fSLisandro Dalcin Vec vec; 503066ea43fSLisandro Dalcin Mat mat; 504066ea43fSLisandro Dalcin char name[32], fec_type[64]; 5050c2bc6bfSStefano Zampini IS perm = NULL; 506066ea43fSLisandro Dalcin 5070c2bc6bfSStefano Zampini PetscCall(GLVisCreateFE(disc, name, &hodisc, &perm)); 5089566063dSJacob Faibussowitsch PetscCall(DMClone(cdm, &hocdm)); 5099566063dSJacob Faibussowitsch PetscCall(DMSetField(hocdm, 0, NULL, (PetscObject)hodisc)); 5109566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&hodisc)); 5119566063dSJacob Faibussowitsch PetscCall(DMCreateDS(hocdm)); 512066ea43fSLisandro Dalcin 5139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &vec)); 5149566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(hocdm, &hovec)); 5159566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdm, hocdm, &mat, NULL)); 5169566063dSJacob Faibussowitsch PetscCall(MatInterpolate(mat, vec, hovec)); 5179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 5180c2bc6bfSStefano Zampini PetscCall(DMGetLocalSection(hocdm, &hoSection)); 5190c2bc6bfSStefano Zampini PetscCall(PetscSectionSetClosurePermutation(hoSection, (PetscObject)hocdm, depth, perm)); 5200c2bc6bfSStefano Zampini PetscCall(ISDestroy(&perm)); 5219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&hocdm)); 5229566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec_type, sizeof(fec_type), "FiniteElementCollection: %s", name)); 5239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)hovec, fec_type)); 524066ea43fSLisandro Dalcin } 525066ea43fSLisandro Dalcin } 52621414b21SStefano Zampini 5279566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 5282827ebadSStefano Zampini PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &p, NULL)); 529c3c203b2SStefano Zampini if (p >= 0) cEnd = p; 5309566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 5319566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 5329566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 5339566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &sdim)); 5349566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 53508401ef6SPierre Jolivet PetscCheck(coordinates || hovec, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Missing local coordinates vector"); 5368135c375SStefano Zampini 5378135c375SStefano Zampini /* 5388135c375SStefano Zampini a couple of sections of the mesh specification are disabled 5393e96840aSStefano 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) 54077eacf09SStefano Zampini - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 5413e6c54aaSStefano Zampini and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) 5428135c375SStefano Zampini */ 5433e96840aSStefano Zampini enable_boundary = PETSC_FALSE; 5448135c375SStefano Zampini enable_ncmesh = PETSC_FALSE; 5453e6c54aaSStefano Zampini enable_mfem = PETSC_FALSE; 546f86f7544SStefano Zampini enable_emark = PETSC_FALSE; 547f86f7544SStefano Zampini enable_bmark = PETSC_FALSE; 5487bf4dd16SStefano Zampini /* I'm tired of problems with negative values in the markers, disable them */ 549d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->prefix, "GLVis PetscViewer DMPlex Options", "PetscViewer"); 5509566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary", "Enable boundary section in mesh representation", NULL, enable_boundary, &enable_boundary, NULL)); 5519566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh", "Enable vertex_parents section in mesh representation (allows derefinement)", NULL, enable_ncmesh, &enable_ncmesh, NULL)); 5529566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem", "Dump a mesh that can be used with MFEM's FiniteElementSpaces", NULL, enable_mfem, &enable_mfem, NULL)); 5539566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_overlap", "Include overlap region in local meshes", NULL, view_ovl, &view_ovl, NULL)); 5549566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_emarker", "String for the material id label", NULL, emark, emark, sizeof(emark), &enable_emark)); 5559566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_bmarker", "String for the boundary id label", NULL, bmark, bmark, sizeof(bmark), &enable_bmark)); 556d0609cedSBarry Smith PetscOptionsEnd(); 557f86f7544SStefano Zampini if (enable_bmark) enable_boundary = PETSC_TRUE; 558f86f7544SStefano Zampini 5599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 5601dca8a05SBarry Smith PetscCheck(!enable_ncmesh || size == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not supported in parallel"); 5619371c9d4SSatish Balay PetscCheck(!enable_boundary || depth < 0 || dim == depth, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, 5629371c9d4SSatish Balay "Mesh must be interpolated. " 5637e1aca4eSStefano Zampini "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); 5649371c9d4SSatish Balay PetscCheck(!enable_ncmesh || depth < 0 || dim == depth, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, 5659371c9d4SSatish Balay "Mesh must be interpolated. " 5667e1aca4eSStefano Zampini "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0"); 567044a5661SStefano Zampini if (depth >= 0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 56863a3b9bcSJacob Faibussowitsch PetscCheck(depth == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported depth %" PetscInt_FMT ". You should interpolate the mesh first", depth); 569044a5661SStefano Zampini cellvertex = PETSC_TRUE; 570044a5661SStefano Zampini } 5718135c375SStefano Zampini 5728135c375SStefano Zampini /* Identify possible cells in the overlap */ 5738135c375SStefano Zampini novl = 0; 5748135c375SStefano Zampini pown = NULL; 5753924b612SStefano Zampini if (size > 1) { 5763e6c54aaSStefano Zampini IS globalNum = NULL; 5773e6c54aaSStefano Zampini const PetscInt *gNum; 5783e6c54aaSStefano Zampini PetscBool ovl = PETSC_FALSE; 5793e6c54aaSStefano Zampini 5809566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "_glvis_plex_gnum", (PetscObject *)&globalNum)); 581b135d7daSStefano Zampini if (!globalNum) { 58228d58a37SPierre Jolivet if (view_ovl) { 5839566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), cEnd - cStart, 0, 1, &globalNum)); 58428d58a37SPierre Jolivet } else { 585*e2b8d0fcSMatthew G. Knepley PetscCall(DMPlexCreateCellNumbering(dm, PETSC_TRUE, &globalNum)); 58628d58a37SPierre Jolivet } 5879566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "_glvis_plex_gnum", (PetscObject)globalNum)); 5889566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)globalNum)); 589b135d7daSStefano Zampini } 5909566063dSJacob Faibussowitsch PetscCall(ISGetIndices(globalNum, &gNum)); 5918135c375SStefano Zampini for (p = cStart; p < cEnd; p++) { 5928135c375SStefano Zampini if (gNum[p - cStart] < 0) { 5938135c375SStefano Zampini ovl = PETSC_TRUE; 5948135c375SStefano Zampini novl++; 5958135c375SStefano Zampini } 5968135c375SStefano Zampini } 5978135c375SStefano Zampini if (ovl) { 5988135c375SStefano Zampini /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 5998135c375SStefano Zampini TODO: garbage collector? attach pown to dm? */ 6009566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(cEnd - cStart, &pown)); 6013e6c54aaSStefano Zampini for (p = cStart; p < cEnd; p++) { 6023e6c54aaSStefano Zampini if (gNum[p - cStart] < 0) continue; 60348a46eb9SPierre Jolivet else PetscCall(PetscBTSet(pown, p - cStart)); 6048135c375SStefano Zampini } 6053e6c54aaSStefano Zampini } 6069566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(globalNum, &gNum)); 6073e6c54aaSStefano Zampini } 6088135c375SStefano Zampini 6091d4815f0SStefano Zampini /* vertex_parents (Non-conforming meshes) */ 6101d4815f0SStefano Zampini parentSection = NULL; 6111d4815f0SStefano Zampini if (enable_ncmesh) { 6129566063dSJacob Faibussowitsch PetscCall(DMPlexGetTree(dm, &parentSection, NULL, NULL, NULL, NULL)); 6131d4815f0SStefano Zampini enable_ncmesh = (PetscBool)(enable_ncmesh && parentSection); 6141d4815f0SStefano Zampini } 6153e6c54aaSStefano Zampini /* return if this process is disabled */ 6168135c375SStefano Zampini if (!enabled) { 6179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh %s\n", enable_ncmesh ? "v1.1" : "v1.0")); 6189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n")); 61963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim)); 6209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n")); 62163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n")); 6229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n")); 62363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n")); 6249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n")); 62563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n")); 62663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim)); 6279566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&pown)); 6289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&hovec)); 6293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6308135c375SStefano Zampini } 6318135c375SStefano Zampini 6323e6c54aaSStefano Zampini if (enable_mfem) { 6336858538eSMatthew G. Knepley if (localized && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */ 6343e6c54aaSStefano Zampini PetscInt vpc = 0; 6353e6c54aaSStefano Zampini char fec[64]; 63696ca5757SLisandro Dalcin PetscInt vids[8] = {0, 1, 2, 3, 4, 5, 6, 7}; 63796ca5757SLisandro Dalcin PetscInt hexv[8] = {0, 1, 3, 2, 4, 5, 7, 6}, tetv[4] = {0, 1, 2, 3}; 63896ca5757SLisandro Dalcin PetscInt quadv[8] = {0, 1, 3, 2}, triv[3] = {0, 1, 2}; 63996ca5757SLisandro Dalcin PetscInt *dof = NULL; 6403e6c54aaSStefano Zampini PetscScalar *array, *ptr; 6413e6c54aaSStefano Zampini 64263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec, sizeof(fec), "FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P1", dim)); 6433e6c54aaSStefano Zampini if (cEnd - cStart) { 6443e6c54aaSStefano Zampini PetscInt fpc; 6453e6c54aaSStefano Zampini 6469566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cStart, &fpc)); 6473e6c54aaSStefano Zampini switch (dim) { 6483e6c54aaSStefano Zampini case 1: 6493e6c54aaSStefano Zampini vpc = 2; 6503e6c54aaSStefano Zampini dof = hexv; 6513e6c54aaSStefano Zampini break; 6523e6c54aaSStefano Zampini case 2: 6533e6c54aaSStefano Zampini switch (fpc) { 6543e6c54aaSStefano Zampini case 3: 6553e6c54aaSStefano Zampini vpc = 3; 656044a5661SStefano Zampini dof = triv; 6573e6c54aaSStefano Zampini break; 6583e6c54aaSStefano Zampini case 4: 6593e6c54aaSStefano Zampini vpc = 4; 660044a5661SStefano Zampini dof = quadv; 6613e6c54aaSStefano Zampini break; 662d71ae5a4SJacob Faibussowitsch default: 663d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: faces per cell %" PetscInt_FMT, fpc); 6643e6c54aaSStefano Zampini } 6653e6c54aaSStefano Zampini break; 6663e6c54aaSStefano Zampini case 3: 6673e6c54aaSStefano Zampini switch (fpc) { 668044a5661SStefano Zampini case 4: /* TODO: still need to understand L2 ordering for tets */ 6693e6c54aaSStefano Zampini vpc = 4; 670044a5661SStefano Zampini dof = tetv; 671044a5661SStefano Zampini SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled tethraedral case"); 6723e6c54aaSStefano Zampini case 6: 67363a3b9bcSJacob Faibussowitsch PetscCheck(!cellvertex, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: vertices per cell %" PetscInt_FMT, fpc); 674044a5661SStefano Zampini vpc = 8; 675044a5661SStefano Zampini dof = hexv; 676044a5661SStefano Zampini break; 677044a5661SStefano Zampini case 8: 67863a3b9bcSJacob Faibussowitsch PetscCheck(cellvertex, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: faces per cell %" PetscInt_FMT, fpc); 6793e6c54aaSStefano Zampini vpc = 8; 6803e6c54aaSStefano Zampini dof = hexv; 6813e6c54aaSStefano Zampini break; 682d71ae5a4SJacob Faibussowitsch default: 683d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unhandled case: faces per cell %" PetscInt_FMT, fpc); 6843e6c54aaSStefano Zampini } 6853e6c54aaSStefano Zampini break; 686d71ae5a4SJacob Faibussowitsch default: 687d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unhandled dim"); 6883e6c54aaSStefano Zampini } 6899566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(dm, cStart, vids)); 6903e6c54aaSStefano Zampini } 69128b400f6SJacob Faibussowitsch PetscCheck(dof, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing dofs"); 6929566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart - novl) * vpc * sdim, &hovec)); 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)hovec, fec)); 6949566063dSJacob Faibussowitsch PetscCall(VecGetArray(hovec, &array)); 6953e6c54aaSStefano Zampini ptr = array; 6963e6c54aaSStefano Zampini for (p = cStart; p < cEnd; p++) { 6973e6c54aaSStefano Zampini PetscInt csize, v, d; 6983e6c54aaSStefano Zampini PetscScalar *vals = NULL; 6993e6c54aaSStefano Zampini 7003e6c54aaSStefano Zampini if (PetscUnlikely(pown && !PetscBTLookup(pown, p - cStart))) continue; 7019566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, p, &csize, &vals)); 7021dca8a05SBarry Smith PetscCheck(csize == vpc * sdim || csize == vpc * sdim * 2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported closure size %" PetscInt_FMT " (vpc %" PetscInt_FMT ", sdim %" PetscInt_FMT ")", csize, vpc, sdim); 7033e6c54aaSStefano Zampini for (v = 0; v < vpc; v++) { 704ad540459SPierre Jolivet for (d = 0; d < sdim; d++) ptr[sdim * dof[v] + d] = vals[sdim * vids[v] + d]; 7053e6c54aaSStefano Zampini } 7063e6c54aaSStefano Zampini ptr += vpc * sdim; 7079566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, p, &csize, &vals)); 7083e6c54aaSStefano Zampini } 7099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(hovec, &array)); 7103e6c54aaSStefano Zampini } 7113e6c54aaSStefano Zampini } 7123e96840aSStefano Zampini /* if we have high-order coordinates in 3D, we need to specify the boundary */ 7133e96840aSStefano Zampini if (hovec && dim == 3) enable_boundary = PETSC_TRUE; 7143e6c54aaSStefano Zampini 7158135c375SStefano Zampini /* header */ 7169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh %s\n", enable_ncmesh ? "v1.1" : "v1.0")); 7178135c375SStefano Zampini 7188135c375SStefano Zampini /* topological dimension */ 7199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n")); 72063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim)); 7218135c375SStefano Zampini 7228135c375SStefano Zampini /* elements */ 723f86f7544SStefano Zampini minl = 1; 724f86f7544SStefano Zampini label = NULL; 725f86f7544SStefano Zampini if (enable_emark) { 726f86f7544SStefano Zampini PetscInt lminl = PETSC_MAX_INT; 727f86f7544SStefano Zampini 7289566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, emark, &label)); 729f86f7544SStefano Zampini if (label) { 730f86f7544SStefano Zampini IS vals; 731f86f7544SStefano Zampini PetscInt ldef; 732f86f7544SStefano Zampini 7339566063dSJacob Faibussowitsch PetscCall(DMLabelGetDefaultValue(label, &ldef)); 7349566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vals)); 7359566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(vals, &lminl, NULL)); 7369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vals)); 737f86f7544SStefano Zampini lminl = PetscMin(ldef, lminl); 738f86f7544SStefano Zampini } 7391c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&lminl, &minl, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm))); 740f86f7544SStefano Zampini if (minl == PETSC_MAX_INT) minl = 1; 741f86f7544SStefano Zampini } 7429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n")); 74363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", cEnd - cStart - novl)); 7448135c375SStefano Zampini for (p = cStart; p < cEnd; p++) { 74596ca5757SLisandro Dalcin PetscInt vids[8]; 74611a4995dSStefano Zampini PetscInt i, nv = 0, cid = -1, mid = 1; 7478135c375SStefano Zampini 7483e6c54aaSStefano Zampini if (PetscUnlikely(pown && !PetscBTLookup(pown, p - cStart))) continue; 7499566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointMFEMCellID_Internal(dm, label, minl, p, &mid, &cid)); 7509566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm, p, (localized && !hovec) ? coordSection : NULL, &nv, vids)); 7519566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(dm, p, vids)); 75263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT, mid, cid)); 75348a46eb9SPierre Jolivet for (i = 0; i < nv; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, vids[i])); 7549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 7558135c375SStefano Zampini } 7568135c375SStefano Zampini 757cc0d3ed7SStefano Zampini /* boundary */ 7589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n")); 759cc0d3ed7SStefano Zampini if (!enable_boundary) { 76063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n")); 761cc0d3ed7SStefano Zampini } else { 76277eacf09SStefano Zampini DMLabel perLabel; 76377eacf09SStefano Zampini PetscBT bfaces; 764b135d7daSStefano Zampini PetscInt fStart, fEnd, *fcells; 765cc0d3ed7SStefano Zampini 7669566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 7679566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(fEnd - fStart, &bfaces)); 7689566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxSizes(dm, NULL, &p)); 7699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p, &fcells)); 7709566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "glvis_periodic_cut", &perLabel)); 7716858538eSMatthew G. Knepley if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 7729566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "glvis_periodic_cut")); 7739566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "glvis_periodic_cut", &perLabel)); 7749566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(perLabel, 1)); 7756858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateSection(dm, &coordSectionCell)); 7766858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &coordinatesCell)); 77777eacf09SStefano Zampini for (p = cStart; p < cEnd; p++) { 778c3c203b2SStefano Zampini DMPolytopeType cellType; 779c3c203b2SStefano Zampini PetscInt dof; 780b135d7daSStefano Zampini 7819566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, p, &cellType)); 7826858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(coordSectionCell, p, &dof)); 78377eacf09SStefano Zampini if (dof) { 7846858538eSMatthew G. Knepley PetscInt uvpc, v, csize, csizeCell, cellClosureSize, *cellClosure = NULL, *vidxs = NULL; 7856858538eSMatthew G. Knepley PetscScalar *vals = NULL, *valsCell = NULL; 786c3c203b2SStefano Zampini 787c3c203b2SStefano Zampini uvpc = DMPolytopeTypeGetNumVertices(cellType); 78863a3b9bcSJacob Faibussowitsch PetscCheck(dof % sdim == 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Incompatible number of cell dofs %" PetscInt_FMT " and space dimension %" PetscInt_FMT, dof, sdim); 7899566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, p, &csize, &vals)); 7906858538eSMatthew G. Knepley PetscCall(DMPlexVecGetClosure(dm, coordSectionCell, coordinatesCell, p, &csizeCell, &valsCell)); 7916858538eSMatthew G. Knepley PetscCheck(csize == csizeCell, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cell %" PetscInt_FMT " has invalid localized coordinates", p); 7929566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &cellClosureSize, &cellClosure)); 79377eacf09SStefano Zampini for (v = 0; v < cellClosureSize; v++) 79477eacf09SStefano Zampini if (cellClosure[2 * v] >= vStart && cellClosure[2 * v] < vEnd) { 79577eacf09SStefano Zampini vidxs = cellClosure + 2 * v; 79677eacf09SStefano Zampini break; 79777eacf09SStefano Zampini } 79828b400f6SJacob Faibussowitsch PetscCheck(vidxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing vertices"); 799b135d7daSStefano Zampini for (v = 0; v < uvpc; v++) { 80077eacf09SStefano Zampini PetscInt s; 801044a5661SStefano Zampini 80277eacf09SStefano Zampini for (s = 0; s < sdim; s++) { 80348a46eb9SPierre Jolivet if (PetscAbsScalar(vals[v * sdim + s] - valsCell[v * sdim + s]) > PETSC_MACHINE_EPSILON) PetscCall(DMLabelSetValue(perLabel, vidxs[2 * v], 2)); 80477eacf09SStefano Zampini } 80577eacf09SStefano Zampini } 8069566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &cellClosureSize, &cellClosure)); 8079566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, p, &csize, &vals)); 8086858538eSMatthew G. Knepley PetscCall(DMPlexVecRestoreClosure(dm, coordSectionCell, coordinatesCell, p, &csizeCell, &valsCell)); 80977eacf09SStefano Zampini } 81077eacf09SStefano Zampini } 81177eacf09SStefano Zampini if (dim > 1) { 812b135d7daSStefano Zampini PetscInt eEnd, eStart; 813044a5661SStefano Zampini 8149566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 81577eacf09SStefano Zampini for (p = eStart; p < eEnd; p++) { 81677eacf09SStefano Zampini const PetscInt *cone; 81777eacf09SStefano Zampini PetscInt coneSize, i; 81877eacf09SStefano Zampini PetscBool ispe = PETSC_TRUE; 81977eacf09SStefano Zampini 8209566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 8219566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 82277eacf09SStefano Zampini for (i = 0; i < coneSize; i++) { 82377eacf09SStefano Zampini PetscInt v; 82477eacf09SStefano Zampini 8259566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(perLabel, cone[i], &v)); 82677eacf09SStefano Zampini ispe = (PetscBool)(ispe && (v == 2)); 82777eacf09SStefano Zampini } 82877eacf09SStefano Zampini if (ispe && coneSize) { 8293e96840aSStefano Zampini PetscInt ch, numChildren; 8303e96840aSStefano Zampini const PetscInt *children; 8313e96840aSStefano Zampini 8329566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(perLabel, p, 2)); 8339566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeChildren(dm, p, &numChildren, &children)); 83448a46eb9SPierre Jolivet for (ch = 0; ch < numChildren; ch++) PetscCall(DMLabelSetValue(perLabel, children[ch], 2)); 83577eacf09SStefano Zampini } 83677eacf09SStefano Zampini } 83777eacf09SStefano Zampini if (dim > 2) { 83877eacf09SStefano Zampini for (p = fStart; p < fEnd; p++) { 83977eacf09SStefano Zampini const PetscInt *cone; 84077eacf09SStefano Zampini PetscInt coneSize, i; 84177eacf09SStefano Zampini PetscBool ispe = PETSC_TRUE; 84277eacf09SStefano Zampini 8439566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 8449566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 84577eacf09SStefano Zampini for (i = 0; i < coneSize; i++) { 84677eacf09SStefano Zampini PetscInt v; 84777eacf09SStefano Zampini 8489566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(perLabel, cone[i], &v)); 84977eacf09SStefano Zampini ispe = (PetscBool)(ispe && (v == 2)); 85077eacf09SStefano Zampini } 85177eacf09SStefano Zampini if (ispe && coneSize) { 8523e96840aSStefano Zampini PetscInt ch, numChildren; 8533e96840aSStefano Zampini const PetscInt *children; 8543e96840aSStefano Zampini 8559566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(perLabel, p, 2)); 8569566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeChildren(dm, p, &numChildren, &children)); 85748a46eb9SPierre Jolivet for (ch = 0; ch < numChildren; ch++) PetscCall(DMLabelSetValue(perLabel, children[ch], 2)); 85877eacf09SStefano Zampini } 85977eacf09SStefano Zampini } 86077eacf09SStefano Zampini } 86177eacf09SStefano Zampini } 86277eacf09SStefano Zampini } 86377eacf09SStefano Zampini for (p = fStart; p < fEnd; p++) { 86477eacf09SStefano Zampini const PetscInt *support; 8658135c375SStefano Zampini PetscInt supportSize; 86677eacf09SStefano Zampini PetscBool isbf = PETSC_FALSE; 8678135c375SStefano Zampini 8689566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 8693e6c54aaSStefano Zampini if (pown) { 8708135c375SStefano Zampini PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 87177eacf09SStefano Zampini PetscInt i; 87277eacf09SStefano Zampini 8739566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, p, &support)); 87477eacf09SStefano Zampini for (i = 0; i < supportSize; i++) { 87577eacf09SStefano Zampini if (PetscLikely(PetscBTLookup(pown, support[i] - cStart))) has_owned = PETSC_TRUE; 87677eacf09SStefano Zampini else has_ghost = PETSC_TRUE; 87777eacf09SStefano Zampini } 87877eacf09SStefano Zampini isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 87977eacf09SStefano Zampini } else { 88077eacf09SStefano Zampini isbf = (PetscBool)(supportSize == 1); 88177eacf09SStefano Zampini } 88277eacf09SStefano Zampini if (!isbf && perLabel) { 88377eacf09SStefano Zampini const PetscInt *cone; 88477eacf09SStefano Zampini PetscInt coneSize, i; 88577eacf09SStefano Zampini 8869566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 8879566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 88877eacf09SStefano Zampini isbf = PETSC_TRUE; 88977eacf09SStefano Zampini for (i = 0; i < coneSize; i++) { 89077eacf09SStefano Zampini PetscInt v, d; 89177eacf09SStefano Zampini 8929566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(perLabel, cone[i], &v)); 8939566063dSJacob Faibussowitsch PetscCall(DMLabelGetDefaultValue(perLabel, &d)); 89477eacf09SStefano Zampini isbf = (PetscBool)(isbf && v != d); 89577eacf09SStefano Zampini } 89677eacf09SStefano Zampini } 8971baa6e33SBarry Smith if (isbf) PetscCall(PetscBTSet(bfaces, p - fStart)); 89877eacf09SStefano Zampini } 89977eacf09SStefano Zampini /* count boundary faces */ 90077eacf09SStefano Zampini for (p = fStart, bf = 0; p < fEnd; p++) { 90177eacf09SStefano Zampini if (PetscUnlikely(PetscBTLookup(bfaces, p - fStart))) { 90277eacf09SStefano Zampini const PetscInt *support; 90377eacf09SStefano Zampini PetscInt supportSize, c; 9048135c375SStefano Zampini 9059566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 9069566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, p, &support)); 90777eacf09SStefano Zampini for (c = 0; c < supportSize; c++) { 9083e96840aSStefano Zampini const PetscInt *cone; 909b135d7daSStefano Zampini PetscInt cell, cl, coneSize; 9103e96840aSStefano Zampini 9113e96840aSStefano Zampini cell = support[c]; 9123e96840aSStefano Zampini if (pown && PetscUnlikely(!PetscBTLookup(pown, cell - cStart))) continue; 9139566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 9149566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 915b135d7daSStefano Zampini for (cl = 0; cl < coneSize; cl++) { 9163e96840aSStefano Zampini if (cone[cl] == p) { 9173e96840aSStefano Zampini bf += 1; 9183e96840aSStefano Zampini break; 9198135c375SStefano Zampini } 92077eacf09SStefano Zampini } 9213e96840aSStefano Zampini } 9228135c375SStefano Zampini } 9238135c375SStefano Zampini } 924f86f7544SStefano Zampini minl = 1; 925f86f7544SStefano Zampini label = NULL; 926f86f7544SStefano Zampini if (enable_bmark) { 927f86f7544SStefano Zampini PetscInt lminl = PETSC_MAX_INT; 928f86f7544SStefano Zampini 9299566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, bmark, &label)); 930f86f7544SStefano Zampini if (label) { 931f86f7544SStefano Zampini IS vals; 932f86f7544SStefano Zampini PetscInt ldef; 933f86f7544SStefano Zampini 9349566063dSJacob Faibussowitsch PetscCall(DMLabelGetDefaultValue(label, &ldef)); 9359566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vals)); 9369566063dSJacob Faibussowitsch PetscCall(ISGetMinMax(vals, &lminl, NULL)); 9379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vals)); 938f86f7544SStefano Zampini lminl = PetscMin(ldef, lminl); 939f86f7544SStefano Zampini } 9401c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&lminl, &minl, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm))); 941f86f7544SStefano Zampini if (minl == PETSC_MAX_INT) minl = 1; 942f86f7544SStefano Zampini } 94363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", bf)); 9448135c375SStefano Zampini for (p = fStart; p < fEnd; p++) { 94577eacf09SStefano Zampini if (PetscUnlikely(PetscBTLookup(bfaces, p - fStart))) { 9468135c375SStefano Zampini const PetscInt *support; 94777eacf09SStefano Zampini PetscInt supportSize, c, nc = 0; 9488135c375SStefano Zampini 9499566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &supportSize)); 9509566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, p, &support)); 9513e6c54aaSStefano Zampini if (pown) { 95277eacf09SStefano Zampini for (c = 0; c < supportSize; c++) { 953ad540459SPierre Jolivet if (PetscLikely(PetscBTLookup(pown, support[c] - cStart))) fcells[nc++] = support[c]; 9548135c375SStefano Zampini } 9559371c9d4SSatish Balay } else 9569371c9d4SSatish Balay for (c = 0; c < supportSize; c++) fcells[nc++] = support[c]; 95777eacf09SStefano Zampini for (c = 0; c < nc; c++) { 958c3c203b2SStefano Zampini const DMPolytopeType *faceTypes; 959c3c203b2SStefano Zampini DMPolytopeType cellType; 960c3c203b2SStefano Zampini const PetscInt *faceSizes, *cone; 961c3c203b2SStefano Zampini PetscInt vids[8], *faces, st, i, coneSize, cell, cl, nv, cid = -1, mid = -1; 9628135c375SStefano Zampini 96377eacf09SStefano Zampini cell = fcells[c]; 9649566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 9659566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 966b135d7daSStefano Zampini for (cl = 0; cl < coneSize; cl++) 9679371c9d4SSatish Balay if (cone[cl] == p) break; 968b135d7daSStefano Zampini if (cl == coneSize) continue; 9698135c375SStefano Zampini 97077eacf09SStefano Zampini /* face material id and type */ 9719566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointMFEMCellID_Internal(dm, label, minl, p, &mid, &cid)); 97263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT, mid, cid)); 97377eacf09SStefano Zampini /* vertex ids */ 9749566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cell, &cellType)); 9759566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm, cell, (localized && !hovec) ? coordSection : NULL, &nv, vids)); 9769566063dSJacob Faibussowitsch PetscCall(DMPlexGetRawFaces_Internal(dm, cellType, vids, NULL, &faceTypes, &faceSizes, (const PetscInt **)&faces)); 977c3c203b2SStefano Zampini st = 0; 978c3c203b2SStefano Zampini for (i = 0; i < cl; i++) st += faceSizes[i]; 9799566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCell(faceTypes[cl], faces + st)); 98048a46eb9SPierre Jolivet for (i = 0; i < faceSizes[cl]; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, faces[st + i])); 9819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 9829566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreRawFaces_Internal(dm, cellType, vids, NULL, &faceTypes, &faceSizes, (const PetscInt **)&faces)); 9833e96840aSStefano Zampini bf -= 1; 98477eacf09SStefano Zampini } 9858135c375SStefano Zampini } 9868135c375SStefano Zampini } 98763a3b9bcSJacob Faibussowitsch PetscCheck(!bf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Remaining boundary faces %" PetscInt_FMT, bf); 9889566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&bfaces)); 9899566063dSJacob Faibussowitsch PetscCall(PetscFree(fcells)); 9908135c375SStefano Zampini } 9918135c375SStefano Zampini 9928135c375SStefano Zampini /* mark owned vertices */ 9933e6c54aaSStefano Zampini vown = NULL; 9943e6c54aaSStefano Zampini if (pown) { 9959566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(vEnd - vStart, &vown)); 9968135c375SStefano Zampini for (p = cStart; p < cEnd; p++) { 9978135c375SStefano Zampini PetscInt i, closureSize, *closure = NULL; 9988135c375SStefano Zampini 9993e6c54aaSStefano Zampini if (PetscUnlikely(!PetscBTLookup(pown, p - cStart))) continue; 10009566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure)); 10018135c375SStefano Zampini for (i = 0; i < closureSize; i++) { 10028135c375SStefano Zampini const PetscInt pp = closure[2 * i]; 10038135c375SStefano Zampini 100448a46eb9SPierre Jolivet if (pp >= vStart && pp < vEnd) PetscCall(PetscBTSet(vown, pp - vStart)); 10058135c375SStefano Zampini } 10069566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure)); 10078135c375SStefano Zampini } 10088135c375SStefano Zampini } 10098135c375SStefano Zampini 10108135c375SStefano Zampini if (parentSection) { 10118135c375SStefano Zampini PetscInt vp, gvp; 10128135c375SStefano Zampini 10138135c375SStefano Zampini for (vp = 0, p = vStart; p < vEnd; p++) { 10148135c375SStefano Zampini DMLabel dlabel; 10158135c375SStefano Zampini PetscInt parent, depth; 10168135c375SStefano Zampini 10173e6c54aaSStefano Zampini if (PetscUnlikely(vown && !PetscBTLookup(vown, p - vStart))) continue; 10189566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &dlabel)); 10199566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(dlabel, p, &depth)); 10209566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeParent(dm, p, &parent, NULL)); 10218135c375SStefano Zampini if (parent != p) vp++; 10228135c375SStefano Zampini } 10231c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&vp, &gvp, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 10248135c375SStefano Zampini if (gvp) { 10258135c375SStefano Zampini PetscInt maxsupp; 10268135c375SStefano Zampini PetscBool *skip = NULL; 10278135c375SStefano Zampini 10289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertex_parents\n")); 102963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", vp)); 10309566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxSizes(dm, NULL, &maxsupp)); 10319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxsupp, &skip)); 10328135c375SStefano Zampini for (p = vStart; p < vEnd; p++) { 10338135c375SStefano Zampini DMLabel dlabel; 10348135c375SStefano Zampini PetscInt parent; 10358135c375SStefano Zampini 10363e6c54aaSStefano Zampini if (PetscUnlikely(vown && !PetscBTLookup(vown, p - vStart))) continue; 10379566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &dlabel)); 10389566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeParent(dm, p, &parent, NULL)); 10398135c375SStefano Zampini if (parent != p) { 104096ca5757SLisandro Dalcin PetscInt vids[8] = {-1, -1, -1, -1, -1, -1, -1, -1}; /* silent overzealous clang static analyzer */ 10413924b612SStefano Zampini PetscInt i, nv, ssize, n, numChildren, depth = -1; 10428135c375SStefano Zampini const PetscInt *children; 10433924b612SStefano Zampini 10449566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, parent, &ssize)); 10453924b612SStefano Zampini switch (ssize) { 10468135c375SStefano Zampini case 2: /* edge */ 10478135c375SStefano Zampini nv = 0; 10489566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointMFEMVertexIDs_Internal(dm, parent, localized ? coordSection : NULL, &nv, vids)); 104963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, p - vStart)); 105048a46eb9SPierre Jolivet for (i = 0; i < nv; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, vids[i])); 10519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 10528135c375SStefano Zampini vp--; 10538135c375SStefano Zampini break; 10548135c375SStefano Zampini case 4: /* face */ 10559566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeChildren(dm, parent, &numChildren, &children)); 10568135c375SStefano Zampini for (n = 0; n < numChildren; n++) { 10579566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(dlabel, children[n], &depth)); 10588135c375SStefano Zampini if (!depth) { 10598135c375SStefano Zampini const PetscInt *hvsupp, *hesupp, *cone; 10608135c375SStefano Zampini PetscInt hvsuppSize, hesuppSize, coneSize; 1061451a39c7SStefano Zampini PetscInt hv = children[n], he = -1, f; 10628135c375SStefano Zampini 10639566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(skip, maxsupp)); 10649566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, hv, &hvsuppSize)); 10659566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, hv, &hvsupp)); 10668135c375SStefano Zampini for (i = 0; i < hvsuppSize; i++) { 10678135c375SStefano Zampini PetscInt ep; 10689566063dSJacob Faibussowitsch PetscCall(DMPlexGetTreeParent(dm, hvsupp[i], &ep, NULL)); 10698135c375SStefano Zampini if (ep != hvsupp[i]) { 10708135c375SStefano Zampini he = hvsupp[i]; 10718135c375SStefano Zampini } else { 10728135c375SStefano Zampini skip[i] = PETSC_TRUE; 10738135c375SStefano Zampini } 10748135c375SStefano Zampini } 107563a3b9bcSJacob Faibussowitsch PetscCheck(he != -1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Vertex %" PetscInt_FMT " support size %" PetscInt_FMT ": hanging edge not found", hv, hvsuppSize); 10769566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, he, &cone)); 107796ca5757SLisandro Dalcin vids[0] = (cone[0] == hv) ? cone[1] : cone[0]; 10789566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, he, &hesuppSize)); 10799566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, he, &hesupp)); 10808135c375SStefano Zampini for (f = 0; f < hesuppSize; f++) { 10818135c375SStefano Zampini PetscInt j; 10828135c375SStefano Zampini 10839566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, hesupp[f], &cone)); 10849566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, hesupp[f], &coneSize)); 10858135c375SStefano Zampini for (j = 0; j < coneSize; j++) { 10868135c375SStefano Zampini PetscInt k; 10878135c375SStefano Zampini for (k = 0; k < hvsuppSize; k++) { 10888135c375SStefano Zampini if (hvsupp[k] == cone[j]) { 10898135c375SStefano Zampini skip[k] = PETSC_TRUE; 10908135c375SStefano Zampini break; 10918135c375SStefano Zampini } 10928135c375SStefano Zampini } 10938135c375SStefano Zampini } 10948135c375SStefano Zampini } 10958135c375SStefano Zampini for (i = 0; i < hvsuppSize; i++) { 10968135c375SStefano Zampini if (!skip[i]) { 10979566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, hvsupp[i], &cone)); 109896ca5757SLisandro Dalcin vids[1] = (cone[0] == hv) ? cone[1] : cone[0]; 10998135c375SStefano Zampini } 11008135c375SStefano Zampini } 110163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT, hv - vStart)); 110248a46eb9SPierre Jolivet for (i = 0; i < 2; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, vids[i] - vStart)); 11039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 11048135c375SStefano Zampini vp--; 11058135c375SStefano Zampini } 11068135c375SStefano Zampini } 11078135c375SStefano Zampini break; 1108d71ae5a4SJacob Faibussowitsch default: 1109d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Don't know how to deal with support size %" PetscInt_FMT, ssize); 11108135c375SStefano Zampini } 11118135c375SStefano Zampini } 11128135c375SStefano Zampini } 11139566063dSJacob Faibussowitsch PetscCall(PetscFree(skip)); 11148135c375SStefano Zampini } 111563a3b9bcSJacob Faibussowitsch PetscCheck(!vp, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected %" PetscInt_FMT " hanging vertices", vp); 11168135c375SStefano Zampini } 11179566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&vown)); 11188135c375SStefano Zampini 11198135c375SStefano Zampini /* vertices */ 112077eacf09SStefano Zampini if (hovec) { /* higher-order meshes */ 112177eacf09SStefano Zampini const char *fec; 11220286d493SLisandro Dalcin PetscInt i, n, s; 11239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n")); 112463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", vEnd - vStart)); 11259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nodes\n")); 11269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)hovec, &fec)); 11279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n")); 11289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", fec)); 112963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", sdim)); 11309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: 1\n\n")); /*Ordering::byVDIM*/ 11310c2bc6bfSStefano Zampini if (hoSection) { 11320c2bc6bfSStefano Zampini DM cdm; 11330c2bc6bfSStefano Zampini 11340c2bc6bfSStefano Zampini PetscCall(VecGetDM(hovec, &cdm)); 11350c2bc6bfSStefano Zampini for (p = cStart; p < cEnd; p++) { 11360c2bc6bfSStefano Zampini PetscScalar *vals = NULL; 11370c2bc6bfSStefano Zampini PetscInt csize; 11380c2bc6bfSStefano Zampini 11390c2bc6bfSStefano Zampini if (PetscUnlikely(pown && !PetscBTLookup(pown, p - cStart))) continue; 11400c2bc6bfSStefano Zampini PetscCall(DMPlexVecGetClosure(cdm, hoSection, hovec, p, &csize, &vals)); 114163a3b9bcSJacob Faibussowitsch PetscCheck(csize % sdim == 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Size of closure %" PetscInt_FMT " incompatible with space dimension %" PetscInt_FMT, csize, sdim); 11420c2bc6bfSStefano Zampini for (i = 0; i < csize / sdim; i++) { 114348a46eb9SPierre Jolivet for (s = 0; s < sdim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, (double)PetscRealPart(vals[i * sdim + s]))); 11440c2bc6bfSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 11450c2bc6bfSStefano Zampini } 11460c2bc6bfSStefano Zampini PetscCall(DMPlexVecRestoreClosure(cdm, hoSection, hovec, p, &csize, &vals)); 11470c2bc6bfSStefano Zampini } 11480c2bc6bfSStefano Zampini } else { 11499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(hovec, &array)); 11509566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(hovec, &n)); 115163a3b9bcSJacob Faibussowitsch PetscCheck(n % sdim == 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Size of local coordinate vector %" PetscInt_FMT " incompatible with space dimension %" PetscInt_FMT, n, sdim); 115277eacf09SStefano Zampini for (i = 0; i < n / sdim; i++) { 115348a46eb9SPierre Jolivet for (s = 0; s < sdim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, (double)PetscRealPart(array[i * sdim + s]))); 11549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 115577eacf09SStefano Zampini } 11569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(hovec, &array)); 11570c2bc6bfSStefano Zampini } 115877eacf09SStefano Zampini } else { 11599566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &nvert)); 11609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n")); 116163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", nvert / sdim)); 116263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim)); 11639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &array)); 1164cc0d3ed7SStefano Zampini for (p = 0; p < nvert / sdim; p++) { 1165cc0d3ed7SStefano Zampini PetscInt s; 1166cc0d3ed7SStefano Zampini for (s = 0; s < sdim; s++) { 11673e96840aSStefano Zampini PetscReal v = PetscRealPart(array[p * sdim + s]); 11683e96840aSStefano Zampini 11699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, fmt, PetscIsInfOrNanReal(v) ? 0.0 : (double)v)); 11708135c375SStefano Zampini } 11719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 11728135c375SStefano Zampini } 11739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &array)); 117477eacf09SStefano Zampini } 11750c2bc6bfSStefano Zampini PetscCall(PetscBTDestroy(&pown)); 11769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&hovec)); 11773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11788135c375SStefano Zampini } 11798135c375SStefano Zampini 1180d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 1181d71ae5a4SJacob Faibussowitsch { 11828135c375SStefano Zampini PetscFunctionBegin; 11839566063dSJacob Faibussowitsch PetscCall(DMView_GLVis(dm, viewer, DMPlexView_GLVis_ASCII)); 11843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11858135c375SStefano Zampini } 1186