15cca0b87SVaclav Hapla static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n"; 25cca0b87SVaclav Hapla 35cca0b87SVaclav Hapla #include <petscdmplex.h> 45cca0b87SVaclav Hapla #include <petscviewerhdf5.h> 55cca0b87SVaclav Hapla #include <petscsf.h> 65cca0b87SVaclav Hapla 75cca0b87SVaclav Hapla typedef struct { 85cca0b87SVaclav Hapla PetscBool compare; /* Compare the meshes using DMPlexEqual() */ 938883f1bSVaclav Hapla PetscBool compare_labels; /* Compare labels in the meshes using DMCompareLabels() */ 105cca0b87SVaclav Hapla PetscBool distribute; /* Distribute the mesh */ 115cca0b87SVaclav Hapla PetscBool interpolate; /* Generate intermediate mesh elements */ 125cca0b87SVaclav Hapla char filename[PETSC_MAX_PATH_LEN]; /* Mesh filename */ 13ed5e4e85SVaclav Hapla char meshname[PETSC_MAX_PATH_LEN]; /* Mesh name */ 145cca0b87SVaclav Hapla PetscViewerFormat format; /* Format to write and read */ 155cca0b87SVaclav Hapla PetscBool second_write_read; /* Write and read for the 2nd time */ 16806c51deSksagiyam PetscBool use_low_level_functions; /* Use low level functions for viewing and loading */ 175cca0b87SVaclav Hapla } AppCtx; 185cca0b87SVaclav Hapla 195cca0b87SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 205cca0b87SVaclav Hapla { 215cca0b87SVaclav Hapla PetscErrorCode ierr; 225cca0b87SVaclav Hapla 235cca0b87SVaclav Hapla PetscFunctionBeginUser; 245cca0b87SVaclav Hapla options->compare = PETSC_FALSE; 2538883f1bSVaclav Hapla options->compare_labels = PETSC_FALSE; 265cca0b87SVaclav Hapla options->distribute = PETSC_TRUE; 275cca0b87SVaclav Hapla options->interpolate = PETSC_FALSE; 285cca0b87SVaclav Hapla options->filename[0] = '\0'; 29ed5e4e85SVaclav Hapla options->meshname[0] = '\0'; 305cca0b87SVaclav Hapla options->format = PETSC_VIEWER_DEFAULT; 315cca0b87SVaclav Hapla options->second_write_read = PETSC_FALSE; 32806c51deSksagiyam options->use_low_level_functions = PETSC_FALSE; 335cca0b87SVaclav Hapla 345cca0b87SVaclav Hapla ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex55.c", options->interpolate, &options->interpolate, NULL)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-filename", "The mesh file", "ex55.c", options->filename, options->filename, sizeof(options->filename), NULL)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-meshname", "The mesh file", "ex55.c", options->meshname, options->meshname, sizeof(options->meshname), NULL)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex55.c", options->second_write_read, &options->second_write_read, NULL)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-use_low_level_functions", "Use low level functions for viewing and loading", "ex55.c", options->use_low_level_functions, &options->use_low_level_functions, NULL)); 441e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 455cca0b87SVaclav Hapla PetscFunctionReturn(0); 465cca0b87SVaclav Hapla }; 475cca0b87SVaclav Hapla 48806c51deSksagiyam static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], const char prefix[], AppCtx user, DM *dm_new) 495cca0b87SVaclav Hapla { 505cca0b87SVaclav Hapla DM dmnew; 51ed5e4e85SVaclav Hapla const char savedName[] = "Mesh"; 52ed5e4e85SVaclav Hapla const char loadedName[] = "Mesh_new"; 535cca0b87SVaclav Hapla PetscViewer v; 545cca0b87SVaclav Hapla 555cca0b87SVaclav Hapla PetscFunctionBeginUser; 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(v, user.format)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, savedName)); 59806c51deSksagiyam if (user.use_low_level_functions) { 605f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTopologyView(dm, v)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCoordinatesView(dm, v)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexLabelsView(dm, v)); 63806c51deSksagiyam } else { 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(dm, v)); 65806c51deSksagiyam } 665cca0b87SVaclav Hapla 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(v, FILE_MODE_READ)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dmnew)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dmnew, DMPLEX)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dmnew, savedName)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(dmnew, prefix)); 72806c51deSksagiyam if (user.use_low_level_functions) { 73c9ad657eSksagiyam PetscSF sfXC; 74c9ad657eSksagiyam 755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTopologyLoad(dmnew, v, &sfXC)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCoordinatesLoad(dmnew, v, sfXC)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexLabelsLoad(dmnew, v, sfXC)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sfXC)); 79806c51deSksagiyam } else { 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMLoad(dmnew, v)); 81806c51deSksagiyam } 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)dmnew,loadedName)); 835cca0b87SVaclav Hapla 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(v)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&v)); 865cca0b87SVaclav Hapla *dm_new = dmnew; 875cca0b87SVaclav Hapla PetscFunctionReturn(0); 885cca0b87SVaclav Hapla } 895cca0b87SVaclav Hapla 905cca0b87SVaclav Hapla int main(int argc, char **argv) 915cca0b87SVaclav Hapla { 925cca0b87SVaclav Hapla DM dm, dmnew; 935cca0b87SVaclav Hapla PetscPartitioner part; 945cca0b87SVaclav Hapla AppCtx user; 955cca0b87SVaclav Hapla PetscBool flg; 965cca0b87SVaclav Hapla PetscErrorCode ierr; 975cca0b87SVaclav Hapla 985cca0b87SVaclav Hapla ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 995f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.meshname, user.interpolate, &dm)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(dm,"orig_")); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 1035cca0b87SVaclav Hapla 1045cca0b87SVaclav Hapla if (user.distribute) { 1055cca0b87SVaclav Hapla DM dmdist; 1065cca0b87SVaclav Hapla 1075f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm, &part)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm, 0, NULL, &dmdist)); 1105cca0b87SVaclav Hapla if (dmdist) { 1115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1125cca0b87SVaclav Hapla dm = dmdist; 1135cca0b87SVaclav Hapla } 1145cca0b87SVaclav Hapla } 1155cca0b87SVaclav Hapla 1165f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOptionsPrefix(dm,NULL)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 1205cca0b87SVaclav Hapla 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew)); 1225cca0b87SVaclav Hapla 1235cca0b87SVaclav Hapla if (user.second_write_read) { 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1255cca0b87SVaclav Hapla dm = dmnew; 1265f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew)); 1275cca0b87SVaclav Hapla } 1285cca0b87SVaclav Hapla 1295f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dmnew, NULL, "-dm_view")); 1305cca0b87SVaclav Hapla /* TODO: Is it still true? */ 1315cca0b87SVaclav Hapla /* The NATIVE format for coordiante viewing is killing parallel output, since we have a local vector. Map it to global, and it will work. */ 1325cca0b87SVaclav Hapla 1335cca0b87SVaclav Hapla /* This currently makes sense only for sequential meshes. */ 1345cca0b87SVaclav Hapla if (user.compare) { 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexEqual(dmnew, dm, &flg)); 136*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal"); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n")); 13838883f1bSVaclav Hapla } 13938883f1bSVaclav Hapla if (user.compare_labels) { 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompareLabels(dmnew, dm, NULL, NULL)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"DMLabels equal\n")); 1425cca0b87SVaclav Hapla } 1435cca0b87SVaclav Hapla 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmnew)); 1465cca0b87SVaclav Hapla ierr = PetscFinalize(); 1475cca0b87SVaclav Hapla return ierr; 1485cca0b87SVaclav Hapla } 1495cca0b87SVaclav Hapla 1505cca0b87SVaclav Hapla /*TEST 1515cca0b87SVaclav Hapla build: 1525cca0b87SVaclav Hapla requires: hdf5 1535cca0b87SVaclav Hapla # Idempotence of saving/loading 1545cca0b87SVaclav Hapla # Have to replace Exodus file, which is creating uninterpolated edges 1555cca0b87SVaclav Hapla test: 1565cca0b87SVaclav Hapla suffix: 0 15738883f1bSVaclav Hapla TODO: broken 15838883f1bSVaclav Hapla requires: exodusii 1595cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1605cca0b87SVaclav Hapla args: -format hdf5_petsc -compare 1615cca0b87SVaclav Hapla test: 1625cca0b87SVaclav Hapla suffix: 1 16338883f1bSVaclav Hapla TODO: broken 16438883f1bSVaclav Hapla requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES) 1655cca0b87SVaclav Hapla nsize: 2 1665cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1675cca0b87SVaclav Hapla args: -petscpartitioner_type parmetis 1685cca0b87SVaclav Hapla args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail 1695cca0b87SVaclav Hapla testset: 1705cca0b87SVaclav Hapla requires: exodusii 1715cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 1725cca0b87SVaclav Hapla args: -petscpartitioner_type simple 1735cca0b87SVaclav Hapla args: -dm_view ascii::ascii_info_detail 1745cca0b87SVaclav Hapla args: -new_dm_view ascii::ascii_info_detail 1755cca0b87SVaclav Hapla test: 1765cca0b87SVaclav Hapla suffix: 2 1775cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1785cca0b87SVaclav Hapla args: -format {{default hdf5_petsc}separate output} 1795cca0b87SVaclav Hapla args: -interpolate {{0 1}separate output} 1805cca0b87SVaclav Hapla test: 1815cca0b87SVaclav Hapla suffix: 2a 1825cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1835cca0b87SVaclav Hapla args: -format {{hdf5_xdmf hdf5_viz}separate output} 1845cca0b87SVaclav Hapla test: 1855cca0b87SVaclav Hapla suffix: 3 1865cca0b87SVaclav Hapla requires: exodusii 18738883f1bSVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels 1885cca0b87SVaclav Hapla 1895cca0b87SVaclav Hapla # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2 19038883f1bSVaclav Hapla testset: 1915cca0b87SVaclav Hapla suffix: 4 1925cca0b87SVaclav Hapla requires: !complex 19338883f1bSVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 19438883f1bSVaclav Hapla args: -distribute 0 -second_write_read -compare 19538883f1bSVaclav Hapla test: 19638883f1bSVaclav Hapla suffix: hdf5_petsc 19738883f1bSVaclav Hapla nsize: {{1 2}} 19838883f1bSVaclav Hapla args: -format hdf5_petsc -compare_labels 19938883f1bSVaclav Hapla test: 20038883f1bSVaclav Hapla suffix: hdf5_xdmf 20138883f1bSVaclav Hapla nsize: {{1 3 8}} 20238883f1bSVaclav Hapla args: -format hdf5_xdmf 2035cca0b87SVaclav Hapla 204806c51deSksagiyam # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load() 205806c51deSksagiyam # Output must be the same as ex55_2_nsize-2_format-hdf5_petsc_interpolate-0.out 206806c51deSksagiyam test: 207806c51deSksagiyam suffix: 5 208806c51deSksagiyam requires: exodusii 209806c51deSksagiyam nsize: 2 210806c51deSksagiyam args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 211806c51deSksagiyam args: -petscpartitioner_type simple 212806c51deSksagiyam args: -dm_view ascii::ascii_info_detail 213806c51deSksagiyam args: -new_dm_view ascii::ascii_info_detail 214806c51deSksagiyam args: -format hdf5_petsc -use_low_level_functions 215806c51deSksagiyam 2165cca0b87SVaclav Hapla testset: 21738883f1bSVaclav Hapla suffix: 6 21838883f1bSVaclav Hapla requires: hdf5 !complex datafilespath 21938883f1bSVaclav Hapla nsize: {{1 3}} 22038883f1bSVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 22138883f1bSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates 22238883f1bSVaclav Hapla args: -format hdf5_petsc -second_write_read -compare -compare_labels 22338883f1bSVaclav Hapla args: -interpolate {{0 1}} -distribute {{0 1}} -petscpartitioner_type simple 22438883f1bSVaclav Hapla 22538883f1bSVaclav Hapla testset: 2265cca0b87SVaclav Hapla # the same data and settings as dm_impls_plex_tests-ex18_9% 22738883f1bSVaclav Hapla suffix: 9 2285cca0b87SVaclav Hapla requires: hdf5 !complex datafilespath 2294d056bb6SVaclav Hapla nsize: {{1 2 4}} 2305cca0b87SVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 2315cca0b87SVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates 2325cca0b87SVaclav Hapla args: -format hdf5_xdmf -second_write_read -compare 2335cca0b87SVaclav Hapla test: 23438883f1bSVaclav Hapla suffix: hdf5_seqload 2355cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type simple 2365cca0b87SVaclav Hapla args: -interpolate {{0 1}} 2375cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 2385cca0b87SVaclav Hapla test: 23938883f1bSVaclav Hapla suffix: hdf5_seqload_metis 2405cca0b87SVaclav Hapla requires: parmetis 2415cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2425cca0b87SVaclav Hapla args: -interpolate 1 2435cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 2445cca0b87SVaclav Hapla test: 24538883f1bSVaclav Hapla suffix: hdf5 2465cca0b87SVaclav Hapla args: -interpolate 1 -petscpartitioner_type simple 2475cca0b87SVaclav Hapla test: 24838883f1bSVaclav Hapla suffix: hdf5_repart 2495cca0b87SVaclav Hapla requires: parmetis 2505cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2515cca0b87SVaclav Hapla args: -interpolate 1 2525cca0b87SVaclav Hapla test: 2535cca0b87SVaclav Hapla TODO: Parallel partitioning of uninterpolated meshes not supported 25438883f1bSVaclav Hapla suffix: hdf5_repart_ppu 2555cca0b87SVaclav Hapla requires: parmetis 2565cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2575cca0b87SVaclav Hapla args: -interpolate 0 2585cca0b87SVaclav Hapla 2595cca0b87SVaclav Hapla # reproduce PetscSFView() crash - fixed, left as regression test 2605cca0b87SVaclav Hapla test: 2615cca0b87SVaclav Hapla suffix: new_dm_view 2625cca0b87SVaclav Hapla requires: exodusii 2635cca0b87SVaclav Hapla nsize: 2 2645cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail 26596b7d01aSVaclav Hapla 26696b7d01aSVaclav Hapla # test backward compatibility of petsc_hdf5 format 26796b7d01aSVaclav Hapla testset: 26896b7d01aSVaclav Hapla suffix: 10-v3.16.0-v1.0.0 26996b7d01aSVaclav Hapla requires: hdf5 !complex datafilespath 270e6ef0f20SVaclav Hapla args: -dm_plex_check_all -compare -compare_labels 27196b7d01aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}} -use_low_level_functions {{0 1}} 27296b7d01aSVaclav Hapla test: 27396b7d01aSVaclav Hapla suffix: a 27496b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 27596b7d01aSVaclav Hapla test: 27696b7d01aSVaclav Hapla suffix: b 27796b7d01aSVaclav Hapla TODO: broken 27896b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 27996b7d01aSVaclav Hapla test: 28096b7d01aSVaclav Hapla suffix: c 28196b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 28296b7d01aSVaclav Hapla test: 28396b7d01aSVaclav Hapla suffix: d 28496b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 28596b7d01aSVaclav Hapla test: 28696b7d01aSVaclav Hapla suffix: e 28796b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 28896b7d01aSVaclav Hapla test: 28996b7d01aSVaclav Hapla suffix: f 29096b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 29196b7d01aSVaclav Hapla 2925cca0b87SVaclav Hapla TEST*/ 293