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