xref: /petsc/src/dm/impls/plex/tests/ex56.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
381d5c9e23SVaclav Hapla   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-verbose", "Verbose printing", EX, options->verbose, &options->verbose, NULL));
471d5c9e23SVaclav Hapla   ierr = PetscOptionsEnd();CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(options->comm, &dm));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetName((PetscObject)dm, &options->meshname));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), options->outfile, FILE_MODE_WRITE, &v));
711d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
72*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTopologyView(dm, v));
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCoordinatesView(dm, v));
74*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexLabelsView(dm, v));
751d5c9e23SVaclav Hapla   } else {
76*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMView(dm, v));
771d5c9e23SVaclav Hapla   }
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(options->comm, &dm));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dm, options->meshname));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexTopologyLoad(dm, v, &sfXC));
941d5c9e23SVaclav Hapla   if (mode == PRE_DIST) {
95*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCoordinatesLoad(dm, v, sfXC));
96*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexLabelsLoad(dm, v, sfXC));
971d5c9e23SVaclav Hapla   }
981d5c9e23SVaclav Hapla   if (explicitDistribute) {
991d5c9e23SVaclav Hapla     DM      dmdist;
1001d5c9e23SVaclav Hapla     PetscSF sfXB = sfXC, sfBC;
1011d5c9e23SVaclav Hapla 
102*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistribute(dm, 0, &sfBC, &dmdist));
1031d5c9e23SVaclav Hapla     if (dmdist) {
1041d5c9e23SVaclav Hapla       const char *name;
1051d5c9e23SVaclav Hapla 
106*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectGetName((PetscObject) dm, &name));
107*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) dmdist, name));
108*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFCompose(sfXB, sfBC, &sfXC));
109*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFDestroy(&sfXB));
110*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFDestroy(&sfBC));
111*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&dm));
1121d5c9e23SVaclav Hapla       dm   = dmdist;
1131d5c9e23SVaclav Hapla     }
1141d5c9e23SVaclav Hapla   }
1151d5c9e23SVaclav Hapla   if (mode == POST_DIST) {
116*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCoordinatesLoad(dm, v, sfXC));
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexLabelsLoad(dm, v, sfXC));
1181d5c9e23SVaclav Hapla   }
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch       CHKERRQ(LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0));
136*5f80ce2aSJacob Faibussowitsch       CHKERRQ(LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm));
137*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCompareLabels(dm0, dm, NULL, NULL));
138*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&dm0));
1391d5c9e23SVaclav Hapla     } else {
140*5f80ce2aSJacob Faibussowitsch       CHKERRQ(LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm));
1411d5c9e23SVaclav Hapla     }
1421d5c9e23SVaclav Hapla   } else {
143*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreate(options->comm, &dm));
144*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetType(dm, DMPLEX));
145*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) dm, options->meshname));
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLoad(dm, v));
1471d5c9e23SVaclav Hapla   }
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&v));
1491d5c9e23SVaclav Hapla 
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOptionsPrefix(dm, "load_"));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexEqual(dm0, dm1, &flg));
1642c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
165*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(options->comm,"DMs equal\n"));
1661d5c9e23SVaclav Hapla   }
1671d5c9e23SVaclav Hapla   if (options->compare_labels) {
168*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCompareLabels(dm0, dm1, NULL, NULL));
169*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PetscObjectComm((PetscObject)dm), LABEL_NAME, &l));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMarkBoundaryFaces(dm, value, l));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLabelComplete(dm, l));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetStratumIS(l, value, &points));
1861d5c9e23SVaclav Hapla 
187*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(points, &n));
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPointDepth(dm, p, &d));
1941d5c9e23SVaclav Hapla     if (d != 0) {
195*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMLabelClearValue(l, p, value));
1961d5c9e23SVaclav Hapla     }
1971d5c9e23SVaclav Hapla   }
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(points, &idx));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm));
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocalSetUp(dm));
2131d5c9e23SVaclav Hapla   if (vertices) {
214*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords));
2151d5c9e23SVaclav Hapla   } else {
216*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF, 0, &coords));
2171d5c9e23SVaclav Hapla   }
2181d5c9e23SVaclav Hapla   {
2191d5c9e23SVaclav Hapla     PetscInt  n;
2201d5c9e23SVaclav Hapla     Vec       mpivec;
2211d5c9e23SVaclav Hapla 
222*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(coords, &n));
223*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateMPI(comm, n, PETSC_DECIDE, &mpivec));
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(coords, mpivec));
225*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&coords));
2261d5c9e23SVaclav Hapla     coords = mpivec;
2271d5c9e23SVaclav Hapla   }
2281d5c9e23SVaclav Hapla 
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreateToAll(coords, &sc, &allCoords_));
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD));
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&sc));
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MarkBoundaryVertices(dm, LABEL_VALUE, &label));
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetStratumIS(label, LABEL_VALUE, &vertices));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VertexCoordinatesToAll(dm, vertices, allCoords));
248*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAddLabel(dm, label));
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&label));
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm));
266*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, LABEL_NAME, &label));
2682c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!label,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label \"%s\" was not loaded", LABEL_NAME);
2691d5c9e23SVaclav Hapla   {
2701d5c9e23SVaclav Hapla     PetscInt pStart, pEnd;
2711d5c9e23SVaclav Hapla 
272*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetChart(dm, &pStart, &pEnd));
273*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelCreateIndex(label, pStart, pEnd));
2741d5c9e23SVaclav Hapla   }
275*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS));
276*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(pointsIS, &points));
277*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(pointsIS, &n));
278*5f80ce2aSJacob Faibussowitsch   if (user->verbose) CHKERRQ(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*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasPoint(label, p, &has));
2861d5c9e23SVaclav Hapla     if (!has) {
287*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetValue(label, p, &v));
2921d5c9e23SVaclav Hapla     if (v != LABEL_VALUE) {
293*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %D has bad value %D", p, v));
2941d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
2951d5c9e23SVaclav Hapla       continue;
2961d5c9e23SVaclav Hapla     }
297*5f80ce2aSJacob Faibussowitsch     if (user->verbose) CHKERRQ(PetscSynchronizedPrintf(comm, "[%d] OK point %D\n", rank, p));
2981d5c9e23SVaclav Hapla   }
299*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT));
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDERR));
301*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(pointsIS, &points));
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&pointsIS));
303*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
3042c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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   PetscErrorCode ierr;
3141d5c9e23SVaclav Hapla 
3151d5c9e23SVaclav Hapla   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
316*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
317*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(&user, &dm));
3181d5c9e23SVaclav Hapla   if (user.compare_boundary) {
319*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundaryLabel_GetCoordinateRepresentation(dm, &allCoords));
3201d5c9e23SVaclav Hapla   }
321*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SaveMesh(&user, dm));
322*5f80ce2aSJacob Faibussowitsch   CHKERRQ(LoadMesh(&user, &dmnew));
323*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CompareMeshes(&user, dm, dmnew));
3241d5c9e23SVaclav Hapla   if (user.compare_boundary) {
325*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetBoundaryLabel_CompareWithCoordinateRepresentation(&user, dmnew, allCoords));
3261d5c9e23SVaclav Hapla   }
327*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmnew));
329*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&allCoords));
3301d5c9e23SVaclav Hapla   ierr = PetscFinalize();
3311d5c9e23SVaclav Hapla   return ierr;
3321d5c9e23SVaclav Hapla }
3331d5c9e23SVaclav Hapla 
3341d5c9e23SVaclav Hapla //TODO we can -compare once the new parallel topology format is in place
3351d5c9e23SVaclav Hapla /*TEST
3361d5c9e23SVaclav Hapla   build:
3371d5c9e23SVaclav Hapla     requires: hdf5
3381d5c9e23SVaclav Hapla 
3391d5c9e23SVaclav Hapla   # load old format, save in new format, reload, distribute
3401d5c9e23SVaclav Hapla   testset:
3411d5c9e23SVaclav Hapla     suffix: 1
3421d5c9e23SVaclav Hapla     requires: !complex datafilespath
3431d5c9e23SVaclav Hapla     args: -dm_plex_name plex
3441d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
345e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
3461d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
347e600fa54SMatthew G. Knepley     args: -use_low_level_functions {{0 1}} -compare_boundary
3481d5c9e23SVaclav Hapla     args: -outfile ex56_1.h5
3491d5c9e23SVaclav Hapla     nsize: {{1 3}}
3501d5c9e23SVaclav Hapla     test:
3511d5c9e23SVaclav Hapla       suffix: a
3521d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
3531d5c9e23SVaclav Hapla     test:
3541d5c9e23SVaclav Hapla       suffix: b
3551d5c9e23SVaclav Hapla       TODO: broken
3561d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
3571d5c9e23SVaclav Hapla     test:
3581d5c9e23SVaclav Hapla       suffix: c
3591d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
3601d5c9e23SVaclav Hapla     test:
3611d5c9e23SVaclav Hapla       suffix: d
3621d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
3631d5c9e23SVaclav Hapla     test:
3641d5c9e23SVaclav Hapla       suffix: e
3651d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
3661d5c9e23SVaclav Hapla     test:
3671d5c9e23SVaclav Hapla       suffix: f
3681d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
3691d5c9e23SVaclav Hapla 
3701d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
3711d5c9e23SVaclav Hapla   testset:
3721d5c9e23SVaclav Hapla     suffix: 2
3731d5c9e23SVaclav Hapla     requires: !complex datafilespath
3741d5c9e23SVaclav Hapla     args: -dm_plex_name plex
3751d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
376e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate
3771d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
378e600fa54SMatthew G. Knepley     args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary
3791d5c9e23SVaclav Hapla     args: -outfile ex56_2.h5
3801d5c9e23SVaclav Hapla     nsize: 3
3811d5c9e23SVaclav Hapla     test:
3821d5c9e23SVaclav Hapla       suffix: a
3831d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
3841d5c9e23SVaclav Hapla     test:
3851d5c9e23SVaclav Hapla       suffix: b
3861d5c9e23SVaclav Hapla       TODO: broken
3871d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
3881d5c9e23SVaclav Hapla     test:
3891d5c9e23SVaclav Hapla       suffix: c
3901d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
3911d5c9e23SVaclav Hapla     test:
3921d5c9e23SVaclav Hapla       suffix: d
3931d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
3941d5c9e23SVaclav Hapla     test:
3951d5c9e23SVaclav Hapla       suffix: e
3961d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
3971d5c9e23SVaclav Hapla     test:
3981d5c9e23SVaclav Hapla       suffix: f
3991d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4001d5c9e23SVaclav Hapla 
4011d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
4021d5c9e23SVaclav Hapla   testset:
4031d5c9e23SVaclav Hapla     suffix: 3
4041d5c9e23SVaclav Hapla     requires: !complex datafilespath
4051d5c9e23SVaclav Hapla     args: -dm_plex_name plex
4061d5c9e23SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 2.0.0
407e600fa54SMatthew G. Knepley     args: -dm_plex_interpolate -load_dm_distribute 0
4081d5c9e23SVaclav Hapla     args: -use_low_level_functions -compare_pre_post
4091d5c9e23SVaclav Hapla     args: -outfile ex56_3.h5
4101d5c9e23SVaclav Hapla     nsize: 3
4111d5c9e23SVaclav Hapla     test:
4121d5c9e23SVaclav Hapla       suffix: a
4131d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
4141d5c9e23SVaclav Hapla     test:
4151d5c9e23SVaclav Hapla       suffix: b
4161d5c9e23SVaclav Hapla       TODO: broken
4171d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
4181d5c9e23SVaclav Hapla     test:
4191d5c9e23SVaclav Hapla       suffix: c
4201d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
4211d5c9e23SVaclav Hapla     test:
4221d5c9e23SVaclav Hapla       suffix: d
4231d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
4241d5c9e23SVaclav Hapla     test:
4251d5c9e23SVaclav Hapla       suffix: e
4261d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
4271d5c9e23SVaclav Hapla     test:
4281d5c9e23SVaclav Hapla       suffix: f
4291d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
4301d5c9e23SVaclav Hapla TEST*/
431