xref: /petsc/src/dm/impls/plex/tests/ex56.c (revision 5c820da0a866dfd2d2d382a86bfaa6df0df955c6)
11d5c9e23SVaclav Hapla #include <petscdmplex.h>
21d5c9e23SVaclav Hapla #include <petscviewerhdf5.h>
31d5c9e23SVaclav Hapla #include <petscsf.h>
41d5c9e23SVaclav Hapla 
51d5c9e23SVaclav Hapla static const char help[] = "Load and save the mesh to the native HDF5 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 */
10*5c820da0SVaclav 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 */
19*5c820da0SVaclav Hapla   PetscInt    verbose;
201d5c9e23SVaclav Hapla } AppCtx;
211d5c9e23SVaclav Hapla 
221d5c9e23SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
231d5c9e23SVaclav Hapla {
241d5c9e23SVaclav Hapla   PetscErrorCode ierr;
251d5c9e23SVaclav Hapla 
261d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
271d5c9e23SVaclav Hapla   options->comm                       = comm;
28*5c820da0SVaclav Hapla   options->num_labels                 = -1;
291d5c9e23SVaclav Hapla   options->compare                    = PETSC_FALSE;
301d5c9e23SVaclav Hapla   options->compare_labels             = PETSC_FALSE;
311d5c9e23SVaclav Hapla   options->compare_boundary           = PETSC_FALSE;
321d5c9e23SVaclav Hapla   options->compare_pre_post           = PETSC_FALSE;
331d5c9e23SVaclav Hapla   options->outfile[0]                 = '\0';
341d5c9e23SVaclav Hapla   options->use_low_level_functions    = PETSC_FALSE;
351d5c9e23SVaclav Hapla   options->distribute_after_topo_load = PETSC_FALSE;
36*5c820da0SVaclav Hapla   options->verbose                    = 0;
371d5c9e23SVaclav Hapla 
389566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");PetscCall(ierr);
39*5c820da0SVaclav 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));
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL));
439566063dSJacob 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));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL));
459566063dSJacob 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));
469566063dSJacob 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));
47*5c820da0SVaclav Hapla   PetscCall(PetscOptionsInt("-verbose", "Verbosity level", EX, options->verbose, &options->verbose, NULL));
489566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
491d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
501d5c9e23SVaclav Hapla };
511d5c9e23SVaclav Hapla 
521d5c9e23SVaclav Hapla static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm)
531d5c9e23SVaclav Hapla {
541d5c9e23SVaclav Hapla   DM             dm;
551d5c9e23SVaclav Hapla 
561d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
579566063dSJacob Faibussowitsch   PetscCall(DMCreate(options->comm, &dm));
589566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
599566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &options->meshname));
619566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
621d5c9e23SVaclav Hapla   *newdm = dm;
631d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
641d5c9e23SVaclav Hapla }
651d5c9e23SVaclav Hapla 
661d5c9e23SVaclav Hapla static PetscErrorCode SaveMesh(AppCtx *options, DM dm)
671d5c9e23SVaclav Hapla {
681d5c9e23SVaclav Hapla   PetscViewer    v;
691d5c9e23SVaclav Hapla 
701d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), options->outfile, FILE_MODE_WRITE, &v));
721d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
739566063dSJacob Faibussowitsch     PetscCall(DMPlexTopologyView(dm, v));
749566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesView(dm, v));
759566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsView(dm, v));
761d5c9e23SVaclav Hapla   } else {
779566063dSJacob Faibussowitsch     PetscCall(DMView(dm, v));
781d5c9e23SVaclav Hapla   }
799566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
801d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
811d5c9e23SVaclav Hapla }
821d5c9e23SVaclav Hapla 
831d5c9e23SVaclav Hapla typedef enum {NONE=0, PRE_DIST=1, POST_DIST=2} AuxObjLoadMode;
841d5c9e23SVaclav Hapla 
851d5c9e23SVaclav Hapla static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm)
861d5c9e23SVaclav Hapla {
871d5c9e23SVaclav Hapla   DM              dm;
881d5c9e23SVaclav Hapla   PetscSF         sfXC;
891d5c9e23SVaclav Hapla 
901d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
919566063dSJacob Faibussowitsch   PetscCall(DMCreate(options->comm, &dm));
929566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname));
949566063dSJacob Faibussowitsch   PetscCall(DMPlexTopologyLoad(dm, v, &sfXC));
951d5c9e23SVaclav Hapla   if (mode == PRE_DIST) {
969566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC));
979566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsLoad(dm, v, sfXC));
981d5c9e23SVaclav Hapla   }
991d5c9e23SVaclav Hapla   if (explicitDistribute) {
1001d5c9e23SVaclav Hapla     DM      dmdist;
1011d5c9e23SVaclav Hapla     PetscSF sfXB = sfXC, sfBC;
1021d5c9e23SVaclav Hapla 
1039566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(dm, 0, &sfBC, &dmdist));
1041d5c9e23SVaclav Hapla     if (dmdist) {
1051d5c9e23SVaclav Hapla       const char *name;
1061d5c9e23SVaclav Hapla 
1079566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) dm, &name));
1089566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) dmdist, name));
1099566063dSJacob Faibussowitsch       PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC));
1109566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&sfXB));
1119566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&sfBC));
1129566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
1131d5c9e23SVaclav Hapla       dm   = dmdist;
1141d5c9e23SVaclav Hapla     }
1151d5c9e23SVaclav Hapla   }
1161d5c9e23SVaclav Hapla   if (mode == POST_DIST) {
1179566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC));
1189566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsLoad(dm, v, sfXC));
1191d5c9e23SVaclav Hapla   }
1209566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfXC));
1211d5c9e23SVaclav Hapla   *newdm = dm;
1221d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1231d5c9e23SVaclav Hapla }
1241d5c9e23SVaclav Hapla 
1251d5c9e23SVaclav Hapla static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew)
1261d5c9e23SVaclav Hapla {
1271d5c9e23SVaclav Hapla   DM             dm;
1281d5c9e23SVaclav Hapla   PetscViewer    v;
1291d5c9e23SVaclav Hapla 
1301d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
1319566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(options->comm, options->outfile, FILE_MODE_READ, &v));
1321d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
1331d5c9e23SVaclav Hapla     if (options->compare_pre_post) {
1341d5c9e23SVaclav Hapla       DM dm0;
1351d5c9e23SVaclav Hapla 
1369566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0));
1379566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm));
1389566063dSJacob Faibussowitsch       PetscCall(DMCompareLabels(dm0, dm, NULL, NULL));
1399566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm0));
1401d5c9e23SVaclav Hapla     } else {
1419566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm));
1421d5c9e23SVaclav Hapla     }
1431d5c9e23SVaclav Hapla   } else {
1449566063dSJacob Faibussowitsch     PetscCall(DMCreate(options->comm, &dm));
1459566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm, DMPLEX));
1469566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname));
1479566063dSJacob Faibussowitsch     PetscCall(DMLoad(dm, v));
1481d5c9e23SVaclav Hapla   }
1499566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
1501d5c9e23SVaclav Hapla 
1519566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dm, "load_"));
1529566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1539566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1541d5c9e23SVaclav Hapla   *dmnew = dm;
1551d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1561d5c9e23SVaclav Hapla }
1571d5c9e23SVaclav Hapla 
1581d5c9e23SVaclav Hapla static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1)
1591d5c9e23SVaclav Hapla {
1601d5c9e23SVaclav Hapla   PetscBool       flg;
1611d5c9e23SVaclav Hapla 
1621d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
1631d5c9e23SVaclav Hapla   if (options->compare) {
1649566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dm0, dm1, &flg));
16528b400f6SJacob Faibussowitsch     PetscCheck(flg,options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
1669566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(options->comm,"DMs equal\n"));
1671d5c9e23SVaclav Hapla   }
1681d5c9e23SVaclav Hapla   if (options->compare_labels) {
1699566063dSJacob Faibussowitsch     PetscCall(DMCompareLabels(dm0, dm1, NULL, NULL));
1709566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(options->comm,"DMLabels equal\n"));
1711d5c9e23SVaclav Hapla   }
1721d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1731d5c9e23SVaclav Hapla }
1741d5c9e23SVaclav Hapla 
175*5c820da0SVaclav Hapla static PetscErrorCode MarkBoundary(DM dm, const char name[], PetscInt value, PetscBool verticesOnly, DMLabel *label)
1761d5c9e23SVaclav Hapla {
1771d5c9e23SVaclav Hapla   DMLabel         l;
178*5c820da0SVaclav Hapla 
179*5c820da0SVaclav Hapla   PetscFunctionBeginUser;
180*5c820da0SVaclav Hapla   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)dm), name, &l));
181*5c820da0SVaclav Hapla   PetscCall(DMAddLabel(dm, l));
182*5c820da0SVaclav Hapla   PetscCall(DMPlexMarkBoundaryFaces(dm, value, l));
183*5c820da0SVaclav Hapla   PetscCall(DMPlexLabelComplete(dm, l));
184*5c820da0SVaclav Hapla   if (verticesOnly) {
1851d5c9e23SVaclav Hapla     IS              points;
1861d5c9e23SVaclav Hapla     const PetscInt *idx;
1871d5c9e23SVaclav Hapla     PetscInt        i, n;
1881d5c9e23SVaclav Hapla 
1899566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(l, value, &points));
1909566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(points, &n));
1919566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(points, &idx));
1921d5c9e23SVaclav Hapla     for (i=0; i<n; i++) {
1931d5c9e23SVaclav Hapla       const PetscInt p = idx[i];
1941d5c9e23SVaclav Hapla       PetscInt       d;
1951d5c9e23SVaclav Hapla 
1969566063dSJacob Faibussowitsch       PetscCall(DMPlexGetPointDepth(dm, p, &d));
1971d5c9e23SVaclav Hapla       if (d != 0) {
1989566063dSJacob Faibussowitsch         PetscCall(DMLabelClearValue(l, p, value));
1991d5c9e23SVaclav Hapla       }
2001d5c9e23SVaclav Hapla     }
2019566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(points, &idx));
2029566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&points));
203*5c820da0SVaclav Hapla   }
204*5c820da0SVaclav Hapla   if (label) *label = l;
205*5c820da0SVaclav Hapla   else PetscCall(DMLabelDestroy(&l));
2061d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2071d5c9e23SVaclav Hapla }
2081d5c9e23SVaclav Hapla 
2091d5c9e23SVaclav Hapla static PetscErrorCode VertexCoordinatesToAll(DM dm, IS vertices, Vec *allCoords)
2101d5c9e23SVaclav Hapla {
2111d5c9e23SVaclav Hapla   Vec             coords, allCoords_;
2121d5c9e23SVaclav Hapla   VecScatter      sc;
2131d5c9e23SVaclav Hapla   MPI_Comm        comm;
2141d5c9e23SVaclav Hapla 
2151d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
2169566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2179566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
2181d5c9e23SVaclav Hapla   if (vertices) {
2199566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords));
2201d5c9e23SVaclav Hapla   } else {
2219566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &coords));
2221d5c9e23SVaclav Hapla   }
2231d5c9e23SVaclav Hapla   {
2241d5c9e23SVaclav Hapla     PetscInt  n;
2251d5c9e23SVaclav Hapla     Vec       mpivec;
2261d5c9e23SVaclav Hapla 
2279566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coords, &n));
2289566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(comm, n, PETSC_DECIDE, &mpivec));
2299566063dSJacob Faibussowitsch     PetscCall(VecCopy(coords, mpivec));
2309566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&coords));
2311d5c9e23SVaclav Hapla     coords = mpivec;
2321d5c9e23SVaclav Hapla   }
2331d5c9e23SVaclav Hapla 
2349566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToAll(coords, &sc, &allCoords_));
2359566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD));
2369566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD));
2379566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sc));
2389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coords));
2391d5c9e23SVaclav Hapla   *allCoords = allCoords_;
2401d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2411d5c9e23SVaclav Hapla }
2421d5c9e23SVaclav Hapla 
243*5c820da0SVaclav Hapla static PetscErrorCode DMLabelGetCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec *allCoords)
2441d5c9e23SVaclav Hapla {
2451d5c9e23SVaclav Hapla   IS vertices;
2461d5c9e23SVaclav Hapla 
2471d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
248*5c820da0SVaclav Hapla   PetscCall(DMLabelGetStratumIS(label, value, &vertices));
2499566063dSJacob Faibussowitsch   PetscCall(VertexCoordinatesToAll(dm, vertices, allCoords));
2509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&vertices));
2511d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2521d5c9e23SVaclav Hapla }
2531d5c9e23SVaclav Hapla 
254*5c820da0SVaclav Hapla static PetscErrorCode DMLabelCompareWithCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec allCoords, PetscInt verbose)
2551d5c9e23SVaclav Hapla {
256*5c820da0SVaclav Hapla   const char     *labelName;
2571d5c9e23SVaclav Hapla   IS              pointsIS;
2581d5c9e23SVaclav Hapla   const PetscInt *points;
2591d5c9e23SVaclav Hapla   PetscInt        i, n;
2601d5c9e23SVaclav Hapla   PetscBool       fail = PETSC_FALSE;
2611d5c9e23SVaclav Hapla   MPI_Comm        comm;
2621d5c9e23SVaclav Hapla   PetscMPIInt     rank;
2631d5c9e23SVaclav Hapla 
2641d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
265*5c820da0SVaclav Hapla   PetscCheck(label,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label not loaded");
2669566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
267*5c820da0SVaclav Hapla   PetscCall(PetscObjectGetName((PetscObject)label, &labelName));
2689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2691d5c9e23SVaclav Hapla   {
2701d5c9e23SVaclav Hapla     PetscInt pStart, pEnd;
2711d5c9e23SVaclav Hapla 
2729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2739566063dSJacob Faibussowitsch     PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
2741d5c9e23SVaclav Hapla   }
2759566063dSJacob Faibussowitsch   PetscCall(DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS));
2769566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointsIS, &points));
2779566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointsIS, &n));
278*5c820da0SVaclav Hapla   if (verbose > 1) PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_(comm)));
2791d5c9e23SVaclav Hapla   for (i=0; i<n; i++) {
2801d5c9e23SVaclav Hapla     const PetscInt  p = points[i];
2811d5c9e23SVaclav Hapla     PetscBool       has;
2821d5c9e23SVaclav Hapla     PetscInt        v;
2831d5c9e23SVaclav Hapla 
2841d5c9e23SVaclav Hapla     if (p < 0) continue;
2859566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(label, p, &has));
2861d5c9e23SVaclav Hapla     if (!has) {
287*5c820da0SVaclav Hapla       if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label \"%s\" does not have point %D\n", rank, labelName, p));
2881d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2891d5c9e23SVaclav Hapla       continue;
2901d5c9e23SVaclav Hapla     }
2919566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(label, p, &v));
292*5c820da0SVaclav Hapla     if (v != value) {
293*5c820da0SVaclav Hapla       if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %D has bad value %D in label \"%s\"", p, v, labelName));
2941d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2951d5c9e23SVaclav Hapla       continue;
2961d5c9e23SVaclav Hapla     }
297*5c820da0SVaclav Hapla     if (verbose > 1) PetscCall(PetscSynchronizedPrintf(comm, "[%d] OK point %D\n", rank, p));
2981d5c9e23SVaclav Hapla   }
2999566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
3009566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
3019566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointsIS, &points));
3029566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointsIS));
3039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
304*5c820da0SVaclav Hapla   PetscCheck(!fail,comm, PETSC_ERR_PLIB, "Label \"%s\" was not loaded correctly%s", labelName, verbose ? " - see details above" : "");
305*5c820da0SVaclav Hapla   PetscFunctionReturn(0);
306*5c820da0SVaclav Hapla }
307*5c820da0SVaclav Hapla 
308*5c820da0SVaclav Hapla static PetscErrorCode CheckNumLabels(DM dm, AppCtx *ctx)
309*5c820da0SVaclav Hapla {
310*5c820da0SVaclav Hapla   PetscInt    actualNum;
311*5c820da0SVaclav Hapla   PetscBool   fail = PETSC_FALSE;
312*5c820da0SVaclav Hapla   MPI_Comm    comm;
313*5c820da0SVaclav Hapla   PetscMPIInt rank;
314*5c820da0SVaclav Hapla 
315*5c820da0SVaclav Hapla   PetscFunctionBeginUser;
316*5c820da0SVaclav Hapla   if (ctx->num_labels < 0) PetscFunctionReturn(0);
317*5c820da0SVaclav Hapla   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
318*5c820da0SVaclav Hapla   PetscCallMPI(MPI_Comm_rank(comm, &rank));
319*5c820da0SVaclav Hapla   PetscCall(DMGetNumLabels(dm, &actualNum));
320*5c820da0SVaclav Hapla   if (ctx->num_labels != actualNum) {
321*5c820da0SVaclav Hapla     fail = PETSC_TRUE;
322*5c820da0SVaclav Hapla     if (ctx->verbose) {
323*5c820da0SVaclav Hapla       PetscInt    i;
324*5c820da0SVaclav Hapla 
325*5c820da0SVaclav Hapla       PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Asserted number of labels: %D, actual: %D\n", rank, ctx->num_labels, actualNum));
326*5c820da0SVaclav Hapla       for (i=0; i<actualNum; i++) {
327*5c820da0SVaclav Hapla         DMLabel     label;
328*5c820da0SVaclav Hapla         const char *name;
329*5c820da0SVaclav Hapla 
330*5c820da0SVaclav Hapla         PetscCall(DMGetLabelByNum(dm, i, &label));
331*5c820da0SVaclav Hapla         PetscCall(PetscObjectGetName((PetscObject)label, &name));
332*5c820da0SVaclav Hapla         PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label %D \"%s\"\n", rank, i, name));
333*5c820da0SVaclav Hapla       }
334*5c820da0SVaclav Hapla       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
335*5c820da0SVaclav Hapla     }
336*5c820da0SVaclav Hapla   }
337*5c820da0SVaclav Hapla   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
338*5c820da0SVaclav Hapla   PetscCheck(!fail, comm, PETSC_ERR_PLIB, "Wrong number of labels%s", ctx->verbose ? " - see details above" : "");
339*5c820da0SVaclav Hapla   PetscFunctionReturn(0);
340*5c820da0SVaclav Hapla }
341*5c820da0SVaclav Hapla 
342*5c820da0SVaclav Hapla static inline PetscErrorCode IncrementNumLabels(AppCtx *ctx)
343*5c820da0SVaclav Hapla {
344*5c820da0SVaclav Hapla   PetscFunctionBeginUser;
345*5c820da0SVaclav Hapla   if (ctx->num_labels >= 0) ctx->num_labels++;
3461d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
3471d5c9e23SVaclav Hapla }
3481d5c9e23SVaclav Hapla 
3491d5c9e23SVaclav Hapla int main(int argc, char **argv)
3501d5c9e23SVaclav Hapla {
351*5c820da0SVaclav Hapla   const char     BOUNDARY_NAME[]          = "Boundary";
352*5c820da0SVaclav Hapla   const char     BOUNDARY_VERTICES_NAME[] = "BoundaryVertices";
353*5c820da0SVaclav Hapla   const PetscInt BOUNDARY_VALUE           = 12345;
354*5c820da0SVaclav Hapla   const PetscInt BOUNDARY_VERTICES_VALUE  = 6789;
3551d5c9e23SVaclav Hapla   DM             dm, dmnew;
3561d5c9e23SVaclav Hapla   AppCtx         user;
3571d5c9e23SVaclav Hapla   Vec            allCoords = NULL;
3581d5c9e23SVaclav Hapla 
3599566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3609566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3619566063dSJacob Faibussowitsch   PetscCall(CreateMesh(&user, &dm));
362*5c820da0SVaclav Hapla   PetscCall(MarkBoundary(dm, BOUNDARY_NAME, BOUNDARY_VALUE, PETSC_FALSE, NULL));
363*5c820da0SVaclav Hapla   PetscCall(IncrementNumLabels(&user));
3641d5c9e23SVaclav Hapla   if (user.compare_boundary) {
365*5c820da0SVaclav Hapla     DMLabel label;
366*5c820da0SVaclav Hapla 
367*5c820da0SVaclav Hapla     PetscCall(MarkBoundary(dm, BOUNDARY_VERTICES_NAME, BOUNDARY_VERTICES_VALUE, PETSC_TRUE, &label));
368*5c820da0SVaclav Hapla     PetscCall(IncrementNumLabels(&user));
369*5c820da0SVaclav Hapla     PetscCall(DMLabelGetCoordinateRepresentation(dm, label, BOUNDARY_VERTICES_VALUE, &allCoords));
370*5c820da0SVaclav Hapla     PetscCall(DMLabelDestroy(&label));
3711d5c9e23SVaclav Hapla   }
3729566063dSJacob Faibussowitsch   PetscCall(SaveMesh(&user, dm));
373*5c820da0SVaclav Hapla 
3749566063dSJacob Faibussowitsch   PetscCall(LoadMesh(&user, &dmnew));
375*5c820da0SVaclav Hapla   PetscCall(IncrementNumLabels(&user)); /* account for depth label */
376*5c820da0SVaclav Hapla   PetscCall(IncrementNumLabels(&user)); /* account for celltype label */
377*5c820da0SVaclav Hapla   PetscCall(CheckNumLabels(dm, &user));
3789566063dSJacob Faibussowitsch   PetscCall(CompareMeshes(&user, dm, dmnew));
3791d5c9e23SVaclav Hapla   if (user.compare_boundary) {
380*5c820da0SVaclav Hapla     DMLabel label;
381*5c820da0SVaclav Hapla 
382*5c820da0SVaclav Hapla     PetscCall(DMGetLabel(dmnew, BOUNDARY_VERTICES_NAME, &label));
383*5c820da0SVaclav Hapla     PetscCall(DMLabelCompareWithCoordinateRepresentation(dmnew, label, BOUNDARY_VERTICES_VALUE, allCoords, user.verbose));
3841d5c9e23SVaclav Hapla   }
3859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmnew));
3879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&allCoords));
3889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
389b122ec5aSJacob Faibussowitsch   return 0;
3901d5c9e23SVaclav Hapla }
3911d5c9e23SVaclav Hapla 
3921d5c9e23SVaclav Hapla //TODO we can -compare once the new parallel topology format is in place
3931d5c9e23SVaclav Hapla /*TEST
3941d5c9e23SVaclav Hapla   build:
3951d5c9e23SVaclav Hapla     requires: hdf5
3961d5c9e23SVaclav Hapla 
3971d5c9e23SVaclav Hapla   # load old format, save in new format, reload, distribute
3981d5c9e23SVaclav Hapla   testset:
3991d5c9e23SVaclav Hapla     suffix: 1
4001d5c9e23SVaclav Hapla     requires: !complex datafilespath
4011d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4021d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
403e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
4041d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
405e600fa54SMatthew G. Knepley     args: -use_low_level_functions {{0 1}} -compare_boundary
406*5c820da0SVaclav Hapla     args: -num_labels 1
4071d5c9e23SVaclav Hapla     args: -outfile ex56_1.h5
4081d5c9e23SVaclav Hapla     nsize: {{1 3}}
4091d5c9e23SVaclav Hapla     test:
4101d5c9e23SVaclav Hapla       suffix: a
4111d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4121d5c9e23SVaclav Hapla     test:
4131d5c9e23SVaclav Hapla       suffix: b
4141d5c9e23SVaclav Hapla       TODO: broken
4151d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4161d5c9e23SVaclav Hapla     test:
4171d5c9e23SVaclav Hapla       suffix: c
4181d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4191d5c9e23SVaclav Hapla     test:
4201d5c9e23SVaclav Hapla       suffix: d
421*5c820da0SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
4221d5c9e23SVaclav Hapla     test:
4231d5c9e23SVaclav Hapla       suffix: e
4241d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4251d5c9e23SVaclav Hapla     test:
4261d5c9e23SVaclav Hapla       suffix: f
4271d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4281d5c9e23SVaclav Hapla 
4291d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
4301d5c9e23SVaclav Hapla   testset:
4311d5c9e23SVaclav Hapla     suffix: 2
4321d5c9e23SVaclav Hapla     requires: !complex datafilespath
4331d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4341d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
435e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
4361d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
437e600fa54SMatthew G. Knepley     args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary
438*5c820da0SVaclav Hapla     args: -num_labels 1
4391d5c9e23SVaclav Hapla     args: -outfile ex56_2.h5
4401d5c9e23SVaclav Hapla     nsize: 3
4411d5c9e23SVaclav Hapla     test:
4421d5c9e23SVaclav Hapla       suffix: a
4431d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4441d5c9e23SVaclav Hapla     test:
4451d5c9e23SVaclav Hapla       suffix: b
4461d5c9e23SVaclav Hapla       TODO: broken
4471d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4481d5c9e23SVaclav Hapla     test:
4491d5c9e23SVaclav Hapla       suffix: c
4501d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4511d5c9e23SVaclav Hapla     test:
4521d5c9e23SVaclav Hapla       suffix: d
453*5c820da0SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
4541d5c9e23SVaclav Hapla     test:
4551d5c9e23SVaclav Hapla       suffix: e
4561d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4571d5c9e23SVaclav Hapla     test:
4581d5c9e23SVaclav Hapla       suffix: f
4591d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4601d5c9e23SVaclav Hapla 
4611d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
4621d5c9e23SVaclav Hapla   testset:
4631d5c9e23SVaclav Hapla     suffix: 3
4641d5c9e23SVaclav Hapla     requires: !complex datafilespath
4651d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4661d5c9e23SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 2.0.0
467e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate -load_dm_distribute 0
4681d5c9e23SVaclav Hapla     args: -use_low_level_functions -compare_pre_post
469*5c820da0SVaclav Hapla     args: -num_labels 1
4701d5c9e23SVaclav Hapla     args: -outfile ex56_3.h5
4711d5c9e23SVaclav Hapla     nsize: 3
4721d5c9e23SVaclav Hapla     test:
4731d5c9e23SVaclav Hapla       suffix: a
4741d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4751d5c9e23SVaclav Hapla     test:
4761d5c9e23SVaclav Hapla       suffix: b
4771d5c9e23SVaclav Hapla       TODO: broken
4781d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4791d5c9e23SVaclav Hapla     test:
4801d5c9e23SVaclav Hapla       suffix: c
4811d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4821d5c9e23SVaclav Hapla     test:
4831d5c9e23SVaclav Hapla       suffix: d
484*5c820da0SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
4851d5c9e23SVaclav Hapla     test:
4861d5c9e23SVaclav Hapla       suffix: e
4871d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4881d5c9e23SVaclav Hapla     test:
4891d5c9e23SVaclav Hapla       suffix: f
4901d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4911d5c9e23SVaclav Hapla TEST*/
492