xref: /petsc/src/dm/impls/plex/tests/ex56.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
11d5c9e23SVaclav Hapla #include <petscdmplex.h>
21d5c9e23SVaclav Hapla #include <petscviewerhdf5.h>
31d5c9e23SVaclav Hapla #include <petscsf.h>
41d5c9e23SVaclav Hapla 
5e22c0602SVaclav 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 
229371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
231d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
241d5c9e23SVaclav Hapla   options->comm                       = comm;
255c820da0SVaclav Hapla   options->num_labels                 = -1;
261d5c9e23SVaclav Hapla   options->compare                    = PETSC_FALSE;
271d5c9e23SVaclav Hapla   options->compare_labels             = PETSC_FALSE;
281d5c9e23SVaclav Hapla   options->compare_boundary           = PETSC_FALSE;
291d5c9e23SVaclav Hapla   options->compare_pre_post           = PETSC_FALSE;
301d5c9e23SVaclav Hapla   options->outfile[0]                 = '\0';
311d5c9e23SVaclav Hapla   options->use_low_level_functions    = PETSC_FALSE;
321d5c9e23SVaclav Hapla   options->distribute_after_topo_load = PETSC_FALSE;
335c820da0SVaclav Hapla   options->verbose                    = 0;
341d5c9e23SVaclav Hapla 
35d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
365c820da0SVaclav 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));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL));
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL));
409566063dSJacob 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));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL));
429566063dSJacob 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));
439566063dSJacob 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));
445c820da0SVaclav Hapla   PetscCall(PetscOptionsInt("-verbose", "Verbosity level", EX, options->verbose, &options->verbose, NULL));
45d0609cedSBarry Smith   PetscOptionsEnd();
461d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
471d5c9e23SVaclav Hapla };
481d5c9e23SVaclav Hapla 
499371c9d4SSatish Balay static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm) {
501d5c9e23SVaclav Hapla   DM dm;
511d5c9e23SVaclav Hapla 
521d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch   PetscCall(DMCreate(options->comm, &dm));
549566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
559566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &options->meshname));
579566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
581d5c9e23SVaclav Hapla   *newdm = dm;
591d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
601d5c9e23SVaclav Hapla }
611d5c9e23SVaclav Hapla 
629371c9d4SSatish Balay static PetscErrorCode SaveMesh(AppCtx *options, DM dm) {
631d5c9e23SVaclav Hapla   PetscViewer v;
641d5c9e23SVaclav Hapla 
651d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
669566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject)dm), options->outfile, FILE_MODE_WRITE, &v));
671d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
689566063dSJacob Faibussowitsch     PetscCall(DMPlexTopologyView(dm, v));
699566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesView(dm, v));
709566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsView(dm, v));
711d5c9e23SVaclav Hapla   } else {
729566063dSJacob Faibussowitsch     PetscCall(DMView(dm, v));
731d5c9e23SVaclav Hapla   }
749566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
751d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
761d5c9e23SVaclav Hapla }
771d5c9e23SVaclav Hapla 
789371c9d4SSatish Balay typedef enum {
799371c9d4SSatish Balay   NONE      = 0,
809371c9d4SSatish Balay   PRE_DIST  = 1,
819371c9d4SSatish Balay   POST_DIST = 2
829371c9d4SSatish Balay } AuxObjLoadMode;
831d5c9e23SVaclav Hapla 
849371c9d4SSatish Balay static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm) {
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 
1239371c9d4SSatish Balay static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew) {
1241d5c9e23SVaclav Hapla   DM          dm;
1251d5c9e23SVaclav Hapla   PetscViewer v;
1261d5c9e23SVaclav Hapla 
1271d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
1289566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(options->comm, options->outfile, FILE_MODE_READ, &v));
1291d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
1301d5c9e23SVaclav Hapla     if (options->compare_pre_post) {
1311d5c9e23SVaclav Hapla       DM dm0;
1321d5c9e23SVaclav Hapla 
1339566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0));
1349566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm));
1359566063dSJacob Faibussowitsch       PetscCall(DMCompareLabels(dm0, dm, NULL, NULL));
1369566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm0));
1371d5c9e23SVaclav Hapla     } else {
1389566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm));
1391d5c9e23SVaclav Hapla     }
1401d5c9e23SVaclav Hapla   } else {
1419566063dSJacob Faibussowitsch     PetscCall(DMCreate(options->comm, &dm));
1429566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm, DMPLEX));
1439566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)dm, options->meshname));
1449566063dSJacob Faibussowitsch     PetscCall(DMLoad(dm, v));
1451d5c9e23SVaclav Hapla   }
1469566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
1471d5c9e23SVaclav Hapla 
1489566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dm, "load_"));
1499566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1509566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1511d5c9e23SVaclav Hapla   *dmnew = dm;
1521d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1531d5c9e23SVaclav Hapla }
1541d5c9e23SVaclav Hapla 
1559371c9d4SSatish Balay static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1) {
1561d5c9e23SVaclav Hapla   PetscBool flg;
1571d5c9e23SVaclav Hapla 
1581d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
1591d5c9e23SVaclav Hapla   if (options->compare) {
1609566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dm0, dm1, &flg));
16128b400f6SJacob Faibussowitsch     PetscCheck(flg, options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
1629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(options->comm, "DMs equal\n"));
1631d5c9e23SVaclav Hapla   }
1641d5c9e23SVaclav Hapla   if (options->compare_labels) {
1659566063dSJacob Faibussowitsch     PetscCall(DMCompareLabels(dm0, dm1, NULL, NULL));
1669566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(options->comm, "DMLabels equal\n"));
1671d5c9e23SVaclav Hapla   }
1681d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1691d5c9e23SVaclav Hapla }
1701d5c9e23SVaclav Hapla 
1719371c9d4SSatish Balay static PetscErrorCode MarkBoundary(DM dm, const char name[], PetscInt value, PetscBool verticesOnly, DMLabel *label) {
1721d5c9e23SVaclav Hapla   DMLabel l;
1735c820da0SVaclav Hapla 
1745c820da0SVaclav Hapla   PetscFunctionBeginUser;
1755c820da0SVaclav Hapla   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)dm), name, &l));
1765c820da0SVaclav Hapla   PetscCall(DMAddLabel(dm, l));
1775c820da0SVaclav Hapla   PetscCall(DMPlexMarkBoundaryFaces(dm, value, l));
1785c820da0SVaclav Hapla   PetscCall(DMPlexLabelComplete(dm, l));
1795c820da0SVaclav Hapla   if (verticesOnly) {
1801d5c9e23SVaclav Hapla     IS              points;
1811d5c9e23SVaclav Hapla     const PetscInt *idx;
1821d5c9e23SVaclav Hapla     PetscInt        i, n;
1831d5c9e23SVaclav Hapla 
1849566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(l, value, &points));
1859566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(points, &n));
1869566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(points, &idx));
1871d5c9e23SVaclav Hapla     for (i = 0; i < n; i++) {
1881d5c9e23SVaclav Hapla       const PetscInt p = idx[i];
1891d5c9e23SVaclav Hapla       PetscInt       d;
1901d5c9e23SVaclav Hapla 
1919566063dSJacob Faibussowitsch       PetscCall(DMPlexGetPointDepth(dm, p, &d));
192*48a46eb9SPierre Jolivet       if (d != 0) PetscCall(DMLabelClearValue(l, p, value));
1931d5c9e23SVaclav Hapla     }
1949566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(points, &idx));
1959566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&points));
1965c820da0SVaclav Hapla   }
1975c820da0SVaclav Hapla   if (label) *label = l;
1985c820da0SVaclav Hapla   else PetscCall(DMLabelDestroy(&l));
1991d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2001d5c9e23SVaclav Hapla }
2011d5c9e23SVaclav Hapla 
2029371c9d4SSatish Balay static PetscErrorCode VertexCoordinatesToAll(DM dm, IS vertices, Vec *allCoords) {
2031d5c9e23SVaclav Hapla   Vec        coords, allCoords_;
2041d5c9e23SVaclav Hapla   VecScatter sc;
2051d5c9e23SVaclav Hapla   MPI_Comm   comm;
2061d5c9e23SVaclav Hapla 
2071d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2099566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
2101d5c9e23SVaclav Hapla   if (vertices) {
2119566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords));
2121d5c9e23SVaclav Hapla   } else {
2139566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &coords));
2141d5c9e23SVaclav Hapla   }
2151d5c9e23SVaclav Hapla   {
2161d5c9e23SVaclav Hapla     PetscInt n;
2171d5c9e23SVaclav Hapla     Vec      mpivec;
2181d5c9e23SVaclav Hapla 
2199566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coords, &n));
2209566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(comm, n, PETSC_DECIDE, &mpivec));
2219566063dSJacob Faibussowitsch     PetscCall(VecCopy(coords, mpivec));
2229566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&coords));
2231d5c9e23SVaclav Hapla     coords = mpivec;
2241d5c9e23SVaclav Hapla   }
2251d5c9e23SVaclav Hapla 
2269566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToAll(coords, &sc, &allCoords_));
2279566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sc, coords, allCoords_, INSERT_VALUES, SCATTER_FORWARD));
2289566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sc, coords, allCoords_, INSERT_VALUES, SCATTER_FORWARD));
2299566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sc));
2309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coords));
2311d5c9e23SVaclav Hapla   *allCoords = allCoords_;
2321d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2331d5c9e23SVaclav Hapla }
2341d5c9e23SVaclav Hapla 
2359371c9d4SSatish Balay static PetscErrorCode DMLabelGetCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec *allCoords) {
2361d5c9e23SVaclav Hapla   IS vertices;
2371d5c9e23SVaclav Hapla 
2381d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
2395c820da0SVaclav Hapla   PetscCall(DMLabelGetStratumIS(label, value, &vertices));
2409566063dSJacob Faibussowitsch   PetscCall(VertexCoordinatesToAll(dm, vertices, allCoords));
2419566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&vertices));
2421d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2431d5c9e23SVaclav Hapla }
2441d5c9e23SVaclav Hapla 
2459371c9d4SSatish Balay static PetscErrorCode DMLabelCompareWithCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec allCoords, PetscInt verbose) {
2465c820da0SVaclav Hapla   const char     *labelName;
2471d5c9e23SVaclav Hapla   IS              pointsIS;
2481d5c9e23SVaclav Hapla   const PetscInt *points;
2491d5c9e23SVaclav Hapla   PetscInt        i, n;
2501d5c9e23SVaclav Hapla   PetscBool       fail = PETSC_FALSE;
2511d5c9e23SVaclav Hapla   MPI_Comm        comm;
2521d5c9e23SVaclav Hapla   PetscMPIInt     rank;
2531d5c9e23SVaclav Hapla 
2541d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
2555c820da0SVaclav Hapla   PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label not loaded");
2569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2575c820da0SVaclav Hapla   PetscCall(PetscObjectGetName((PetscObject)label, &labelName));
2589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2591d5c9e23SVaclav Hapla   {
2601d5c9e23SVaclav Hapla     PetscInt pStart, pEnd;
2611d5c9e23SVaclav Hapla 
2629566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2639566063dSJacob Faibussowitsch     PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
2641d5c9e23SVaclav Hapla   }
2659566063dSJacob Faibussowitsch   PetscCall(DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS));
2669566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointsIS, &points));
2679566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointsIS, &n));
2685c820da0SVaclav Hapla   if (verbose > 1) PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_(comm)));
2691d5c9e23SVaclav Hapla   for (i = 0; i < n; i++) {
2701d5c9e23SVaclav Hapla     const PetscInt p = points[i];
2711d5c9e23SVaclav Hapla     PetscBool      has;
2721d5c9e23SVaclav Hapla     PetscInt       v;
2731d5c9e23SVaclav Hapla 
2741d5c9e23SVaclav Hapla     if (p < 0) continue;
2759566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(label, p, &has));
2761d5c9e23SVaclav Hapla     if (!has) {
27763a3b9bcSJacob Faibussowitsch       if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label \"%s\" does not have point %" PetscInt_FMT "\n", rank, labelName, p));
2781d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2791d5c9e23SVaclav Hapla       continue;
2801d5c9e23SVaclav Hapla     }
2819566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(label, p, &v));
2825c820da0SVaclav Hapla     if (v != value) {
28363a3b9bcSJacob Faibussowitsch       if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %" PetscInt_FMT " has bad value %" PetscInt_FMT " in label \"%s\"", p, v, labelName));
2841d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2851d5c9e23SVaclav Hapla       continue;
2861d5c9e23SVaclav Hapla     }
28763a3b9bcSJacob Faibussowitsch     if (verbose > 1) PetscCall(PetscSynchronizedPrintf(comm, "[%d] OK point %" PetscInt_FMT "\n", rank, p));
2881d5c9e23SVaclav Hapla   }
2899566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
2909566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
2919566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointsIS, &points));
2929566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointsIS));
2939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
2945c820da0SVaclav Hapla   PetscCheck(!fail, comm, PETSC_ERR_PLIB, "Label \"%s\" was not loaded correctly%s", labelName, verbose ? " - see details above" : "");
2955c820da0SVaclav Hapla   PetscFunctionReturn(0);
2965c820da0SVaclav Hapla }
2975c820da0SVaclav Hapla 
2989371c9d4SSatish Balay static PetscErrorCode CheckNumLabels(DM dm, AppCtx *ctx) {
2995c820da0SVaclav Hapla   PetscInt    actualNum;
3005c820da0SVaclav Hapla   PetscBool   fail = PETSC_FALSE;
3015c820da0SVaclav Hapla   MPI_Comm    comm;
3025c820da0SVaclav Hapla   PetscMPIInt rank;
3035c820da0SVaclav Hapla 
3045c820da0SVaclav Hapla   PetscFunctionBeginUser;
3055c820da0SVaclav Hapla   if (ctx->num_labels < 0) PetscFunctionReturn(0);
3065c820da0SVaclav Hapla   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3075c820da0SVaclav Hapla   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3085c820da0SVaclav Hapla   PetscCall(DMGetNumLabels(dm, &actualNum));
3095c820da0SVaclav Hapla   if (ctx->num_labels != actualNum) {
3105c820da0SVaclav Hapla     fail = PETSC_TRUE;
3115c820da0SVaclav Hapla     if (ctx->verbose) {
3125c820da0SVaclav Hapla       PetscInt i;
3135c820da0SVaclav Hapla 
31463a3b9bcSJacob Faibussowitsch       PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Asserted number of labels: %" PetscInt_FMT ", actual: %" PetscInt_FMT "\n", rank, ctx->num_labels, actualNum));
3155c820da0SVaclav Hapla       for (i = 0; i < actualNum; i++) {
3165c820da0SVaclav Hapla         DMLabel     label;
3175c820da0SVaclav Hapla         const char *name;
3185c820da0SVaclav Hapla 
3195c820da0SVaclav Hapla         PetscCall(DMGetLabelByNum(dm, i, &label));
3205c820da0SVaclav Hapla         PetscCall(PetscObjectGetName((PetscObject)label, &name));
32163a3b9bcSJacob Faibussowitsch         PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label %" PetscInt_FMT " \"%s\"\n", rank, i, name));
3225c820da0SVaclav Hapla       }
3235c820da0SVaclav Hapla       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
3245c820da0SVaclav Hapla     }
3255c820da0SVaclav Hapla   }
3265c820da0SVaclav Hapla   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
3275c820da0SVaclav Hapla   PetscCheck(!fail, comm, PETSC_ERR_PLIB, "Wrong number of labels%s", ctx->verbose ? " - see details above" : "");
3285c820da0SVaclav Hapla   PetscFunctionReturn(0);
3295c820da0SVaclav Hapla }
3305c820da0SVaclav Hapla 
3319371c9d4SSatish Balay static inline PetscErrorCode IncrementNumLabels(AppCtx *ctx) {
3325c820da0SVaclav Hapla   PetscFunctionBeginUser;
3335c820da0SVaclav Hapla   if (ctx->num_labels >= 0) ctx->num_labels++;
3341d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
3351d5c9e23SVaclav Hapla }
3361d5c9e23SVaclav Hapla 
3379371c9d4SSatish Balay int main(int argc, char **argv) {
3385c820da0SVaclav Hapla   const char     BOUNDARY_NAME[]          = "Boundary";
3395c820da0SVaclav Hapla   const char     BOUNDARY_VERTICES_NAME[] = "BoundaryVertices";
3405c820da0SVaclav Hapla   const PetscInt BOUNDARY_VALUE           = 12345;
3415c820da0SVaclav Hapla   const PetscInt BOUNDARY_VERTICES_VALUE  = 6789;
3421d5c9e23SVaclav Hapla   DM             dm, dmnew;
3431d5c9e23SVaclav Hapla   AppCtx         user;
3441d5c9e23SVaclav Hapla   Vec            allCoords = NULL;
3451d5c9e23SVaclav Hapla 
346327415f7SBarry Smith   PetscFunctionBeginUser;
3479566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3489566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3499566063dSJacob Faibussowitsch   PetscCall(CreateMesh(&user, &dm));
3505c820da0SVaclav Hapla   PetscCall(MarkBoundary(dm, BOUNDARY_NAME, BOUNDARY_VALUE, PETSC_FALSE, NULL));
3515c820da0SVaclav Hapla   PetscCall(IncrementNumLabels(&user));
3521d5c9e23SVaclav Hapla   if (user.compare_boundary) {
3535c820da0SVaclav Hapla     DMLabel label;
3545c820da0SVaclav Hapla 
3555c820da0SVaclav Hapla     PetscCall(MarkBoundary(dm, BOUNDARY_VERTICES_NAME, BOUNDARY_VERTICES_VALUE, PETSC_TRUE, &label));
3565c820da0SVaclav Hapla     PetscCall(IncrementNumLabels(&user));
3575c820da0SVaclav Hapla     PetscCall(DMLabelGetCoordinateRepresentation(dm, label, BOUNDARY_VERTICES_VALUE, &allCoords));
3585c820da0SVaclav Hapla     PetscCall(DMLabelDestroy(&label));
3591d5c9e23SVaclav Hapla   }
3609566063dSJacob Faibussowitsch   PetscCall(SaveMesh(&user, dm));
3615c820da0SVaclav Hapla 
3629566063dSJacob Faibussowitsch   PetscCall(LoadMesh(&user, &dmnew));
3635c820da0SVaclav Hapla   PetscCall(IncrementNumLabels(&user)); /* account for depth label */
3645c820da0SVaclav Hapla   PetscCall(IncrementNumLabels(&user)); /* account for celltype label */
3655c820da0SVaclav Hapla   PetscCall(CheckNumLabels(dm, &user));
3669566063dSJacob Faibussowitsch   PetscCall(CompareMeshes(&user, dm, dmnew));
3671d5c9e23SVaclav Hapla   if (user.compare_boundary) {
3685c820da0SVaclav Hapla     DMLabel label;
3695c820da0SVaclav Hapla 
3705c820da0SVaclav Hapla     PetscCall(DMGetLabel(dmnew, BOUNDARY_VERTICES_NAME, &label));
3715c820da0SVaclav Hapla     PetscCall(DMLabelCompareWithCoordinateRepresentation(dmnew, label, BOUNDARY_VERTICES_VALUE, allCoords, user.verbose));
3721d5c9e23SVaclav Hapla   }
3739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmnew));
3759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&allCoords));
3769566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
377b122ec5aSJacob Faibussowitsch   return 0;
3781d5c9e23SVaclav Hapla }
3791d5c9e23SVaclav Hapla 
3801d5c9e23SVaclav Hapla //TODO we can -compare once the new parallel topology format is in place
3811d5c9e23SVaclav Hapla /*TEST
3821d5c9e23SVaclav Hapla   build:
3831d5c9e23SVaclav Hapla     requires: hdf5
3841d5c9e23SVaclav Hapla 
3851d5c9e23SVaclav Hapla   # load old format, save in new format, reload, distribute
3861d5c9e23SVaclav Hapla   testset:
3871d5c9e23SVaclav Hapla     suffix: 1
3881d5c9e23SVaclav Hapla     requires: !complex datafilespath
3891d5c9e23SVaclav Hapla     args: -dm_plex_name plex
3901d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
391e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
3921d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
393e600fa54SMatthew G. Knepley     args: -use_low_level_functions {{0 1}} -compare_boundary
3945c820da0SVaclav Hapla     args: -num_labels 1
3951d5c9e23SVaclav Hapla     args: -outfile ex56_1.h5
3961d5c9e23SVaclav Hapla     nsize: {{1 3}}
3971d5c9e23SVaclav Hapla     test:
3981d5c9e23SVaclav Hapla       suffix: a
3991d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4001d5c9e23SVaclav Hapla     test:
4011d5c9e23SVaclav Hapla       suffix: b
4021d5c9e23SVaclav Hapla       TODO: broken
4031d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4041d5c9e23SVaclav Hapla     test:
4051d5c9e23SVaclav Hapla       suffix: c
4061d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4071d5c9e23SVaclav Hapla     test:
4081d5c9e23SVaclav Hapla       suffix: d
4095c820da0SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
4101d5c9e23SVaclav Hapla     test:
4111d5c9e23SVaclav Hapla       suffix: e
4121d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4131d5c9e23SVaclav Hapla     test:
4141d5c9e23SVaclav Hapla       suffix: f
4151d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4161d5c9e23SVaclav Hapla 
4171d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
4181d5c9e23SVaclav Hapla   testset:
4191d5c9e23SVaclav Hapla     suffix: 2
4201d5c9e23SVaclav Hapla     requires: !complex datafilespath
4211d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4221d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
423e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
4241d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
425e600fa54SMatthew G. Knepley     args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary
4265c820da0SVaclav Hapla     args: -num_labels 1
4271d5c9e23SVaclav Hapla     args: -outfile ex56_2.h5
4281d5c9e23SVaclav Hapla     nsize: 3
4291d5c9e23SVaclav Hapla     test:
4301d5c9e23SVaclav Hapla       suffix: a
4311d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4321d5c9e23SVaclav Hapla     test:
4331d5c9e23SVaclav Hapla       suffix: b
4341d5c9e23SVaclav Hapla       TODO: broken
4351d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4361d5c9e23SVaclav Hapla     test:
4371d5c9e23SVaclav Hapla       suffix: c
4381d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4391d5c9e23SVaclav Hapla     test:
4401d5c9e23SVaclav Hapla       suffix: d
4415c820da0SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
4421d5c9e23SVaclav Hapla     test:
4431d5c9e23SVaclav Hapla       suffix: e
4441d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4451d5c9e23SVaclav Hapla     test:
4461d5c9e23SVaclav Hapla       suffix: f
4471d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4481d5c9e23SVaclav Hapla 
4491d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
4501d5c9e23SVaclav Hapla   testset:
4511d5c9e23SVaclav Hapla     suffix: 3
4521d5c9e23SVaclav Hapla     requires: !complex datafilespath
4531d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4541d5c9e23SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 2.0.0
455e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate -load_dm_distribute 0
4561d5c9e23SVaclav Hapla     args: -use_low_level_functions -compare_pre_post
4575c820da0SVaclav Hapla     args: -num_labels 1
4581d5c9e23SVaclav Hapla     args: -outfile ex56_3.h5
4591d5c9e23SVaclav Hapla     nsize: 3
4601d5c9e23SVaclav Hapla     test:
4611d5c9e23SVaclav Hapla       suffix: a
4621d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4631d5c9e23SVaclav Hapla     test:
4641d5c9e23SVaclav Hapla       suffix: b
4651d5c9e23SVaclav Hapla       TODO: broken
4661d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4671d5c9e23SVaclav Hapla     test:
4681d5c9e23SVaclav Hapla       suffix: c
4691d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4701d5c9e23SVaclav Hapla     test:
4711d5c9e23SVaclav Hapla       suffix: d
4725c820da0SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
4731d5c9e23SVaclav Hapla     test:
4741d5c9e23SVaclav Hapla       suffix: e
4751d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4761d5c9e23SVaclav Hapla     test:
4771d5c9e23SVaclav Hapla       suffix: f
4781d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4791d5c9e23SVaclav Hapla TEST*/
480