xref: /petsc/src/dm/impls/plex/tests/ex56.c (revision 1d5c9e23f9f2b4e998e7bacb840d646f45de35ef)
1*1d5c9e23SVaclav Hapla #include <petscdmplex.h>
2*1d5c9e23SVaclav Hapla #include <petscviewerhdf5.h>
3*1d5c9e23SVaclav Hapla #include <petscsf.h>
4*1d5c9e23SVaclav Hapla 
5*1d5c9e23SVaclav Hapla static const char help[] = "Load and save the mesh to the native HDF5 format\n\n";
6*1d5c9e23SVaclav Hapla static const char EX[] = "ex56.c";
7*1d5c9e23SVaclav Hapla static const char LABEL_NAME[] = "BoundaryVertices";
8*1d5c9e23SVaclav Hapla static const PetscInt LABEL_VALUE = 12345;
9*1d5c9e23SVaclav Hapla typedef struct {
10*1d5c9e23SVaclav Hapla   MPI_Comm    comm;
11*1d5c9e23SVaclav Hapla   const char *meshname;                     /* Mesh name */
12*1d5c9e23SVaclav Hapla   PetscBool   compare;                      /* Compare the meshes using DMPlexEqual() and DMCompareLabels() */
13*1d5c9e23SVaclav Hapla   PetscBool   compare_labels;               /* Compare labels in the meshes using DMCompareLabels() */
14*1d5c9e23SVaclav Hapla   PetscBool   compare_boundary;             /* Check label I/O via boundary vertex coordinates */
15*1d5c9e23SVaclav Hapla   PetscBool   compare_pre_post;             /* Compare labels loaded before distribution with those loaded after distribution */
16*1d5c9e23SVaclav Hapla   char        outfile[PETSC_MAX_PATH_LEN];  /* Output file */
17*1d5c9e23SVaclav Hapla   PetscBool   use_low_level_functions;      /* Use low level functions for viewing and loading */
18*1d5c9e23SVaclav Hapla   //TODO This is meant as temporary option; can be removed once we have full parallel loading in place
19*1d5c9e23SVaclav Hapla   PetscBool   distribute_after_topo_load;   /* Distribute topology right after DMPlexTopologyLoad(), if use_low_level_functions=true */
20*1d5c9e23SVaclav Hapla   PetscBool   verbose;
21*1d5c9e23SVaclav Hapla } AppCtx;
22*1d5c9e23SVaclav Hapla 
23*1d5c9e23SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
24*1d5c9e23SVaclav Hapla {
25*1d5c9e23SVaclav Hapla   PetscErrorCode ierr;
26*1d5c9e23SVaclav Hapla 
27*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
28*1d5c9e23SVaclav Hapla   options->comm                       = comm;
29*1d5c9e23SVaclav Hapla   options->compare                    = PETSC_FALSE;
30*1d5c9e23SVaclav Hapla   options->compare_labels             = PETSC_FALSE;
31*1d5c9e23SVaclav Hapla   options->compare_boundary           = PETSC_FALSE;
32*1d5c9e23SVaclav Hapla   options->compare_pre_post           = PETSC_FALSE;
33*1d5c9e23SVaclav Hapla   options->outfile[0]                 = '\0';
34*1d5c9e23SVaclav Hapla   options->use_low_level_functions    = PETSC_FALSE;
35*1d5c9e23SVaclav Hapla   options->distribute_after_topo_load = PETSC_FALSE;
36*1d5c9e23SVaclav Hapla   options->verbose                    = PETSC_FALSE;
37*1d5c9e23SVaclav Hapla 
38*1d5c9e23SVaclav Hapla   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
39*1d5c9e23SVaclav Hapla   ierr = PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL);CHKERRQ(ierr);
40*1d5c9e23SVaclav Hapla   ierr = PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL);CHKERRQ(ierr);
41*1d5c9e23SVaclav Hapla   ierr = PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL);CHKERRQ(ierr);
42*1d5c9e23SVaclav Hapla   ierr = 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);CHKERRQ(ierr);
43*1d5c9e23SVaclav Hapla   ierr = PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL);CHKERRQ(ierr);
44*1d5c9e23SVaclav Hapla   ierr = 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);CHKERRQ(ierr);
45*1d5c9e23SVaclav Hapla   ierr = 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);CHKERRQ(ierr);
46*1d5c9e23SVaclav Hapla   ierr = PetscOptionsBool("-verbose", "Verbose printing", EX, options->verbose, &options->verbose, NULL);CHKERRQ(ierr);
47*1d5c9e23SVaclav Hapla   ierr = PetscOptionsEnd();CHKERRQ(ierr);
48*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
49*1d5c9e23SVaclav Hapla };
50*1d5c9e23SVaclav Hapla 
51*1d5c9e23SVaclav Hapla static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm)
52*1d5c9e23SVaclav Hapla {
53*1d5c9e23SVaclav Hapla   DM             dm;
54*1d5c9e23SVaclav Hapla   PetscErrorCode ierr;
55*1d5c9e23SVaclav Hapla 
56*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
57*1d5c9e23SVaclav Hapla   ierr = DMCreate(options->comm, &dm);CHKERRQ(ierr);
58*1d5c9e23SVaclav Hapla   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
59*1d5c9e23SVaclav Hapla   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
60*1d5c9e23SVaclav Hapla   ierr = PetscObjectGetName((PetscObject)dm, &options->meshname);CHKERRQ(ierr);
61*1d5c9e23SVaclav Hapla   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
62*1d5c9e23SVaclav Hapla   *newdm = dm;
63*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
64*1d5c9e23SVaclav Hapla }
65*1d5c9e23SVaclav Hapla 
66*1d5c9e23SVaclav Hapla static PetscErrorCode SaveMesh(AppCtx *options, DM dm)
67*1d5c9e23SVaclav Hapla {
68*1d5c9e23SVaclav Hapla   PetscViewer    v;
69*1d5c9e23SVaclav Hapla   PetscErrorCode ierr;
70*1d5c9e23SVaclav Hapla 
71*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
72*1d5c9e23SVaclav Hapla   ierr = PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), options->outfile, FILE_MODE_WRITE, &v);CHKERRQ(ierr);
73*1d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
74*1d5c9e23SVaclav Hapla     ierr = DMPlexTopologyView(dm, v);CHKERRQ(ierr);
75*1d5c9e23SVaclav Hapla     ierr = DMPlexCoordinatesView(dm, v);CHKERRQ(ierr);
76*1d5c9e23SVaclav Hapla     ierr = DMPlexLabelsView(dm, v);CHKERRQ(ierr);
77*1d5c9e23SVaclav Hapla   } else {
78*1d5c9e23SVaclav Hapla     ierr = DMView(dm, v);CHKERRQ(ierr);
79*1d5c9e23SVaclav Hapla   }
80*1d5c9e23SVaclav Hapla   ierr = PetscViewerDestroy(&v);CHKERRQ(ierr);
81*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
82*1d5c9e23SVaclav Hapla }
83*1d5c9e23SVaclav Hapla 
84*1d5c9e23SVaclav Hapla typedef enum {NONE=0, PRE_DIST=1, POST_DIST=2} AuxObjLoadMode;
85*1d5c9e23SVaclav Hapla 
86*1d5c9e23SVaclav Hapla static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm)
87*1d5c9e23SVaclav Hapla {
88*1d5c9e23SVaclav Hapla   DM              dm;
89*1d5c9e23SVaclav Hapla   PetscSF         sfXC;
90*1d5c9e23SVaclav Hapla   PetscErrorCode  ierr;
91*1d5c9e23SVaclav Hapla 
92*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
93*1d5c9e23SVaclav Hapla   ierr = DMCreate(options->comm, &dm);CHKERRQ(ierr);
94*1d5c9e23SVaclav Hapla   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
95*1d5c9e23SVaclav Hapla   ierr = PetscObjectSetName((PetscObject) dm, options->meshname);CHKERRQ(ierr);
96*1d5c9e23SVaclav Hapla   ierr = DMPlexTopologyLoad(dm, v, &sfXC);CHKERRQ(ierr);
97*1d5c9e23SVaclav Hapla   if (mode == PRE_DIST) {
98*1d5c9e23SVaclav Hapla     ierr = DMPlexCoordinatesLoad(dm, v, sfXC);CHKERRQ(ierr);
99*1d5c9e23SVaclav Hapla     ierr = DMPlexLabelsLoad(dm, v, sfXC);CHKERRQ(ierr);
100*1d5c9e23SVaclav Hapla   }
101*1d5c9e23SVaclav Hapla   if (explicitDistribute) {
102*1d5c9e23SVaclav Hapla     DM      dmdist;
103*1d5c9e23SVaclav Hapla     PetscSF sfXB = sfXC, sfBC;
104*1d5c9e23SVaclav Hapla 
105*1d5c9e23SVaclav Hapla     ierr = DMPlexDistribute(dm, 0, &sfBC, &dmdist);CHKERRQ(ierr);
106*1d5c9e23SVaclav Hapla     if (dmdist) {
107*1d5c9e23SVaclav Hapla       const char *name;
108*1d5c9e23SVaclav Hapla 
109*1d5c9e23SVaclav Hapla       ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr);
110*1d5c9e23SVaclav Hapla       ierr = PetscObjectSetName((PetscObject) dmdist, name);CHKERRQ(ierr);
111*1d5c9e23SVaclav Hapla       ierr = PetscSFCompose(sfXB, sfBC, &sfXC);CHKERRQ(ierr);
112*1d5c9e23SVaclav Hapla       ierr = PetscSFDestroy(&sfXB);CHKERRQ(ierr);
113*1d5c9e23SVaclav Hapla       ierr = PetscSFDestroy(&sfBC);CHKERRQ(ierr);
114*1d5c9e23SVaclav Hapla       ierr = DMDestroy(&dm);CHKERRQ(ierr);
115*1d5c9e23SVaclav Hapla       dm   = dmdist;
116*1d5c9e23SVaclav Hapla     }
117*1d5c9e23SVaclav Hapla   }
118*1d5c9e23SVaclav Hapla   if (mode == POST_DIST) {
119*1d5c9e23SVaclav Hapla     ierr = DMPlexCoordinatesLoad(dm, v, sfXC);CHKERRQ(ierr);
120*1d5c9e23SVaclav Hapla     ierr = DMPlexLabelsLoad(dm, v, sfXC);CHKERRQ(ierr);
121*1d5c9e23SVaclav Hapla   }
122*1d5c9e23SVaclav Hapla   ierr = PetscSFDestroy(&sfXC);CHKERRQ(ierr);
123*1d5c9e23SVaclav Hapla   *newdm = dm;
124*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
125*1d5c9e23SVaclav Hapla }
126*1d5c9e23SVaclav Hapla 
127*1d5c9e23SVaclav Hapla static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew)
128*1d5c9e23SVaclav Hapla {
129*1d5c9e23SVaclav Hapla   DM             dm;
130*1d5c9e23SVaclav Hapla   PetscViewer    v;
131*1d5c9e23SVaclav Hapla   PetscErrorCode ierr;
132*1d5c9e23SVaclav Hapla 
133*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
134*1d5c9e23SVaclav Hapla   ierr = PetscViewerHDF5Open(options->comm, options->outfile, FILE_MODE_READ, &v);CHKERRQ(ierr);
135*1d5c9e23SVaclav Hapla   if (options->use_low_level_functions) {
136*1d5c9e23SVaclav Hapla     if (options->compare_pre_post) {
137*1d5c9e23SVaclav Hapla       DM dm0;
138*1d5c9e23SVaclav Hapla 
139*1d5c9e23SVaclav Hapla       ierr = LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0);CHKERRQ(ierr);
140*1d5c9e23SVaclav Hapla       ierr = LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm);CHKERRQ(ierr);
141*1d5c9e23SVaclav Hapla       ierr = DMCompareLabels(dm0, dm, NULL, NULL);CHKERRQ(ierr);
142*1d5c9e23SVaclav Hapla       ierr = DMDestroy(&dm0);CHKERRQ(ierr);
143*1d5c9e23SVaclav Hapla     } else {
144*1d5c9e23SVaclav Hapla       ierr = LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm);CHKERRQ(ierr);
145*1d5c9e23SVaclav Hapla     }
146*1d5c9e23SVaclav Hapla   } else {
147*1d5c9e23SVaclav Hapla     ierr = DMCreate(options->comm, &dm);CHKERRQ(ierr);
148*1d5c9e23SVaclav Hapla     ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
149*1d5c9e23SVaclav Hapla     ierr = PetscObjectSetName((PetscObject) dm, options->meshname);CHKERRQ(ierr);
150*1d5c9e23SVaclav Hapla     ierr = DMLoad(dm, v);CHKERRQ(ierr);
151*1d5c9e23SVaclav Hapla   }
152*1d5c9e23SVaclav Hapla   ierr = PetscViewerDestroy(&v);CHKERRQ(ierr);
153*1d5c9e23SVaclav Hapla 
154*1d5c9e23SVaclav Hapla   ierr = DMSetOptionsPrefix(dm, "load_");CHKERRQ(ierr);
155*1d5c9e23SVaclav Hapla   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
156*1d5c9e23SVaclav Hapla   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
157*1d5c9e23SVaclav Hapla   *dmnew = dm;
158*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
159*1d5c9e23SVaclav Hapla }
160*1d5c9e23SVaclav Hapla 
161*1d5c9e23SVaclav Hapla static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1)
162*1d5c9e23SVaclav Hapla {
163*1d5c9e23SVaclav Hapla   PetscBool       flg;
164*1d5c9e23SVaclav Hapla   PetscErrorCode  ierr;
165*1d5c9e23SVaclav Hapla 
166*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
167*1d5c9e23SVaclav Hapla   if (options->compare) {
168*1d5c9e23SVaclav Hapla     ierr = DMPlexEqual(dm0, dm1, &flg);CHKERRQ(ierr);
169*1d5c9e23SVaclav Hapla     if (!flg) SETERRQ(options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
170*1d5c9e23SVaclav Hapla     ierr = PetscPrintf(options->comm,"DMs equal\n");CHKERRQ(ierr);
171*1d5c9e23SVaclav Hapla   }
172*1d5c9e23SVaclav Hapla   if (options->compare_labels) {
173*1d5c9e23SVaclav Hapla     ierr = DMCompareLabels(dm0, dm1, NULL, NULL);CHKERRQ(ierr);
174*1d5c9e23SVaclav Hapla     ierr = PetscPrintf(options->comm,"DMLabels equal\n");CHKERRQ(ierr);
175*1d5c9e23SVaclav Hapla   }
176*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
177*1d5c9e23SVaclav Hapla }
178*1d5c9e23SVaclav Hapla 
179*1d5c9e23SVaclav Hapla static PetscErrorCode MarkBoundaryVertices(DM dm, PetscInt value, DMLabel *label)
180*1d5c9e23SVaclav Hapla {
181*1d5c9e23SVaclav Hapla   DMLabel         l;
182*1d5c9e23SVaclav Hapla   IS              points;
183*1d5c9e23SVaclav Hapla   const PetscInt *idx;
184*1d5c9e23SVaclav Hapla   PetscInt        i, n;
185*1d5c9e23SVaclav Hapla   PetscErrorCode  ierr;
186*1d5c9e23SVaclav Hapla 
187*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
188*1d5c9e23SVaclav Hapla   ierr = DMLabelCreate(PetscObjectComm((PetscObject)dm), LABEL_NAME, &l);CHKERRQ(ierr);
189*1d5c9e23SVaclav Hapla   ierr = DMPlexMarkBoundaryFaces(dm, value, l);CHKERRQ(ierr);
190*1d5c9e23SVaclav Hapla   ierr = DMPlexLabelComplete(dm, l);CHKERRQ(ierr);
191*1d5c9e23SVaclav Hapla   ierr = DMLabelGetStratumIS(l, value, &points);CHKERRQ(ierr);
192*1d5c9e23SVaclav Hapla 
193*1d5c9e23SVaclav Hapla   ierr = ISGetLocalSize(points, &n);CHKERRQ(ierr);
194*1d5c9e23SVaclav Hapla   ierr = ISGetIndices(points, &idx);CHKERRQ(ierr);
195*1d5c9e23SVaclav Hapla   for (i=0; i<n; i++) {
196*1d5c9e23SVaclav Hapla     const PetscInt p = idx[i];
197*1d5c9e23SVaclav Hapla     PetscInt       d;
198*1d5c9e23SVaclav Hapla 
199*1d5c9e23SVaclav Hapla     ierr = DMPlexGetPointDepth(dm, p, &d);CHKERRQ(ierr);
200*1d5c9e23SVaclav Hapla     if (d != 0) {
201*1d5c9e23SVaclav Hapla       ierr = DMLabelClearValue(l, p, value);CHKERRQ(ierr);
202*1d5c9e23SVaclav Hapla     }
203*1d5c9e23SVaclav Hapla   }
204*1d5c9e23SVaclav Hapla   ierr = ISRestoreIndices(points, &idx);CHKERRQ(ierr);
205*1d5c9e23SVaclav Hapla   ierr = ISDestroy(&points);CHKERRQ(ierr);
206*1d5c9e23SVaclav Hapla   *label = l;
207*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
208*1d5c9e23SVaclav Hapla }
209*1d5c9e23SVaclav Hapla 
210*1d5c9e23SVaclav Hapla static PetscErrorCode VertexCoordinatesToAll(DM dm, IS vertices, Vec *allCoords)
211*1d5c9e23SVaclav Hapla {
212*1d5c9e23SVaclav Hapla   Vec             coords, allCoords_;
213*1d5c9e23SVaclav Hapla   VecScatter      sc;
214*1d5c9e23SVaclav Hapla   MPI_Comm        comm;
215*1d5c9e23SVaclav Hapla   PetscErrorCode  ierr;
216*1d5c9e23SVaclav Hapla 
217*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
218*1d5c9e23SVaclav Hapla   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
219*1d5c9e23SVaclav Hapla   ierr = DMGetCoordinatesLocalSetUp(dm);CHKERRQ(ierr);
220*1d5c9e23SVaclav Hapla   if (vertices) {
221*1d5c9e23SVaclav Hapla     ierr = DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords);CHKERRQ(ierr);
222*1d5c9e23SVaclav Hapla   } else {
223*1d5c9e23SVaclav Hapla     ierr = VecCreateSeq(PETSC_COMM_SELF, 0, &coords);CHKERRQ(ierr);
224*1d5c9e23SVaclav Hapla   }
225*1d5c9e23SVaclav Hapla   {
226*1d5c9e23SVaclav Hapla     PetscInt  n;
227*1d5c9e23SVaclav Hapla     Vec       mpivec;
228*1d5c9e23SVaclav Hapla 
229*1d5c9e23SVaclav Hapla     ierr = VecGetLocalSize(coords, &n);CHKERRQ(ierr);
230*1d5c9e23SVaclav Hapla     ierr = VecCreateMPI(comm, n, PETSC_DECIDE, &mpivec);CHKERRQ(ierr);
231*1d5c9e23SVaclav Hapla     ierr = VecCopy(coords, mpivec);CHKERRQ(ierr);
232*1d5c9e23SVaclav Hapla     ierr = VecDestroy(&coords);CHKERRQ(ierr);
233*1d5c9e23SVaclav Hapla     coords = mpivec;
234*1d5c9e23SVaclav Hapla   }
235*1d5c9e23SVaclav Hapla 
236*1d5c9e23SVaclav Hapla   ierr = VecScatterCreateToAll(coords, &sc, &allCoords_);CHKERRQ(ierr);
237*1d5c9e23SVaclav Hapla   ierr = VecScatterBegin(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
238*1d5c9e23SVaclav Hapla   ierr = VecScatterEnd(sc,coords,allCoords_,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
239*1d5c9e23SVaclav Hapla   ierr = VecScatterDestroy(&sc);CHKERRQ(ierr);
240*1d5c9e23SVaclav Hapla   ierr = VecDestroy(&coords);CHKERRQ(ierr);
241*1d5c9e23SVaclav Hapla   *allCoords = allCoords_;
242*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
243*1d5c9e23SVaclav Hapla }
244*1d5c9e23SVaclav Hapla 
245*1d5c9e23SVaclav Hapla /* Compute boundary label, remember boundary vertices using coordinates, save and load label, check it is defined on the original boundary vertices */
246*1d5c9e23SVaclav Hapla static PetscErrorCode DMAddBoundaryLabel_GetCoordinateRepresentation(DM dm, Vec *allCoords)
247*1d5c9e23SVaclav Hapla {
248*1d5c9e23SVaclav Hapla   DMLabel         label;
249*1d5c9e23SVaclav Hapla   IS              vertices;
250*1d5c9e23SVaclav Hapla   PetscErrorCode  ierr;
251*1d5c9e23SVaclav Hapla 
252*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
253*1d5c9e23SVaclav Hapla   ierr = MarkBoundaryVertices(dm, LABEL_VALUE, &label);CHKERRQ(ierr);
254*1d5c9e23SVaclav Hapla   ierr = DMLabelGetStratumIS(label, LABEL_VALUE, &vertices);CHKERRQ(ierr);
255*1d5c9e23SVaclav Hapla   ierr = VertexCoordinatesToAll(dm, vertices, allCoords);CHKERRQ(ierr);
256*1d5c9e23SVaclav Hapla   ierr = DMAddLabel(dm, label);CHKERRQ(ierr);
257*1d5c9e23SVaclav Hapla   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
258*1d5c9e23SVaclav Hapla   ierr = ISDestroy(&vertices);CHKERRQ(ierr);
259*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
260*1d5c9e23SVaclav Hapla }
261*1d5c9e23SVaclav Hapla 
262*1d5c9e23SVaclav Hapla static PetscErrorCode DMGetBoundaryLabel_CompareWithCoordinateRepresentation(AppCtx *user, DM dm, Vec allCoords)
263*1d5c9e23SVaclav Hapla {
264*1d5c9e23SVaclav Hapla   DMLabel         label;
265*1d5c9e23SVaclav Hapla   IS              pointsIS;
266*1d5c9e23SVaclav Hapla   const PetscInt *points;
267*1d5c9e23SVaclav Hapla   PetscInt        i, n;
268*1d5c9e23SVaclav Hapla   PetscBool       fail = PETSC_FALSE;
269*1d5c9e23SVaclav Hapla   MPI_Comm        comm;
270*1d5c9e23SVaclav Hapla   PetscMPIInt     rank;
271*1d5c9e23SVaclav Hapla   PetscErrorCode  ierr;
272*1d5c9e23SVaclav Hapla 
273*1d5c9e23SVaclav Hapla   PetscFunctionBeginUser;
274*1d5c9e23SVaclav Hapla   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
275*1d5c9e23SVaclav Hapla   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
276*1d5c9e23SVaclav Hapla   ierr = DMGetLabel(dm, LABEL_NAME, &label);CHKERRQ(ierr);
277*1d5c9e23SVaclav Hapla   if (!label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label \"%s\" was not loaded", LABEL_NAME);
278*1d5c9e23SVaclav Hapla   {
279*1d5c9e23SVaclav Hapla     PetscInt pStart, pEnd;
280*1d5c9e23SVaclav Hapla 
281*1d5c9e23SVaclav Hapla     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
282*1d5c9e23SVaclav Hapla     ierr = DMLabelCreateIndex(label, pStart, pEnd);CHKERRQ(ierr);
283*1d5c9e23SVaclav Hapla   }
284*1d5c9e23SVaclav Hapla   ierr = DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS);CHKERRQ(ierr);
285*1d5c9e23SVaclav Hapla   ierr = ISGetIndices(pointsIS, &points);CHKERRQ(ierr);
286*1d5c9e23SVaclav Hapla   ierr = ISGetLocalSize(pointsIS, &n);CHKERRQ(ierr);
287*1d5c9e23SVaclav Hapla   if (user->verbose) {ierr = DMLabelView(label, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);}
288*1d5c9e23SVaclav Hapla   for (i=0; i<n; i++) {
289*1d5c9e23SVaclav Hapla     const PetscInt  p = points[i];
290*1d5c9e23SVaclav Hapla     PetscBool       has;
291*1d5c9e23SVaclav Hapla     PetscInt        v;
292*1d5c9e23SVaclav Hapla 
293*1d5c9e23SVaclav Hapla     if (p < 0) continue;
294*1d5c9e23SVaclav Hapla     ierr = DMLabelHasPoint(label, p, &has);CHKERRQ(ierr);
295*1d5c9e23SVaclav Hapla     if (!has) {
296*1d5c9e23SVaclav Hapla       ierr = PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label does not have point %D\n", rank, p);CHKERRQ(ierr);
297*1d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
298*1d5c9e23SVaclav Hapla       continue;
299*1d5c9e23SVaclav Hapla     }
300*1d5c9e23SVaclav Hapla     ierr = DMLabelGetValue(label, p, &v);CHKERRQ(ierr);
301*1d5c9e23SVaclav Hapla     if (v != LABEL_VALUE) {
302*1d5c9e23SVaclav Hapla       ierr = PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %D has bad value %D", p, v);CHKERRQ(ierr);
303*1d5c9e23SVaclav Hapla       fail = PETSC_TRUE;
304*1d5c9e23SVaclav Hapla       continue;
305*1d5c9e23SVaclav Hapla     }
306*1d5c9e23SVaclav Hapla     if (user->verbose) {ierr = PetscSynchronizedPrintf(comm, "[%d] OK point %D\n", rank, p);CHKERRQ(ierr);}
307*1d5c9e23SVaclav Hapla   }
308*1d5c9e23SVaclav Hapla   ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
309*1d5c9e23SVaclav Hapla   ierr = PetscSynchronizedFlush(comm, PETSC_STDERR);CHKERRQ(ierr);
310*1d5c9e23SVaclav Hapla   ierr = ISRestoreIndices(pointsIS, &points);CHKERRQ(ierr);
311*1d5c9e23SVaclav Hapla   ierr = ISDestroy(&pointsIS);CHKERRQ(ierr);
312*1d5c9e23SVaclav Hapla   ierr = MPI_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm);CHKERRMPI(ierr);
313*1d5c9e23SVaclav Hapla   if (fail) SETERRQ1(comm, PETSC_ERR_PLIB, "Label \"%s\" was not loaded correctly - see details above", LABEL_NAME);
314*1d5c9e23SVaclav Hapla   PetscFunctionReturn(0);
315*1d5c9e23SVaclav Hapla }
316*1d5c9e23SVaclav Hapla 
317*1d5c9e23SVaclav Hapla int main(int argc, char **argv)
318*1d5c9e23SVaclav Hapla {
319*1d5c9e23SVaclav Hapla   DM             dm, dmnew;
320*1d5c9e23SVaclav Hapla   AppCtx         user;
321*1d5c9e23SVaclav Hapla   Vec            allCoords = NULL;
322*1d5c9e23SVaclav Hapla   PetscErrorCode ierr;
323*1d5c9e23SVaclav Hapla 
324*1d5c9e23SVaclav Hapla   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
325*1d5c9e23SVaclav Hapla   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
326*1d5c9e23SVaclav Hapla   ierr = CreateMesh(&user, &dm);CHKERRQ(ierr);
327*1d5c9e23SVaclav Hapla   if (user.compare_boundary) {
328*1d5c9e23SVaclav Hapla     ierr = DMAddBoundaryLabel_GetCoordinateRepresentation(dm, &allCoords);CHKERRQ(ierr);
329*1d5c9e23SVaclav Hapla   }
330*1d5c9e23SVaclav Hapla   ierr = SaveMesh(&user, dm);CHKERRQ(ierr);
331*1d5c9e23SVaclav Hapla   ierr = LoadMesh(&user, &dmnew);CHKERRQ(ierr);
332*1d5c9e23SVaclav Hapla   ierr = CompareMeshes(&user, dm, dmnew);CHKERRQ(ierr);
333*1d5c9e23SVaclav Hapla   if (user.compare_boundary) {
334*1d5c9e23SVaclav Hapla     ierr = DMGetBoundaryLabel_CompareWithCoordinateRepresentation(&user, dmnew, allCoords);CHKERRQ(ierr);
335*1d5c9e23SVaclav Hapla   }
336*1d5c9e23SVaclav Hapla   ierr = DMDestroy(&dm);CHKERRQ(ierr);
337*1d5c9e23SVaclav Hapla   ierr = DMDestroy(&dmnew);CHKERRQ(ierr);
338*1d5c9e23SVaclav Hapla   ierr = VecDestroy(&allCoords);CHKERRQ(ierr);
339*1d5c9e23SVaclav Hapla   ierr = PetscFinalize();
340*1d5c9e23SVaclav Hapla   return ierr;
341*1d5c9e23SVaclav Hapla }
342*1d5c9e23SVaclav Hapla 
343*1d5c9e23SVaclav Hapla //TODO we can -compare once the new parallel topology format is in place
344*1d5c9e23SVaclav Hapla /*TEST
345*1d5c9e23SVaclav Hapla   build:
346*1d5c9e23SVaclav Hapla     requires: hdf5
347*1d5c9e23SVaclav Hapla 
348*1d5c9e23SVaclav Hapla   # load old format, save in new format, reload, distribute
349*1d5c9e23SVaclav Hapla   testset:
350*1d5c9e23SVaclav Hapla     suffix: 1
351*1d5c9e23SVaclav Hapla     requires: !complex datafilespath
352*1d5c9e23SVaclav Hapla     args: -dm_plex_name plex
353*1d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
354*1d5c9e23SVaclav Hapla     args: -dm_distribute -dm_plex_interpolate
355*1d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
356*1d5c9e23SVaclav Hapla     args: -use_low_level_functions {{0 1}} -load_dm_distribute -compare_boundary
357*1d5c9e23SVaclav Hapla     args: -outfile ex56_1.h5
358*1d5c9e23SVaclav Hapla     nsize: {{1 3}}
359*1d5c9e23SVaclav Hapla     test:
360*1d5c9e23SVaclav Hapla       suffix: a
361*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
362*1d5c9e23SVaclav Hapla     test:
363*1d5c9e23SVaclav Hapla       suffix: b
364*1d5c9e23SVaclav Hapla       TODO: broken
365*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
366*1d5c9e23SVaclav Hapla     test:
367*1d5c9e23SVaclav Hapla       suffix: c
368*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
369*1d5c9e23SVaclav Hapla     test:
370*1d5c9e23SVaclav Hapla       suffix: d
371*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
372*1d5c9e23SVaclav Hapla     test:
373*1d5c9e23SVaclav Hapla       suffix: e
374*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
375*1d5c9e23SVaclav Hapla     test:
376*1d5c9e23SVaclav Hapla       suffix: f
377*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
378*1d5c9e23SVaclav Hapla 
379*1d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
380*1d5c9e23SVaclav Hapla   testset:
381*1d5c9e23SVaclav Hapla     suffix: 2
382*1d5c9e23SVaclav Hapla     requires: !complex datafilespath
383*1d5c9e23SVaclav Hapla     args: -dm_plex_name plex
384*1d5c9e23SVaclav Hapla     args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
385*1d5c9e23SVaclav Hapla     args: -dm_distribute -dm_plex_interpolate
386*1d5c9e23SVaclav Hapla     args: -load_dm_plex_check_all
387*1d5c9e23SVaclav Hapla     args: -use_low_level_functions -distribute_after_topo_load -compare_boundary
388*1d5c9e23SVaclav Hapla     args: -outfile ex56_2.h5
389*1d5c9e23SVaclav Hapla     nsize: 3
390*1d5c9e23SVaclav Hapla     test:
391*1d5c9e23SVaclav Hapla       suffix: a
392*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
393*1d5c9e23SVaclav Hapla     test:
394*1d5c9e23SVaclav Hapla       suffix: b
395*1d5c9e23SVaclav Hapla       TODO: broken
396*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
397*1d5c9e23SVaclav Hapla     test:
398*1d5c9e23SVaclav Hapla       suffix: c
399*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
400*1d5c9e23SVaclav Hapla     test:
401*1d5c9e23SVaclav Hapla       suffix: d
402*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
403*1d5c9e23SVaclav Hapla     test:
404*1d5c9e23SVaclav Hapla       suffix: e
405*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
406*1d5c9e23SVaclav Hapla     test:
407*1d5c9e23SVaclav Hapla       suffix: f
408*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
409*1d5c9e23SVaclav Hapla 
410*1d5c9e23SVaclav Hapla   # load old format, save in new format, reload topology, distribute, load geometry and labels
411*1d5c9e23SVaclav Hapla   testset:
412*1d5c9e23SVaclav Hapla     suffix: 3
413*1d5c9e23SVaclav Hapla     requires: !complex datafilespath
414*1d5c9e23SVaclav Hapla     args: -dm_plex_name plex
415*1d5c9e23SVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 2.0.0
416*1d5c9e23SVaclav Hapla     args: -dm_distribute -dm_plex_interpolate
417*1d5c9e23SVaclav Hapla     args: -use_low_level_functions -compare_pre_post
418*1d5c9e23SVaclav Hapla     args: -outfile ex56_3.h5
419*1d5c9e23SVaclav Hapla     nsize: 3
420*1d5c9e23SVaclav Hapla     test:
421*1d5c9e23SVaclav Hapla       suffix: a
422*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
423*1d5c9e23SVaclav Hapla     test:
424*1d5c9e23SVaclav Hapla       suffix: b
425*1d5c9e23SVaclav Hapla       TODO: broken
426*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
427*1d5c9e23SVaclav Hapla     test:
428*1d5c9e23SVaclav Hapla       suffix: c
429*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
430*1d5c9e23SVaclav Hapla     test:
431*1d5c9e23SVaclav Hapla       suffix: d
432*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
433*1d5c9e23SVaclav Hapla     test:
434*1d5c9e23SVaclav Hapla       suffix: e
435*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
436*1d5c9e23SVaclav Hapla     test:
437*1d5c9e23SVaclav Hapla       suffix: f
438*1d5c9e23SVaclav Hapla       args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
439*1d5c9e23SVaclav Hapla TEST*/
440