11d5c9e23SVaclav Hapla #include <petscdmplex.h> 21d5c9e23SVaclav Hapla #include <petscviewerhdf5.h> 31d5c9e23SVaclav Hapla #include <petscsf.h> 41d5c9e23SVaclav Hapla 5*e22c0602SVaclav Hapla static const char help[] = "Test DMLabel I/O with PETSc native HDF5 mesh format\n\n"; 61d5c9e23SVaclav Hapla static const char EX[] = "ex56.c"; 71d5c9e23SVaclav Hapla typedef struct { 81d5c9e23SVaclav Hapla MPI_Comm comm; 91d5c9e23SVaclav Hapla const char *meshname; /* Mesh name */ 105c820da0SVaclav Hapla PetscInt num_labels; /* Asserted number of labels in loaded mesh */ 111d5c9e23SVaclav Hapla PetscBool compare; /* Compare the meshes using DMPlexEqual() and DMCompareLabels() */ 121d5c9e23SVaclav Hapla PetscBool compare_labels; /* Compare labels in the meshes using DMCompareLabels() */ 131d5c9e23SVaclav Hapla PetscBool compare_boundary; /* Check label I/O via boundary vertex coordinates */ 141d5c9e23SVaclav Hapla PetscBool compare_pre_post; /* Compare labels loaded before distribution with those loaded after distribution */ 151d5c9e23SVaclav Hapla char outfile[PETSC_MAX_PATH_LEN]; /* Output file */ 161d5c9e23SVaclav Hapla PetscBool use_low_level_functions; /* Use low level functions for viewing and loading */ 171d5c9e23SVaclav Hapla //TODO This is meant as temporary option; can be removed once we have full parallel loading in place 181d5c9e23SVaclav Hapla PetscBool distribute_after_topo_load; /* Distribute topology right after DMPlexTopologyLoad(), if use_low_level_functions=true */ 195c820da0SVaclav Hapla PetscInt verbose; 201d5c9e23SVaclav Hapla } AppCtx; 211d5c9e23SVaclav Hapla 221d5c9e23SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 231d5c9e23SVaclav Hapla { 241d5c9e23SVaclav Hapla PetscFunctionBeginUser; 251d5c9e23SVaclav Hapla options->comm = comm; 265c820da0SVaclav Hapla options->num_labels = -1; 271d5c9e23SVaclav Hapla options->compare = PETSC_FALSE; 281d5c9e23SVaclav Hapla options->compare_labels = PETSC_FALSE; 291d5c9e23SVaclav Hapla options->compare_boundary = PETSC_FALSE; 301d5c9e23SVaclav Hapla options->compare_pre_post = PETSC_FALSE; 311d5c9e23SVaclav Hapla options->outfile[0] = '\0'; 321d5c9e23SVaclav Hapla options->use_low_level_functions = PETSC_FALSE; 331d5c9e23SVaclav Hapla options->distribute_after_topo_load = PETSC_FALSE; 345c820da0SVaclav Hapla options->verbose = 0; 351d5c9e23SVaclav Hapla 36d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 375c820da0SVaclav Hapla PetscCall(PetscOptionsInt("-num_labels", "Asserted number of labels in meshfile; don't count depth and celltype; -1 to deactivate", EX, options->num_labels, &options->num_labels, NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare_pre_post", "Compare labels loaded before distribution with those loaded after distribution", "ex55.c", options->compare_pre_post, &options->compare_pre_post, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL)); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_low_level_functions", "Use low level functions for viewing and loading", EX, options->use_low_level_functions, &options->use_low_level_functions, NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-distribute_after_topo_load", "Distribute topology right after DMPlexTopologyLoad(), if use_low_level_functions=true", EX, options->distribute_after_topo_load, &options->distribute_after_topo_load, NULL)); 455c820da0SVaclav Hapla PetscCall(PetscOptionsInt("-verbose", "Verbosity level", EX, options->verbose, &options->verbose, NULL)); 46d0609cedSBarry Smith PetscOptionsEnd(); 471d5c9e23SVaclav Hapla PetscFunctionReturn(0); 481d5c9e23SVaclav Hapla }; 491d5c9e23SVaclav Hapla 501d5c9e23SVaclav Hapla static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm) 511d5c9e23SVaclav Hapla { 521d5c9e23SVaclav Hapla DM dm; 531d5c9e23SVaclav Hapla 541d5c9e23SVaclav Hapla PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(DMCreate(options->comm, &dm)); 569566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &options->meshname)); 599566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 601d5c9e23SVaclav Hapla *newdm = dm; 611d5c9e23SVaclav Hapla PetscFunctionReturn(0); 621d5c9e23SVaclav Hapla } 631d5c9e23SVaclav Hapla 641d5c9e23SVaclav Hapla static PetscErrorCode SaveMesh(AppCtx *options, DM dm) 651d5c9e23SVaclav Hapla { 661d5c9e23SVaclav Hapla PetscViewer v; 671d5c9e23SVaclav Hapla 681d5c9e23SVaclav Hapla PetscFunctionBeginUser; 699566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), options->outfile, FILE_MODE_WRITE, &v)); 701d5c9e23SVaclav Hapla if (options->use_low_level_functions) { 719566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyView(dm, v)); 729566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesView(dm, v)); 739566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsView(dm, v)); 741d5c9e23SVaclav Hapla } else { 759566063dSJacob Faibussowitsch PetscCall(DMView(dm, v)); 761d5c9e23SVaclav Hapla } 779566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&v)); 781d5c9e23SVaclav Hapla PetscFunctionReturn(0); 791d5c9e23SVaclav Hapla } 801d5c9e23SVaclav Hapla 811d5c9e23SVaclav Hapla typedef enum {NONE=0, PRE_DIST=1, POST_DIST=2} AuxObjLoadMode; 821d5c9e23SVaclav Hapla 831d5c9e23SVaclav Hapla static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm) 841d5c9e23SVaclav Hapla { 851d5c9e23SVaclav Hapla DM dm; 861d5c9e23SVaclav Hapla PetscSF sfXC; 871d5c9e23SVaclav Hapla 881d5c9e23SVaclav Hapla PetscFunctionBeginUser; 899566063dSJacob Faibussowitsch PetscCall(DMCreate(options->comm, &dm)); 909566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 919566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname)); 929566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyLoad(dm, v, &sfXC)); 931d5c9e23SVaclav Hapla if (mode == PRE_DIST) { 949566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC)); 959566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsLoad(dm, v, sfXC)); 961d5c9e23SVaclav Hapla } 971d5c9e23SVaclav Hapla if (explicitDistribute) { 981d5c9e23SVaclav Hapla DM dmdist; 991d5c9e23SVaclav Hapla PetscSF sfXB = sfXC, sfBC; 1001d5c9e23SVaclav Hapla 1019566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, &sfBC, &dmdist)); 1021d5c9e23SVaclav Hapla if (dmdist) { 1031d5c9e23SVaclav Hapla const char *name; 1041d5c9e23SVaclav Hapla 1059566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) dm, &name)); 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dmdist, name)); 1079566063dSJacob Faibussowitsch PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC)); 1089566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfXB)); 1099566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfBC)); 1109566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1111d5c9e23SVaclav Hapla dm = dmdist; 1121d5c9e23SVaclav Hapla } 1131d5c9e23SVaclav Hapla } 1141d5c9e23SVaclav Hapla if (mode == POST_DIST) { 1159566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC)); 1169566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsLoad(dm, v, sfXC)); 1171d5c9e23SVaclav Hapla } 1189566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfXC)); 1191d5c9e23SVaclav Hapla *newdm = dm; 1201d5c9e23SVaclav Hapla PetscFunctionReturn(0); 1211d5c9e23SVaclav Hapla } 1221d5c9e23SVaclav Hapla 1231d5c9e23SVaclav Hapla static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew) 1241d5c9e23SVaclav Hapla { 1251d5c9e23SVaclav Hapla DM dm; 1261d5c9e23SVaclav Hapla PetscViewer v; 1271d5c9e23SVaclav Hapla 1281d5c9e23SVaclav Hapla PetscFunctionBeginUser; 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(options->comm, options->outfile, FILE_MODE_READ, &v)); 1301d5c9e23SVaclav Hapla if (options->use_low_level_functions) { 1311d5c9e23SVaclav Hapla if (options->compare_pre_post) { 1321d5c9e23SVaclav Hapla DM dm0; 1331d5c9e23SVaclav Hapla 1349566063dSJacob Faibussowitsch PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0)); 1359566063dSJacob Faibussowitsch PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm)); 1369566063dSJacob Faibussowitsch PetscCall(DMCompareLabels(dm0, dm, NULL, NULL)); 1379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm0)); 1381d5c9e23SVaclav Hapla } else { 1399566063dSJacob Faibussowitsch PetscCall(LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm)); 1401d5c9e23SVaclav Hapla } 1411d5c9e23SVaclav Hapla } else { 1429566063dSJacob Faibussowitsch PetscCall(DMCreate(options->comm, &dm)); 1439566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 1449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname)); 1459566063dSJacob Faibussowitsch PetscCall(DMLoad(dm, v)); 1461d5c9e23SVaclav Hapla } 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&v)); 1481d5c9e23SVaclav Hapla 1499566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm, "load_")); 1509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 1519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 1521d5c9e23SVaclav Hapla *dmnew = dm; 1531d5c9e23SVaclav Hapla PetscFunctionReturn(0); 1541d5c9e23SVaclav Hapla } 1551d5c9e23SVaclav Hapla 1561d5c9e23SVaclav Hapla static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1) 1571d5c9e23SVaclav Hapla { 1581d5c9e23SVaclav Hapla PetscBool flg; 1591d5c9e23SVaclav Hapla 1601d5c9e23SVaclav Hapla PetscFunctionBeginUser; 1611d5c9e23SVaclav Hapla if (options->compare) { 1629566063dSJacob Faibussowitsch PetscCall(DMPlexEqual(dm0, dm1, &flg)); 16328b400f6SJacob Faibussowitsch PetscCheck(flg,options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal"); 1649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(options->comm,"DMs equal\n")); 1651d5c9e23SVaclav Hapla } 1661d5c9e23SVaclav Hapla if (options->compare_labels) { 1679566063dSJacob Faibussowitsch PetscCall(DMCompareLabels(dm0, dm1, NULL, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(options->comm,"DMLabels equal\n")); 1691d5c9e23SVaclav Hapla } 1701d5c9e23SVaclav Hapla PetscFunctionReturn(0); 1711d5c9e23SVaclav Hapla } 1721d5c9e23SVaclav Hapla 1735c820da0SVaclav Hapla static PetscErrorCode MarkBoundary(DM dm, const char name[], PetscInt value, PetscBool verticesOnly, DMLabel *label) 1741d5c9e23SVaclav Hapla { 1751d5c9e23SVaclav Hapla DMLabel l; 1765c820da0SVaclav Hapla 1775c820da0SVaclav Hapla PetscFunctionBeginUser; 1785c820da0SVaclav Hapla PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)dm), name, &l)); 1795c820da0SVaclav Hapla PetscCall(DMAddLabel(dm, l)); 1805c820da0SVaclav Hapla PetscCall(DMPlexMarkBoundaryFaces(dm, value, l)); 1815c820da0SVaclav Hapla PetscCall(DMPlexLabelComplete(dm, l)); 1825c820da0SVaclav Hapla if (verticesOnly) { 1831d5c9e23SVaclav Hapla IS points; 1841d5c9e23SVaclav Hapla const PetscInt *idx; 1851d5c9e23SVaclav Hapla PetscInt i, n; 1861d5c9e23SVaclav Hapla 1879566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(l, value, &points)); 1889566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(points, &n)); 1899566063dSJacob Faibussowitsch PetscCall(ISGetIndices(points, &idx)); 1901d5c9e23SVaclav Hapla for (i=0; i<n; i++) { 1911d5c9e23SVaclav Hapla const PetscInt p = idx[i]; 1921d5c9e23SVaclav Hapla PetscInt d; 1931d5c9e23SVaclav Hapla 1949566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointDepth(dm, p, &d)); 1951d5c9e23SVaclav Hapla if (d != 0) { 1969566063dSJacob Faibussowitsch PetscCall(DMLabelClearValue(l, p, value)); 1971d5c9e23SVaclav Hapla } 1981d5c9e23SVaclav Hapla } 1999566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(points, &idx)); 2009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&points)); 2015c820da0SVaclav Hapla } 2025c820da0SVaclav Hapla if (label) *label = l; 2035c820da0SVaclav Hapla else PetscCall(DMLabelDestroy(&l)); 2041d5c9e23SVaclav Hapla PetscFunctionReturn(0); 2051d5c9e23SVaclav Hapla } 2061d5c9e23SVaclav Hapla 2071d5c9e23SVaclav Hapla static PetscErrorCode VertexCoordinatesToAll(DM dm, IS vertices, Vec *allCoords) 2081d5c9e23SVaclav Hapla { 2091d5c9e23SVaclav Hapla Vec coords, allCoords_; 2101d5c9e23SVaclav Hapla VecScatter sc; 2111d5c9e23SVaclav Hapla MPI_Comm comm; 2121d5c9e23SVaclav Hapla 2131d5c9e23SVaclav Hapla PetscFunctionBeginUser; 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(dm)); 2161d5c9e23SVaclav Hapla if (vertices) { 2179566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords)); 2181d5c9e23SVaclav Hapla } else { 2199566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &coords)); 2201d5c9e23SVaclav Hapla } 2211d5c9e23SVaclav Hapla { 2221d5c9e23SVaclav Hapla PetscInt n; 2231d5c9e23SVaclav Hapla Vec mpivec; 2241d5c9e23SVaclav Hapla 2259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coords, &n)); 2269566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, n, PETSC_DECIDE, &mpivec)); 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(coords, mpivec)); 2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 2291d5c9e23SVaclav Hapla coords = mpivec; 2301d5c9e23SVaclav Hapla } 2311d5c9e23SVaclav Hapla 2329566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToAll(coords, &sc, &allCoords_)); 2339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD)); 2349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD)); 2359566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sc)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 2371d5c9e23SVaclav Hapla *allCoords = allCoords_; 2381d5c9e23SVaclav Hapla PetscFunctionReturn(0); 2391d5c9e23SVaclav Hapla } 2401d5c9e23SVaclav Hapla 2415c820da0SVaclav Hapla static PetscErrorCode DMLabelGetCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec *allCoords) 2421d5c9e23SVaclav Hapla { 2431d5c9e23SVaclav Hapla IS vertices; 2441d5c9e23SVaclav Hapla 2451d5c9e23SVaclav Hapla PetscFunctionBeginUser; 2465c820da0SVaclav Hapla PetscCall(DMLabelGetStratumIS(label, value, &vertices)); 2479566063dSJacob Faibussowitsch PetscCall(VertexCoordinatesToAll(dm, vertices, allCoords)); 2489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vertices)); 2491d5c9e23SVaclav Hapla PetscFunctionReturn(0); 2501d5c9e23SVaclav Hapla } 2511d5c9e23SVaclav Hapla 2525c820da0SVaclav Hapla static PetscErrorCode DMLabelCompareWithCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec allCoords, PetscInt verbose) 2531d5c9e23SVaclav Hapla { 2545c820da0SVaclav Hapla const char *labelName; 2551d5c9e23SVaclav Hapla IS pointsIS; 2561d5c9e23SVaclav Hapla const PetscInt *points; 2571d5c9e23SVaclav Hapla PetscInt i, n; 2581d5c9e23SVaclav Hapla PetscBool fail = PETSC_FALSE; 2591d5c9e23SVaclav Hapla MPI_Comm comm; 2601d5c9e23SVaclav Hapla PetscMPIInt rank; 2611d5c9e23SVaclav Hapla 2621d5c9e23SVaclav Hapla PetscFunctionBeginUser; 2635c820da0SVaclav Hapla PetscCheck(label,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label not loaded"); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2655c820da0SVaclav Hapla PetscCall(PetscObjectGetName((PetscObject)label, &labelName)); 2669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2671d5c9e23SVaclav Hapla { 2681d5c9e23SVaclav Hapla PetscInt pStart, pEnd; 2691d5c9e23SVaclav Hapla 2709566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 2719566063dSJacob Faibussowitsch PetscCall(DMLabelCreateIndex(label, pStart, pEnd)); 2721d5c9e23SVaclav Hapla } 2739566063dSJacob Faibussowitsch PetscCall(DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS)); 2749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(pointsIS, &points)); 2759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointsIS, &n)); 2765c820da0SVaclav Hapla if (verbose > 1) PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_(comm))); 2771d5c9e23SVaclav Hapla for (i=0; i<n; i++) { 2781d5c9e23SVaclav Hapla const PetscInt p = points[i]; 2791d5c9e23SVaclav Hapla PetscBool has; 2801d5c9e23SVaclav Hapla PetscInt v; 2811d5c9e23SVaclav Hapla 2821d5c9e23SVaclav Hapla if (p < 0) continue; 2839566063dSJacob Faibussowitsch PetscCall(DMLabelHasPoint(label, p, &has)); 2841d5c9e23SVaclav Hapla if (!has) { 28563a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label \"%s\" does not have point %" PetscInt_FMT "\n", rank, labelName, p)); 2861d5c9e23SVaclav Hapla fail = PETSC_TRUE; 2871d5c9e23SVaclav Hapla continue; 2881d5c9e23SVaclav Hapla } 2899566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, p, &v)); 2905c820da0SVaclav Hapla if (v != value) { 29163a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %" PetscInt_FMT " has bad value %" PetscInt_FMT " in label \"%s\"", p, v, labelName)); 2921d5c9e23SVaclav Hapla fail = PETSC_TRUE; 2931d5c9e23SVaclav Hapla continue; 2941d5c9e23SVaclav Hapla } 29563a3b9bcSJacob Faibussowitsch if (verbose > 1) PetscCall(PetscSynchronizedPrintf(comm, "[%d] OK point %" PetscInt_FMT "\n", rank, p)); 2961d5c9e23SVaclav Hapla } 2979566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 2989566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR)); 2999566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(pointsIS, &points)); 3009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointsIS)); 3019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm)); 3025c820da0SVaclav Hapla PetscCheck(!fail,comm, PETSC_ERR_PLIB, "Label \"%s\" was not loaded correctly%s", labelName, verbose ? " - see details above" : ""); 3035c820da0SVaclav Hapla PetscFunctionReturn(0); 3045c820da0SVaclav Hapla } 3055c820da0SVaclav Hapla 3065c820da0SVaclav Hapla static PetscErrorCode CheckNumLabels(DM dm, AppCtx *ctx) 3075c820da0SVaclav Hapla { 3085c820da0SVaclav Hapla PetscInt actualNum; 3095c820da0SVaclav Hapla PetscBool fail = PETSC_FALSE; 3105c820da0SVaclav Hapla MPI_Comm comm; 3115c820da0SVaclav Hapla PetscMPIInt rank; 3125c820da0SVaclav Hapla 3135c820da0SVaclav Hapla PetscFunctionBeginUser; 3145c820da0SVaclav Hapla if (ctx->num_labels < 0) PetscFunctionReturn(0); 3155c820da0SVaclav Hapla PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3165c820da0SVaclav Hapla PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3175c820da0SVaclav Hapla PetscCall(DMGetNumLabels(dm, &actualNum)); 3185c820da0SVaclav Hapla if (ctx->num_labels != actualNum) { 3195c820da0SVaclav Hapla fail = PETSC_TRUE; 3205c820da0SVaclav Hapla if (ctx->verbose) { 3215c820da0SVaclav Hapla PetscInt i; 3225c820da0SVaclav Hapla 32363a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Asserted number of labels: %" PetscInt_FMT ", actual: %" PetscInt_FMT "\n", rank, ctx->num_labels, actualNum)); 3245c820da0SVaclav Hapla for (i=0; i<actualNum; i++) { 3255c820da0SVaclav Hapla DMLabel label; 3265c820da0SVaclav Hapla const char *name; 3275c820da0SVaclav Hapla 3285c820da0SVaclav Hapla PetscCall(DMGetLabelByNum(dm, i, &label)); 3295c820da0SVaclav Hapla PetscCall(PetscObjectGetName((PetscObject)label, &name)); 33063a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label %" PetscInt_FMT " \"%s\"\n", rank, i, name)); 3315c820da0SVaclav Hapla } 3325c820da0SVaclav Hapla PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR)); 3335c820da0SVaclav Hapla } 3345c820da0SVaclav Hapla } 3355c820da0SVaclav Hapla PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm)); 3365c820da0SVaclav Hapla PetscCheck(!fail, comm, PETSC_ERR_PLIB, "Wrong number of labels%s", ctx->verbose ? " - see details above" : ""); 3375c820da0SVaclav Hapla PetscFunctionReturn(0); 3385c820da0SVaclav Hapla } 3395c820da0SVaclav Hapla 3405c820da0SVaclav Hapla static inline PetscErrorCode IncrementNumLabels(AppCtx *ctx) 3415c820da0SVaclav Hapla { 3425c820da0SVaclav Hapla PetscFunctionBeginUser; 3435c820da0SVaclav Hapla if (ctx->num_labels >= 0) ctx->num_labels++; 3441d5c9e23SVaclav Hapla PetscFunctionReturn(0); 3451d5c9e23SVaclav Hapla } 3461d5c9e23SVaclav Hapla 3471d5c9e23SVaclav Hapla int main(int argc, char **argv) 3481d5c9e23SVaclav Hapla { 3495c820da0SVaclav Hapla const char BOUNDARY_NAME[] = "Boundary"; 3505c820da0SVaclav Hapla const char BOUNDARY_VERTICES_NAME[] = "BoundaryVertices"; 3515c820da0SVaclav Hapla const PetscInt BOUNDARY_VALUE = 12345; 3525c820da0SVaclav Hapla const PetscInt BOUNDARY_VERTICES_VALUE = 6789; 3531d5c9e23SVaclav Hapla DM dm, dmnew; 3541d5c9e23SVaclav Hapla AppCtx user; 3551d5c9e23SVaclav Hapla Vec allCoords = NULL; 3561d5c9e23SVaclav Hapla 3579566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3589566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 3599566063dSJacob Faibussowitsch PetscCall(CreateMesh(&user, &dm)); 3605c820da0SVaclav Hapla PetscCall(MarkBoundary(dm, BOUNDARY_NAME, BOUNDARY_VALUE, PETSC_FALSE, NULL)); 3615c820da0SVaclav Hapla PetscCall(IncrementNumLabels(&user)); 3621d5c9e23SVaclav Hapla if (user.compare_boundary) { 3635c820da0SVaclav Hapla DMLabel label; 3645c820da0SVaclav Hapla 3655c820da0SVaclav Hapla PetscCall(MarkBoundary(dm, BOUNDARY_VERTICES_NAME, BOUNDARY_VERTICES_VALUE, PETSC_TRUE, &label)); 3665c820da0SVaclav Hapla PetscCall(IncrementNumLabels(&user)); 3675c820da0SVaclav Hapla PetscCall(DMLabelGetCoordinateRepresentation(dm, label, BOUNDARY_VERTICES_VALUE, &allCoords)); 3685c820da0SVaclav Hapla PetscCall(DMLabelDestroy(&label)); 3691d5c9e23SVaclav Hapla } 3709566063dSJacob Faibussowitsch PetscCall(SaveMesh(&user, dm)); 3715c820da0SVaclav Hapla 3729566063dSJacob Faibussowitsch PetscCall(LoadMesh(&user, &dmnew)); 3735c820da0SVaclav Hapla PetscCall(IncrementNumLabels(&user)); /* account for depth label */ 3745c820da0SVaclav Hapla PetscCall(IncrementNumLabels(&user)); /* account for celltype label */ 3755c820da0SVaclav Hapla PetscCall(CheckNumLabels(dm, &user)); 3769566063dSJacob Faibussowitsch PetscCall(CompareMeshes(&user, dm, dmnew)); 3771d5c9e23SVaclav Hapla if (user.compare_boundary) { 3785c820da0SVaclav Hapla DMLabel label; 3795c820da0SVaclav Hapla 3805c820da0SVaclav Hapla PetscCall(DMGetLabel(dmnew, BOUNDARY_VERTICES_NAME, &label)); 3815c820da0SVaclav Hapla PetscCall(DMLabelCompareWithCoordinateRepresentation(dmnew, label, BOUNDARY_VERTICES_VALUE, allCoords, user.verbose)); 3821d5c9e23SVaclav Hapla } 3839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmnew)); 3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&allCoords)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 387b122ec5aSJacob Faibussowitsch return 0; 3881d5c9e23SVaclav Hapla } 3891d5c9e23SVaclav Hapla 3901d5c9e23SVaclav Hapla //TODO we can -compare once the new parallel topology format is in place 3911d5c9e23SVaclav Hapla /*TEST 3921d5c9e23SVaclav Hapla build: 3931d5c9e23SVaclav Hapla requires: hdf5 3941d5c9e23SVaclav Hapla 3951d5c9e23SVaclav Hapla # load old format, save in new format, reload, distribute 3961d5c9e23SVaclav Hapla testset: 3971d5c9e23SVaclav Hapla suffix: 1 3981d5c9e23SVaclav Hapla requires: !complex datafilespath 3991d5c9e23SVaclav Hapla args: -dm_plex_name plex 4001d5c9e23SVaclav Hapla args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0 401e600fa54SMatthew G. Knepley args: -dm_plex_interpolate 4021d5c9e23SVaclav Hapla args: -load_dm_plex_check_all 403e600fa54SMatthew G. Knepley args: -use_low_level_functions {{0 1}} -compare_boundary 4045c820da0SVaclav Hapla args: -num_labels 1 4051d5c9e23SVaclav Hapla args: -outfile ex56_1.h5 4061d5c9e23SVaclav Hapla nsize: {{1 3}} 4071d5c9e23SVaclav Hapla test: 4081d5c9e23SVaclav Hapla suffix: a 4091d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 4101d5c9e23SVaclav Hapla test: 4111d5c9e23SVaclav Hapla suffix: b 4121d5c9e23SVaclav Hapla TODO: broken 4131d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 4141d5c9e23SVaclav Hapla test: 4151d5c9e23SVaclav Hapla suffix: c 4161d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 4171d5c9e23SVaclav Hapla test: 4181d5c9e23SVaclav Hapla suffix: d 4195c820da0SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0 4201d5c9e23SVaclav Hapla test: 4211d5c9e23SVaclav Hapla suffix: e 4221d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 4231d5c9e23SVaclav Hapla test: 4241d5c9e23SVaclav Hapla suffix: f 4251d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 4261d5c9e23SVaclav Hapla 4271d5c9e23SVaclav Hapla # load old format, save in new format, reload topology, distribute, load geometry and labels 4281d5c9e23SVaclav Hapla testset: 4291d5c9e23SVaclav Hapla suffix: 2 4301d5c9e23SVaclav Hapla requires: !complex datafilespath 4311d5c9e23SVaclav Hapla args: -dm_plex_name plex 4321d5c9e23SVaclav Hapla args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0 433e600fa54SMatthew G. Knepley args: -dm_plex_interpolate 4341d5c9e23SVaclav Hapla args: -load_dm_plex_check_all 435e600fa54SMatthew G. Knepley args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary 4365c820da0SVaclav Hapla args: -num_labels 1 4371d5c9e23SVaclav Hapla args: -outfile ex56_2.h5 4381d5c9e23SVaclav Hapla nsize: 3 4391d5c9e23SVaclav Hapla test: 4401d5c9e23SVaclav Hapla suffix: a 4411d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 4421d5c9e23SVaclav Hapla test: 4431d5c9e23SVaclav Hapla suffix: b 4441d5c9e23SVaclav Hapla TODO: broken 4451d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 4461d5c9e23SVaclav Hapla test: 4471d5c9e23SVaclav Hapla suffix: c 4481d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 4491d5c9e23SVaclav Hapla test: 4501d5c9e23SVaclav Hapla suffix: d 4515c820da0SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0 4521d5c9e23SVaclav Hapla test: 4531d5c9e23SVaclav Hapla suffix: e 4541d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 4551d5c9e23SVaclav Hapla test: 4561d5c9e23SVaclav Hapla suffix: f 4571d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 4581d5c9e23SVaclav Hapla 4591d5c9e23SVaclav Hapla # load old format, save in new format, reload topology, distribute, load geometry and labels 4601d5c9e23SVaclav Hapla testset: 4611d5c9e23SVaclav Hapla suffix: 3 4621d5c9e23SVaclav Hapla requires: !complex datafilespath 4631d5c9e23SVaclav Hapla args: -dm_plex_name plex 4641d5c9e23SVaclav Hapla args: -dm_plex_view_hdf5_storage_version 2.0.0 465e600fa54SMatthew G. Knepley args: -dm_plex_interpolate -load_dm_distribute 0 4661d5c9e23SVaclav Hapla args: -use_low_level_functions -compare_pre_post 4675c820da0SVaclav Hapla args: -num_labels 1 4681d5c9e23SVaclav Hapla args: -outfile ex56_3.h5 4691d5c9e23SVaclav Hapla nsize: 3 4701d5c9e23SVaclav Hapla test: 4711d5c9e23SVaclav Hapla suffix: a 4721d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 4731d5c9e23SVaclav Hapla test: 4741d5c9e23SVaclav Hapla suffix: b 4751d5c9e23SVaclav Hapla TODO: broken 4761d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 4771d5c9e23SVaclav Hapla test: 4781d5c9e23SVaclav Hapla suffix: c 4791d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 4801d5c9e23SVaclav Hapla test: 4811d5c9e23SVaclav Hapla suffix: d 4825c820da0SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0 4831d5c9e23SVaclav Hapla test: 4841d5c9e23SVaclav Hapla suffix: e 4851d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 4861d5c9e23SVaclav Hapla test: 4871d5c9e23SVaclav Hapla suffix: f 4881d5c9e23SVaclav Hapla args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 4891d5c9e23SVaclav Hapla TEST*/ 490