xref: /petsc/src/dm/impls/plex/tests/ex56.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 static const char LABEL_NAME[] = "BoundaryVertices";
81d5c9e23SVaclav Hapla static const PetscInt LABEL_VALUE = 12345;
91d5c9e23SVaclav Hapla typedef struct {
101d5c9e23SVaclav Hapla   MPI_Comm    comm;
111d5c9e23SVaclav Hapla   const char *meshname;                     /* Mesh name */
121d5c9e23SVaclav Hapla   PetscBool   compare;                      /* Compare the meshes using DMPlexEqual() and DMCompareLabels() */
131d5c9e23SVaclav Hapla   PetscBool   compare_labels;               /* Compare labels in the meshes using DMCompareLabels() */
141d5c9e23SVaclav Hapla   PetscBool   compare_boundary;             /* Check label I/O via boundary vertex coordinates */
151d5c9e23SVaclav Hapla   PetscBool   compare_pre_post;             /* Compare labels loaded before distribution with those loaded after distribution */
161d5c9e23SVaclav Hapla   char        outfile[PETSC_MAX_PATH_LEN];  /* Output file */
171d5c9e23SVaclav Hapla   PetscBool   use_low_level_functions;      /* Use low level functions for viewing and loading */
181d5c9e23SVaclav Hapla   //TODO This is meant as temporary option; can be removed once we have full parallel loading in place
191d5c9e23SVaclav Hapla   PetscBool   distribute_after_topo_load;   /* Distribute topology right after DMPlexTopologyLoad(), if use_low_level_functions=true */
201d5c9e23SVaclav Hapla   PetscBool   verbose;
211d5c9e23SVaclav Hapla } AppCtx;
221d5c9e23SVaclav Hapla 
231d5c9e23SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
241d5c9e23SVaclav Hapla {
251d5c9e23SVaclav Hapla   PetscErrorCode ierr;
261d5c9e23SVaclav Hapla 
271d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
281d5c9e23SVaclav Hapla   options->comm                       = comm;
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;
361d5c9e23SVaclav Hapla   options->verbose                    = PETSC_FALSE;
371d5c9e23SVaclav Hapla 
38*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");PetscCall(ierr);
39*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL));
40*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL));
41*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL));
42*9566063dSJacob 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));
43*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL));
44*9566063dSJacob 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));
45*9566063dSJacob 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));
46*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-verbose", "Verbose printing", EX, options->verbose, &options->verbose, NULL));
47*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
481d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
491d5c9e23SVaclav Hapla };
501d5c9e23SVaclav Hapla 
511d5c9e23SVaclav Hapla static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm)
521d5c9e23SVaclav Hapla {
531d5c9e23SVaclav Hapla   DM             dm;
541d5c9e23SVaclav Hapla 
551d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
56*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(options->comm, &dm));
57*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
58*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
59*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)dm, &options->meshname));
60*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
611d5c9e23SVaclav Hapla   *newdm = dm;
621d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
631d5c9e23SVaclav Hapla }
641d5c9e23SVaclav Hapla 
651d5c9e23SVaclav Hapla static PetscErrorCode SaveMesh(AppCtx *options, DM dm)
661d5c9e23SVaclav Hapla {
671d5c9e23SVaclav Hapla   PetscViewer    v;
681d5c9e23SVaclav Hapla 
691d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
70*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), options->outfile, FILE_MODE_WRITE, &v));
711d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
72*9566063dSJacob Faibussowitsch     PetscCall(DMPlexTopologyView(dm, v));
73*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesView(dm, v));
74*9566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsView(dm, v));
751d5c9e23SVaclav Hapla   } else {
76*9566063dSJacob Faibussowitsch     PetscCall(DMView(dm, v));
771d5c9e23SVaclav Hapla   }
78*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
791d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
801d5c9e23SVaclav Hapla }
811d5c9e23SVaclav Hapla 
821d5c9e23SVaclav Hapla typedef enum {NONE=0, PRE_DIST=1, POST_DIST=2} AuxObjLoadMode;
831d5c9e23SVaclav Hapla 
841d5c9e23SVaclav Hapla static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm)
851d5c9e23SVaclav Hapla {
861d5c9e23SVaclav Hapla   DM              dm;
871d5c9e23SVaclav Hapla   PetscSF         sfXC;
881d5c9e23SVaclav Hapla 
891d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
90*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(options->comm, &dm));
91*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
92*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname));
93*9566063dSJacob Faibussowitsch   PetscCall(DMPlexTopologyLoad(dm, v, &sfXC));
941d5c9e23SVaclav Hapla   if (mode == PRE_DIST) {
95*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC));
96*9566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsLoad(dm, v, sfXC));
971d5c9e23SVaclav Hapla   }
981d5c9e23SVaclav Hapla   if (explicitDistribute) {
991d5c9e23SVaclav Hapla     DM      dmdist;
1001d5c9e23SVaclav Hapla     PetscSF sfXB = sfXC, sfBC;
1011d5c9e23SVaclav Hapla 
102*9566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(dm, 0, &sfBC, &dmdist));
1031d5c9e23SVaclav Hapla     if (dmdist) {
1041d5c9e23SVaclav Hapla       const char *name;
1051d5c9e23SVaclav Hapla 
106*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) dm, &name));
107*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) dmdist, name));
108*9566063dSJacob Faibussowitsch       PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC));
109*9566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&sfXB));
110*9566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&sfBC));
111*9566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
1121d5c9e23SVaclav Hapla       dm   = dmdist;
1131d5c9e23SVaclav Hapla     }
1141d5c9e23SVaclav Hapla   }
1151d5c9e23SVaclav Hapla   if (mode == POST_DIST) {
116*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC));
117*9566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelsLoad(dm, v, sfXC));
1181d5c9e23SVaclav Hapla   }
119*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfXC));
1201d5c9e23SVaclav Hapla   *newdm = dm;
1211d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1221d5c9e23SVaclav Hapla }
1231d5c9e23SVaclav Hapla 
1241d5c9e23SVaclav Hapla static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew)
1251d5c9e23SVaclav Hapla {
1261d5c9e23SVaclav Hapla   DM             dm;
1271d5c9e23SVaclav Hapla   PetscViewer    v;
1281d5c9e23SVaclav Hapla 
1291d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
130*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(options->comm, options->outfile, FILE_MODE_READ, &v));
1311d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
1321d5c9e23SVaclav Hapla     if (options->compare_pre_post) {
1331d5c9e23SVaclav Hapla       DM dm0;
1341d5c9e23SVaclav Hapla 
135*9566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0));
136*9566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm));
137*9566063dSJacob Faibussowitsch       PetscCall(DMCompareLabels(dm0, dm, NULL, NULL));
138*9566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm0));
1391d5c9e23SVaclav Hapla     } else {
140*9566063dSJacob Faibussowitsch       PetscCall(LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm));
1411d5c9e23SVaclav Hapla     }
1421d5c9e23SVaclav Hapla   } else {
143*9566063dSJacob Faibussowitsch     PetscCall(DMCreate(options->comm, &dm));
144*9566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm, DMPLEX));
145*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname));
146*9566063dSJacob Faibussowitsch     PetscCall(DMLoad(dm, v));
1471d5c9e23SVaclav Hapla   }
148*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&v));
1491d5c9e23SVaclav Hapla 
150*9566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(dm, "load_"));
151*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
152*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1531d5c9e23SVaclav Hapla   *dmnew = dm;
1541d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1551d5c9e23SVaclav Hapla }
1561d5c9e23SVaclav Hapla 
1571d5c9e23SVaclav Hapla static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1)
1581d5c9e23SVaclav Hapla {
1591d5c9e23SVaclav Hapla   PetscBool       flg;
1601d5c9e23SVaclav Hapla 
1611d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
1621d5c9e23SVaclav Hapla   if (options->compare) {
163*9566063dSJacob Faibussowitsch     PetscCall(DMPlexEqual(dm0, dm1, &flg));
16428b400f6SJacob Faibussowitsch     PetscCheck(flg,options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
165*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(options->comm,"DMs equal\n"));
1661d5c9e23SVaclav Hapla   }
1671d5c9e23SVaclav Hapla   if (options->compare_labels) {
168*9566063dSJacob Faibussowitsch     PetscCall(DMCompareLabels(dm0, dm1, NULL, NULL));
169*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(options->comm,"DMLabels equal\n"));
1701d5c9e23SVaclav Hapla   }
1711d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
1721d5c9e23SVaclav Hapla }
1731d5c9e23SVaclav Hapla 
1741d5c9e23SVaclav Hapla static PetscErrorCode MarkBoundaryVertices(DM dm, PetscInt value, DMLabel *label)
1751d5c9e23SVaclav Hapla {
1761d5c9e23SVaclav Hapla   DMLabel         l;
1771d5c9e23SVaclav Hapla   IS              points;
1781d5c9e23SVaclav Hapla   const PetscInt *idx;
1791d5c9e23SVaclav Hapla   PetscInt        i, n;
1801d5c9e23SVaclav Hapla 
1811d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
182*9566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)dm), LABEL_NAME, &l));
183*9566063dSJacob Faibussowitsch   PetscCall(DMPlexMarkBoundaryFaces(dm, value, l));
184*9566063dSJacob Faibussowitsch   PetscCall(DMPlexLabelComplete(dm, l));
185*9566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(l, value, &points));
1861d5c9e23SVaclav Hapla 
187*9566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(points, &n));
188*9566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(points, &idx));
1891d5c9e23SVaclav Hapla   for (i=0; i<n; i++) {
1901d5c9e23SVaclav Hapla     const PetscInt p = idx[i];
1911d5c9e23SVaclav Hapla     PetscInt       d;
1921d5c9e23SVaclav Hapla 
193*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPointDepth(dm, p, &d));
1941d5c9e23SVaclav Hapla     if (d != 0) {
195*9566063dSJacob Faibussowitsch       PetscCall(DMLabelClearValue(l, p, value));
1961d5c9e23SVaclav Hapla     }
1971d5c9e23SVaclav Hapla   }
198*9566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(points, &idx));
199*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&points));
2001d5c9e23SVaclav Hapla   *label = l;
2011d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2021d5c9e23SVaclav Hapla }
2031d5c9e23SVaclav Hapla 
2041d5c9e23SVaclav Hapla static PetscErrorCode VertexCoordinatesToAll(DM dm, IS vertices, Vec *allCoords)
2051d5c9e23SVaclav Hapla {
2061d5c9e23SVaclav Hapla   Vec             coords, allCoords_;
2071d5c9e23SVaclav Hapla   VecScatter      sc;
2081d5c9e23SVaclav Hapla   MPI_Comm        comm;
2091d5c9e23SVaclav Hapla 
2101d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
211*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
212*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
2131d5c9e23SVaclav Hapla   if (vertices) {
214*9566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords));
2151d5c9e23SVaclav Hapla   } else {
216*9566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &coords));
2171d5c9e23SVaclav Hapla   }
2181d5c9e23SVaclav Hapla   {
2191d5c9e23SVaclav Hapla     PetscInt  n;
2201d5c9e23SVaclav Hapla     Vec       mpivec;
2211d5c9e23SVaclav Hapla 
222*9566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coords, &n));
223*9566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(comm, n, PETSC_DECIDE, &mpivec));
224*9566063dSJacob Faibussowitsch     PetscCall(VecCopy(coords, mpivec));
225*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&coords));
2261d5c9e23SVaclav Hapla     coords = mpivec;
2271d5c9e23SVaclav Hapla   }
2281d5c9e23SVaclav Hapla 
229*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToAll(coords, &sc, &allCoords_));
230*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD));
231*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD));
232*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sc));
233*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coords));
2341d5c9e23SVaclav Hapla   *allCoords = allCoords_;
2351d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2361d5c9e23SVaclav Hapla }
2371d5c9e23SVaclav Hapla 
2381d5c9e23SVaclav Hapla /* Compute boundary label, remember boundary vertices using coordinates, save and load label, check it is defined on the original boundary vertices */
2391d5c9e23SVaclav Hapla static PetscErrorCode DMAddBoundaryLabel_GetCoordinateRepresentation(DM dm, Vec *allCoords)
2401d5c9e23SVaclav Hapla {
2411d5c9e23SVaclav Hapla   DMLabel         label;
2421d5c9e23SVaclav Hapla   IS              vertices;
2431d5c9e23SVaclav Hapla 
2441d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
245*9566063dSJacob Faibussowitsch   PetscCall(MarkBoundaryVertices(dm, LABEL_VALUE, &label));
246*9566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(label, LABEL_VALUE, &vertices));
247*9566063dSJacob Faibussowitsch   PetscCall(VertexCoordinatesToAll(dm, vertices, allCoords));
248*9566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(dm, label));
249*9566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
250*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&vertices));
2511d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
2521d5c9e23SVaclav Hapla }
2531d5c9e23SVaclav Hapla 
2541d5c9e23SVaclav Hapla static PetscErrorCode DMGetBoundaryLabel_CompareWithCoordinateRepresentation(AppCtx *user, DM dm, Vec allCoords)
2551d5c9e23SVaclav Hapla {
2561d5c9e23SVaclav Hapla   DMLabel         label;
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*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
266*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
267*9566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, LABEL_NAME, &label));
26828b400f6SJacob Faibussowitsch   PetscCheck(label,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label \"%s\" was not loaded", LABEL_NAME);
2691d5c9e23SVaclav Hapla   {
2701d5c9e23SVaclav Hapla     PetscInt pStart, pEnd;
2711d5c9e23SVaclav Hapla 
272*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
273*9566063dSJacob Faibussowitsch     PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
2741d5c9e23SVaclav Hapla   }
275*9566063dSJacob Faibussowitsch   PetscCall(DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS));
276*9566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointsIS, &points));
277*9566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointsIS, &n));
278*9566063dSJacob Faibussowitsch   if (user->verbose) 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;
285*9566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(label, p, &has));
2861d5c9e23SVaclav Hapla     if (!has) {
287*9566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label does not have point %D\n", rank, p));
2881d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2891d5c9e23SVaclav Hapla       continue;
2901d5c9e23SVaclav Hapla     }
291*9566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(label, p, &v));
2921d5c9e23SVaclav Hapla     if (v != LABEL_VALUE) {
293*9566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %D has bad value %D", p, v));
2941d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2951d5c9e23SVaclav Hapla       continue;
2961d5c9e23SVaclav Hapla     }
297*9566063dSJacob Faibussowitsch     if (user->verbose) PetscCall(PetscSynchronizedPrintf(comm, "[%d] OK point %D\n", rank, p));
2981d5c9e23SVaclav Hapla   }
299*9566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
300*9566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
301*9566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointsIS, &points));
302*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointsIS));
303*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
30428b400f6SJacob Faibussowitsch   PetscCheck(!fail,comm, PETSC_ERR_PLIB, "Label \"%s\" was not loaded correctly - see details above", LABEL_NAME);
3051d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
3061d5c9e23SVaclav Hapla }
3071d5c9e23SVaclav Hapla 
3081d5c9e23SVaclav Hapla int main(int argc, char **argv)
3091d5c9e23SVaclav Hapla {
3101d5c9e23SVaclav Hapla   DM             dm, dmnew;
3111d5c9e23SVaclav Hapla   AppCtx         user;
3121d5c9e23SVaclav Hapla   Vec            allCoords = NULL;
3131d5c9e23SVaclav Hapla 
314*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
315*9566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
316*9566063dSJacob Faibussowitsch   PetscCall(CreateMesh(&user, &dm));
3171d5c9e23SVaclav Hapla   if (user.compare_boundary) {
318*9566063dSJacob Faibussowitsch     PetscCall(DMAddBoundaryLabel_GetCoordinateRepresentation(dm, &allCoords));
3191d5c9e23SVaclav Hapla   }
320*9566063dSJacob Faibussowitsch   PetscCall(SaveMesh(&user, dm));
321*9566063dSJacob Faibussowitsch   PetscCall(LoadMesh(&user, &dmnew));
322*9566063dSJacob Faibussowitsch   PetscCall(CompareMeshes(&user, dm, dmnew));
3231d5c9e23SVaclav Hapla   if (user.compare_boundary) {
324*9566063dSJacob Faibussowitsch     PetscCall(DMGetBoundaryLabel_CompareWithCoordinateRepresentation(&user, dmnew, allCoords));
3251d5c9e23SVaclav Hapla   }
326*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
327*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmnew));
328*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&allCoords));
329*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
330b122ec5aSJacob Faibussowitsch   return 0;
3311d5c9e23SVaclav Hapla }
3321d5c9e23SVaclav Hapla 
3331d5c9e23SVaclav Hapla //TODO we can -compare once the new parallel topology format is in place
3341d5c9e23SVaclav Hapla /*TEST
3351d5c9e23SVaclav Hapla   build:
3361d5c9e23SVaclav Hapla     requires: hdf5
3371d5c9e23SVaclav Hapla 
3381d5c9e23SVaclav Hapla   # load old format, save in new format, reload, distribute
3391d5c9e23SVaclav Hapla   testset:
3401d5c9e23SVaclav Hapla     suffix: 1
3411d5c9e23SVaclav Hapla     requires: !complex datafilespath
3421d5c9e23SVaclav Hapla     args: -dm_plex_name plex
3431d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
344e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
3451d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
346e600fa54SMatthew G. Knepley     args: -use_low_level_functions {{0 1}} -compare_boundary
3471d5c9e23SVaclav Hapla     args: -outfile ex56_1.h5
3481d5c9e23SVaclav Hapla     nsize: {{1 3}}
3491d5c9e23SVaclav Hapla     test:
3501d5c9e23SVaclav Hapla       suffix: a
3511d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
3521d5c9e23SVaclav Hapla     test:
3531d5c9e23SVaclav Hapla       suffix: b
3541d5c9e23SVaclav Hapla       TODO: broken
3551d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
3561d5c9e23SVaclav Hapla     test:
3571d5c9e23SVaclav Hapla       suffix: c
3581d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
3591d5c9e23SVaclav Hapla     test:
3601d5c9e23SVaclav Hapla       suffix: d
3611d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
3621d5c9e23SVaclav Hapla     test:
3631d5c9e23SVaclav Hapla       suffix: e
3641d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
3651d5c9e23SVaclav Hapla     test:
3661d5c9e23SVaclav Hapla       suffix: f
3671d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
3681d5c9e23SVaclav Hapla 
3691d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
3701d5c9e23SVaclav Hapla   testset:
3711d5c9e23SVaclav Hapla     suffix: 2
3721d5c9e23SVaclav Hapla     requires: !complex datafilespath
3731d5c9e23SVaclav Hapla     args: -dm_plex_name plex
3741d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
375e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
3761d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
377e600fa54SMatthew G. Knepley     args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary
3781d5c9e23SVaclav Hapla     args: -outfile ex56_2.h5
3791d5c9e23SVaclav Hapla     nsize: 3
3801d5c9e23SVaclav Hapla     test:
3811d5c9e23SVaclav Hapla       suffix: a
3821d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
3831d5c9e23SVaclav Hapla     test:
3841d5c9e23SVaclav Hapla       suffix: b
3851d5c9e23SVaclav Hapla       TODO: broken
3861d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
3871d5c9e23SVaclav Hapla     test:
3881d5c9e23SVaclav Hapla       suffix: c
3891d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
3901d5c9e23SVaclav Hapla     test:
3911d5c9e23SVaclav Hapla       suffix: d
3921d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
3931d5c9e23SVaclav Hapla     test:
3941d5c9e23SVaclav Hapla       suffix: e
3951d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
3961d5c9e23SVaclav Hapla     test:
3971d5c9e23SVaclav Hapla       suffix: f
3981d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
3991d5c9e23SVaclav Hapla 
4001d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
4011d5c9e23SVaclav Hapla   testset:
4021d5c9e23SVaclav Hapla     suffix: 3
4031d5c9e23SVaclav Hapla     requires: !complex datafilespath
4041d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4051d5c9e23SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 2.0.0
406e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate -load_dm_distribute 0
4071d5c9e23SVaclav Hapla     args: -use_low_level_functions -compare_pre_post
4081d5c9e23SVaclav Hapla     args: -outfile ex56_3.h5
4091d5c9e23SVaclav Hapla     nsize: 3
4101d5c9e23SVaclav Hapla     test:
4111d5c9e23SVaclav Hapla       suffix: a
4121d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4131d5c9e23SVaclav Hapla     test:
4141d5c9e23SVaclav Hapla       suffix: b
4151d5c9e23SVaclav Hapla       TODO: broken
4161d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4171d5c9e23SVaclav Hapla     test:
4181d5c9e23SVaclav Hapla       suffix: c
4191d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4201d5c9e23SVaclav Hapla     test:
4211d5c9e23SVaclav Hapla       suffix: d
4221d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
4231d5c9e23SVaclav Hapla     test:
4241d5c9e23SVaclav Hapla       suffix: e
4251d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4261d5c9e23SVaclav Hapla     test:
4271d5c9e23SVaclav Hapla       suffix: f
4281d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4291d5c9e23SVaclav Hapla TEST*/
430