1 #include <petscdmplex.h> 2 #include <petscviewerhdf5.h> 3 #include <petscsf.h> 4 5 static const char help[] = "Load and save the mesh to the native HDF5 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 PetscErrorCode ierr; 25 26 PetscFunctionBeginUser; 27 options->comm = comm; 28 options->num_labels = -1; 29 options->compare = PETSC_FALSE; 30 options->compare_labels = PETSC_FALSE; 31 options->compare_boundary = PETSC_FALSE; 32 options->compare_pre_post = PETSC_FALSE; 33 options->outfile[0] = '\0'; 34 options->use_low_level_functions = PETSC_FALSE; 35 options->distribute_after_topo_load = PETSC_FALSE; 36 options->verbose = 0; 37 38 ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");PetscCall(ierr); 39 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)); 40 PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual() and DMCompareLabels()", EX, options->compare, &options->compare, NULL)); 41 PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL)); 42 PetscCall(PetscOptionsBool("-compare_boundary", "Check label I/O via boundary vertex coordinates", "ex55.c", options->compare_boundary, &options->compare_boundary, NULL)); 43 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)); 44 PetscCall(PetscOptionsString("-outfile", "Output mesh file", EX, options->outfile, options->outfile, sizeof(options->outfile), NULL)); 45 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)); 46 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)); 47 PetscCall(PetscOptionsInt("-verbose", "Verbosity level", EX, options->verbose, &options->verbose, NULL)); 48 ierr = PetscOptionsEnd();PetscCall(ierr); 49 PetscFunctionReturn(0); 50 }; 51 52 static PetscErrorCode CreateMesh(AppCtx *options, DM *newdm) 53 { 54 DM dm; 55 56 PetscFunctionBeginUser; 57 PetscCall(DMCreate(options->comm, &dm)); 58 PetscCall(DMSetType(dm, DMPLEX)); 59 PetscCall(DMSetFromOptions(dm)); 60 PetscCall(PetscObjectGetName((PetscObject)dm, &options->meshname)); 61 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 62 *newdm = dm; 63 PetscFunctionReturn(0); 64 } 65 66 static PetscErrorCode SaveMesh(AppCtx *options, DM dm) 67 { 68 PetscViewer v; 69 70 PetscFunctionBeginUser; 71 PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), options->outfile, FILE_MODE_WRITE, &v)); 72 if (options->use_low_level_functions) { 73 PetscCall(DMPlexTopologyView(dm, v)); 74 PetscCall(DMPlexCoordinatesView(dm, v)); 75 PetscCall(DMPlexLabelsView(dm, v)); 76 } else { 77 PetscCall(DMView(dm, v)); 78 } 79 PetscCall(PetscViewerDestroy(&v)); 80 PetscFunctionReturn(0); 81 } 82 83 typedef enum {NONE=0, PRE_DIST=1, POST_DIST=2} AuxObjLoadMode; 84 85 static PetscErrorCode LoadMeshLowLevel(AppCtx *options, PetscViewer v, PetscBool explicitDistribute, AuxObjLoadMode mode, DM *newdm) 86 { 87 DM dm; 88 PetscSF sfXC; 89 90 PetscFunctionBeginUser; 91 PetscCall(DMCreate(options->comm, &dm)); 92 PetscCall(DMSetType(dm, DMPLEX)); 93 PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname)); 94 PetscCall(DMPlexTopologyLoad(dm, v, &sfXC)); 95 if (mode == PRE_DIST) { 96 PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC)); 97 PetscCall(DMPlexLabelsLoad(dm, v, sfXC)); 98 } 99 if (explicitDistribute) { 100 DM dmdist; 101 PetscSF sfXB = sfXC, sfBC; 102 103 PetscCall(DMPlexDistribute(dm, 0, &sfBC, &dmdist)); 104 if (dmdist) { 105 const char *name; 106 107 PetscCall(PetscObjectGetName((PetscObject) dm, &name)); 108 PetscCall(PetscObjectSetName((PetscObject) dmdist, name)); 109 PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC)); 110 PetscCall(PetscSFDestroy(&sfXB)); 111 PetscCall(PetscSFDestroy(&sfBC)); 112 PetscCall(DMDestroy(&dm)); 113 dm = dmdist; 114 } 115 } 116 if (mode == POST_DIST) { 117 PetscCall(DMPlexCoordinatesLoad(dm, v, sfXC)); 118 PetscCall(DMPlexLabelsLoad(dm, v, sfXC)); 119 } 120 PetscCall(PetscSFDestroy(&sfXC)); 121 *newdm = dm; 122 PetscFunctionReturn(0); 123 } 124 125 static PetscErrorCode LoadMesh(AppCtx *options, DM *dmnew) 126 { 127 DM dm; 128 PetscViewer v; 129 130 PetscFunctionBeginUser; 131 PetscCall(PetscViewerHDF5Open(options->comm, options->outfile, FILE_MODE_READ, &v)); 132 if (options->use_low_level_functions) { 133 if (options->compare_pre_post) { 134 DM dm0; 135 136 PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, PRE_DIST, &dm0)); 137 PetscCall(LoadMeshLowLevel(options, v, PETSC_TRUE, POST_DIST, &dm)); 138 PetscCall(DMCompareLabels(dm0, dm, NULL, NULL)); 139 PetscCall(DMDestroy(&dm0)); 140 } else { 141 PetscCall(LoadMeshLowLevel(options, v, options->distribute_after_topo_load, POST_DIST, &dm)); 142 } 143 } else { 144 PetscCall(DMCreate(options->comm, &dm)); 145 PetscCall(DMSetType(dm, DMPLEX)); 146 PetscCall(PetscObjectSetName((PetscObject) dm, options->meshname)); 147 PetscCall(DMLoad(dm, v)); 148 } 149 PetscCall(PetscViewerDestroy(&v)); 150 151 PetscCall(DMSetOptionsPrefix(dm, "load_")); 152 PetscCall(DMSetFromOptions(dm)); 153 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 154 *dmnew = dm; 155 PetscFunctionReturn(0); 156 } 157 158 static PetscErrorCode CompareMeshes(AppCtx *options, DM dm0, DM dm1) 159 { 160 PetscBool flg; 161 162 PetscFunctionBeginUser; 163 if (options->compare) { 164 PetscCall(DMPlexEqual(dm0, dm1, &flg)); 165 PetscCheck(flg,options->comm, PETSC_ERR_ARG_INCOMP, "DMs are not equal"); 166 PetscCall(PetscPrintf(options->comm,"DMs equal\n")); 167 } 168 if (options->compare_labels) { 169 PetscCall(DMCompareLabels(dm0, dm1, NULL, NULL)); 170 PetscCall(PetscPrintf(options->comm,"DMLabels equal\n")); 171 } 172 PetscFunctionReturn(0); 173 } 174 175 static PetscErrorCode MarkBoundary(DM dm, const char name[], PetscInt value, PetscBool verticesOnly, DMLabel *label) 176 { 177 DMLabel l; 178 179 PetscFunctionBeginUser; 180 PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)dm), name, &l)); 181 PetscCall(DMAddLabel(dm, l)); 182 PetscCall(DMPlexMarkBoundaryFaces(dm, value, l)); 183 PetscCall(DMPlexLabelComplete(dm, l)); 184 if (verticesOnly) { 185 IS points; 186 const PetscInt *idx; 187 PetscInt i, n; 188 189 PetscCall(DMLabelGetStratumIS(l, value, &points)); 190 PetscCall(ISGetLocalSize(points, &n)); 191 PetscCall(ISGetIndices(points, &idx)); 192 for (i=0; i<n; i++) { 193 const PetscInt p = idx[i]; 194 PetscInt d; 195 196 PetscCall(DMPlexGetPointDepth(dm, p, &d)); 197 if (d != 0) { 198 PetscCall(DMLabelClearValue(l, p, value)); 199 } 200 } 201 PetscCall(ISRestoreIndices(points, &idx)); 202 PetscCall(ISDestroy(&points)); 203 } 204 if (label) *label = l; 205 else PetscCall(DMLabelDestroy(&l)); 206 PetscFunctionReturn(0); 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(0); 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(0); 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 %D\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 %D has bad value %D in label \"%s\"", p, v, labelName)); 294 fail = PETSC_TRUE; 295 continue; 296 } 297 if (verbose > 1) PetscCall(PetscSynchronizedPrintf(comm, "[%d] OK point %D\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(0); 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(0); 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: %D, actual: %D\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 %D \"%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(0); 340 } 341 342 static inline PetscErrorCode IncrementNumLabels(AppCtx *ctx) 343 { 344 PetscFunctionBeginUser; 345 if (ctx->num_labels >= 0) ctx->num_labels++; 346 PetscFunctionReturn(0); 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 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 360 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 361 PetscCall(CreateMesh(&user, &dm)); 362 PetscCall(MarkBoundary(dm, BOUNDARY_NAME, BOUNDARY_VALUE, PETSC_FALSE, NULL)); 363 PetscCall(IncrementNumLabels(&user)); 364 if (user.compare_boundary) { 365 DMLabel label; 366 367 PetscCall(MarkBoundary(dm, BOUNDARY_VERTICES_NAME, BOUNDARY_VERTICES_VALUE, PETSC_TRUE, &label)); 368 PetscCall(IncrementNumLabels(&user)); 369 PetscCall(DMLabelGetCoordinateRepresentation(dm, label, BOUNDARY_VERTICES_VALUE, &allCoords)); 370 PetscCall(DMLabelDestroy(&label)); 371 } 372 PetscCall(SaveMesh(&user, dm)); 373 374 PetscCall(LoadMesh(&user, &dmnew)); 375 PetscCall(IncrementNumLabels(&user)); /* account for depth label */ 376 PetscCall(IncrementNumLabels(&user)); /* account for celltype label */ 377 PetscCall(CheckNumLabels(dm, &user)); 378 PetscCall(CompareMeshes(&user, dm, dmnew)); 379 if (user.compare_boundary) { 380 DMLabel label; 381 382 PetscCall(DMGetLabel(dmnew, BOUNDARY_VERTICES_NAME, &label)); 383 PetscCall(DMLabelCompareWithCoordinateRepresentation(dmnew, label, BOUNDARY_VERTICES_VALUE, allCoords, user.verbose)); 384 } 385 PetscCall(DMDestroy(&dm)); 386 PetscCall(DMDestroy(&dmnew)); 387 PetscCall(VecDestroy(&allCoords)); 388 PetscCall(PetscFinalize()); 389 return 0; 390 } 391 392 //TODO we can -compare once the new parallel topology format is in place 393 /*TEST 394 build: 395 requires: hdf5 396 397 # load old format, save in new format, reload, distribute 398 testset: 399 suffix: 1 400 requires: !complex datafilespath 401 args: -dm_plex_name plex 402 args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0 403 args: -dm_plex_interpolate 404 args: -load_dm_plex_check_all 405 args: -use_low_level_functions {{0 1}} -compare_boundary 406 args: -num_labels 1 407 args: -outfile ex56_1.h5 408 nsize: {{1 3}} 409 test: 410 suffix: a 411 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 412 test: 413 suffix: b 414 TODO: broken 415 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 416 test: 417 suffix: c 418 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 419 test: 420 suffix: d 421 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0 422 test: 423 suffix: e 424 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 425 test: 426 suffix: f 427 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 428 429 # load old format, save in new format, reload topology, distribute, load geometry and labels 430 testset: 431 suffix: 2 432 requires: !complex datafilespath 433 args: -dm_plex_name plex 434 args: -dm_plex_check_all -dm_plex_view_hdf5_storage_version 2.0.0 435 args: -dm_plex_interpolate 436 args: -load_dm_plex_check_all 437 args: -use_low_level_functions -load_dm_distribute 0 -distribute_after_topo_load -compare_boundary 438 args: -num_labels 1 439 args: -outfile ex56_2.h5 440 nsize: 3 441 test: 442 suffix: a 443 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 444 test: 445 suffix: b 446 TODO: broken 447 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 448 test: 449 suffix: c 450 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 451 test: 452 suffix: d 453 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0 454 test: 455 suffix: e 456 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 457 test: 458 suffix: f 459 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 460 461 # load old format, save in new format, reload topology, distribute, load geometry and labels 462 testset: 463 suffix: 3 464 requires: !complex datafilespath 465 args: -dm_plex_name plex 466 args: -dm_plex_view_hdf5_storage_version 2.0.0 467 args: -dm_plex_interpolate -load_dm_distribute 0 468 args: -use_low_level_functions -compare_pre_post 469 args: -num_labels 1 470 args: -outfile ex56_3.h5 471 nsize: 3 472 test: 473 suffix: a 474 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 475 test: 476 suffix: b 477 TODO: broken 478 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 479 test: 480 suffix: c 481 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 482 test: 483 suffix: d 484 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 -num_labels 0 485 test: 486 suffix: e 487 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 488 test: 489 suffix: f 490 args: -dm_plex_filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 491 TEST*/ 492