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