Actual source code: ex56.c
1: #include <petscdmplex.h>
2: #include <petscviewerhdf5.h>
3: #include <petscsf.h>
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;
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;
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: };
50: static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm)
51: {
52: DM dm;
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: }
64: static PetscErrorCode SaveMesh(AppCtx *options, DM dm)
65: {
66: PetscViewer v;
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: }
81: typedef enum {
82: NONE = 0,
83: PRE_DIST = 1,
84: POST_DIST = 2
85: } AuxObjLoadMode;
87: static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm)
88: {
89: DM dm;
90: PetscSF sfXC;
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;
105: PetscCall(DMPlexDistribute(dm, 0, &sfBC, &dmdist));
106: if (dmdist) {
107: const char *name;
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(DMPlexLabelsLoad(dm, v, sfXC));
120: PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC));
121: }
122: PetscCall(PetscSFDestroy(&sfXC));
123: *newdm = dm;
124: PetscFunctionReturn(PETSC_SUCCESS);
125: }
127: static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew)
128: {
129: DM dm;
130: PetscViewer v;
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;
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));
153: PetscCall(DMSetOptionsPrefix(dm, "load_"));
154: PetscCall(DMSetFromOptions(dm));
155: PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
156: *dmnew = dm;
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1)
161: {
162: PetscBool flg;
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: }
177: static PetscErrorCode MarkBoundary(DM dm, const char name[], PetscInt value, PetscBool verticesOnly, DMLabel *label)
178: {
179: DMLabel l;
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 = 0;
191: PetscCall(DMLabelGetStratumIS(l, value, &points));
192: if (points) {
193: PetscCall(ISGetLocalSize(points, &n));
194: PetscCall(ISGetIndices(points, &idx));
195: }
196: for (i = 0; i < n; i++) {
197: const PetscInt p = idx[i];
198: PetscInt d;
200: PetscCall(DMPlexGetPointDepth(dm, p, &d));
201: if (d != 0) PetscCall(DMLabelClearValue(l, p, value));
202: }
203: if (points) PetscCall(ISRestoreIndices(points, &idx));
204: PetscCall(ISDestroy(&points));
205: }
206: if (label) *label = l;
207: else PetscCall(DMLabelDestroy(&l));
208: PetscFunctionReturn(PETSC_SUCCESS);
209: }
211: static PetscErrorCode VertexCoordinatesToAll(DM dm, IS vertices, Vec *allCoords)
212: {
213: Vec coords, allCoords_;
214: VecScatter sc;
215: MPI_Comm comm;
217: PetscFunctionBeginUser;
218: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
219: PetscCall(DMGetCoordinatesLocalSetUp(dm));
220: if (vertices) {
221: PetscCall(DMGetCoordinatesLocalTuple(dm, vertices, NULL, &coords));
222: } else {
223: PetscCall(VecCreateFromOptions(PETSC_COMM_SELF, NULL, 1, 0, 0, &coords));
224: }
225: {
226: PetscInt n;
227: Vec mpivec;
229: PetscCall(VecGetLocalSize(coords, &n));
230: PetscCall(VecCreateFromOptions(comm, NULL, 1, n, PETSC_DECIDE, &mpivec));
231: PetscCall(VecCopy(coords, mpivec));
232: PetscCall(VecDestroy(&coords));
233: coords = mpivec;
234: }
236: PetscCall(VecScatterCreateToAll(coords, &sc, &allCoords_));
237: PetscCall(VecScatterBegin(sc, coords, allCoords_, INSERT_VALUES, SCATTER_FORWARD));
238: PetscCall(VecScatterEnd(sc, coords, allCoords_, INSERT_VALUES, SCATTER_FORWARD));
239: PetscCall(VecScatterDestroy(&sc));
240: PetscCall(VecDestroy(&coords));
241: *allCoords = allCoords_;
242: PetscFunctionReturn(PETSC_SUCCESS);
243: }
245: static PetscErrorCode DMLabelGetCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec *allCoords)
246: {
247: IS vertices;
249: PetscFunctionBeginUser;
250: PetscCall(DMLabelGetStratumIS(label, value, &vertices));
251: PetscCall(VertexCoordinatesToAll(dm, vertices, allCoords));
252: PetscCall(ISDestroy(&vertices));
253: PetscFunctionReturn(PETSC_SUCCESS);
254: }
256: static PetscErrorCode DMLabelCompareWithCoordinateRepresentation(DM dm, DMLabel label, PetscInt value, Vec allCoords, PetscInt verbose)
257: {
258: const char *labelName;
259: IS pointsIS;
260: const PetscInt *points;
261: PetscInt i, n;
262: PetscBool fail = PETSC_FALSE;
263: MPI_Comm comm;
264: PetscMPIInt rank;
266: PetscFunctionBeginUser;
267: PetscCheck(label, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label not loaded");
268: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
269: PetscCall(PetscObjectGetName((PetscObject)label, &labelName));
270: PetscCallMPI(MPI_Comm_rank(comm, &rank));
271: {
272: PetscInt pStart, pEnd;
274: PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
275: PetscCall(DMLabelCreateIndex(label, pStart, pEnd));
276: }
277: PetscCall(DMPlexFindVertices(dm, allCoords, 0.0, &pointsIS));
278: PetscCall(ISGetIndices(pointsIS, &points));
279: PetscCall(ISGetLocalSize(pointsIS, &n));
280: if (verbose > 1) PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_(comm)));
281: for (i = 0; i < n; i++) {
282: const PetscInt p = points[i];
283: PetscBool has;
284: PetscInt v;
286: if (p < 0) continue;
287: PetscCall(DMLabelHasPoint(label, p, &has));
288: if (!has) {
289: if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label \"%s\" does not have point %" PetscInt_FMT "\n", rank, labelName, p));
290: fail = PETSC_TRUE;
291: continue;
292: }
293: PetscCall(DMLabelGetValue(label, p, &v));
294: if (v != value) {
295: if (verbose) PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "Point %" PetscInt_FMT " has bad value %" PetscInt_FMT " in label \"%s\"", p, v, labelName));
296: fail = PETSC_TRUE;
297: continue;
298: }
299: if (verbose > 1) PetscCall(PetscSynchronizedPrintf(comm, "[%d] OK point %" PetscInt_FMT "\n", rank, p));
300: }
301: PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
302: PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
303: PetscCall(ISRestoreIndices(pointsIS, &points));
304: PetscCall(ISDestroy(&pointsIS));
305: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
306: PetscCheck(!fail, comm, PETSC_ERR_PLIB, "Label \"%s\" was not loaded correctly%s", labelName, verbose ? " - see details above" : "");
307: PetscFunctionReturn(PETSC_SUCCESS);
308: }
310: static PetscErrorCode CheckNumLabels(DM dm, AppCtx *ctx)
311: {
312: PetscInt actualNum;
313: PetscBool fail = PETSC_FALSE;
314: MPI_Comm comm;
315: PetscMPIInt rank;
317: PetscFunctionBeginUser;
318: if (ctx->num_labels < 0) PetscFunctionReturn(PETSC_SUCCESS);
319: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
320: PetscCallMPI(MPI_Comm_rank(comm, &rank));
321: PetscCall(DMGetNumLabels(dm, &actualNum));
322: if (ctx->num_labels != actualNum) {
323: fail = PETSC_TRUE;
324: if (ctx->verbose) {
325: PetscInt i;
327: PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Asserted number of labels: %" PetscInt_FMT ", actual: %" PetscInt_FMT "\n", rank, ctx->num_labels, actualNum));
328: for (i = 0; i < actualNum; i++) {
329: DMLabel label;
330: const char *name;
332: PetscCall(DMGetLabelByNum(dm, i, &label));
333: PetscCall(PetscObjectGetName((PetscObject)label, &name));
334: PetscCall(PetscSynchronizedFPrintf(comm, PETSC_STDERR, "[%d] Label %" PetscInt_FMT " \"%s\"\n", rank, i, name));
335: }
336: PetscCall(PetscSynchronizedFlush(comm, PETSC_STDERR));
337: }
338: }
339: PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &fail, 1, MPIU_BOOL, MPI_LOR, comm));
340: PetscCheck(!fail, comm, PETSC_ERR_PLIB, "Wrong number of labels%s", ctx->verbose ? " - see details above" : "");
341: PetscFunctionReturn(PETSC_SUCCESS);
342: }
344: static inline PetscErrorCode IncrementNumLabels(AppCtx *ctx)
345: {
346: PetscFunctionBeginUser;
347: if (ctx->num_labels >= 0) ctx->num_labels++;
348: PetscFunctionReturn(PETSC_SUCCESS);
349: }
351: int main(int argc, char **argv)
352: {
353: const char BOUNDARY_NAME[] = "Boundary";
354: const char BOUNDARY_VERTICES_NAME[] = "BoundaryVertices";
355: const PetscInt BOUNDARY_VALUE = 12345;
356: const PetscInt BOUNDARY_VERTICES_VALUE = 6789;
357: DM dm, dmnew;
358: AppCtx user;
359: Vec allCoords = NULL;
361: PetscFunctionBeginUser;
362: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
363: PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
364: PetscCall(CreateMesh(&user, &dm));
365: PetscCall(MarkBoundary(dm, BOUNDARY_NAME, BOUNDARY_VALUE, PETSC_FALSE, NULL));
366: PetscCall(IncrementNumLabels(&user));
367: if (user.compare_boundary) {
368: DMLabel label;
370: PetscCall(MarkBoundary(dm, BOUNDARY_VERTICES_NAME, BOUNDARY_VERTICES_VALUE, PETSC_TRUE, &label));
371: PetscCall(IncrementNumLabels(&user));
372: PetscCall(DMLabelGetCoordinateRepresentation(dm, label, BOUNDARY_VERTICES_VALUE, &allCoords));
373: PetscCall(DMLabelDestroy(&label));
374: }
375: PetscCall(SaveMesh(&user, dm));
377: PetscCall(LoadMesh(&user, &dmnew));
378: PetscCall(IncrementNumLabels(&user)); /* account for depth label */
379: PetscCall(IncrementNumLabels(&user)); /* account for celltype label */
380: PetscCall(CheckNumLabels(dm, &user));
381: PetscCall(CompareMeshes(&user, dm, dmnew));
382: if (user.compare_boundary) {
383: DMLabel label;
385: PetscCall(DMGetLabel(dmnew, BOUNDARY_VERTICES_NAME, &label));
386: PetscCall(DMLabelCompareWithCoordinateRepresentation(dmnew, label, BOUNDARY_VERTICES_VALUE, allCoords, user.verbose));
387: }
388: PetscCall(DMDestroy(&dm));
389: PetscCall(DMDestroy(&dmnew));
390: PetscCall(VecDestroy(&allCoords));
391: PetscCall(PetscFinalize());
392: return 0;
393: }
395: //TODO we can -compare once the new parallel topology format is in place
396: /*TEST
397: build:
398: requires: hdf5
400: # load old format, save in new format, reload, distribute
401: testset:
402: suffix: 1
403: requires: !complex datafilespath
404: args: -dm_plex_name plex
405: args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
406: args: -dm_plex_interpolate -petscpartitioner_type simple
407: args: -load_dm_plex_check_all
408: args: -use_low_level_functions {{0 1}} -compare_boundary
409: args: -num_labels 1
410: args: -outfile ex56_1.h5
411: nsize: {{1 3}}
412: test:
413: suffix: a
414: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
415: test:
416: suffix: b
417: TODO: broken
418: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
419: test:
420: suffix: c
421: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
422: test:
423: suffix: d
424: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
425: test:
426: suffix: e
427: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
428: test:
429: suffix: f
430: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
432: # load old format, save in new format, reload topology, distribute, load geometry and labels
433: testset:
434: suffix: 2
435: requires: !complex datafilespath
436: args: -dm_plex_name plex
437: args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0
438: args: -dm_plex_interpolate -petscpartitioner_type simple
439: args: -load_dm_plex_check_all
440: args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary
441: args: -num_labels 1
442: args: -outfile ex56_2.h5
443: nsize: 3
444: test:
445: suffix: a
446: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
447: test:
448: suffix: b
449: TODO: broken
450: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
451: test:
452: suffix: c
453: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
454: test:
455: suffix: d
456: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
457: test:
458: suffix: e
459: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
460: test:
461: suffix: f
462: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
464: # load old format, save in new format, reload topology, distribute, load geometry and labels
465: testset:
466: suffix: 3
467: requires: !complex datafilespath
468: args: -dm_plex_name plex
469: args: -dm_plex_view_hdf5_storage_version 2.0.0
470: args: -dm_plex_interpolate -load_dm_distribute 0 -petscpartitioner_type simple
471: args: -use_low_level_functions -compare_pre_post
472: args: -num_labels 1
473: args: -outfile ex56_3.h5
474: nsize: 3
475: test:
476: suffix: a
477: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
478: test:
479: suffix: b
480: TODO: broken
481: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
482: test:
483: suffix: c
484: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
485: test:
486: suffix: d
487: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0
488: test:
489: suffix: e
490: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
491: test:
492: suffix: f
493: args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
494: TEST*/