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