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 PetscFunctionBeginUser; 225cca0b87SVaclav Hapla options->compare = PETSC_FALSE; 2338883f1bSVaclav Hapla options->compare_labels = PETSC_FALSE; 245cca0b87SVaclav Hapla options->distribute = PETSC_TRUE; 255cca0b87SVaclav Hapla options->interpolate = PETSC_FALSE; 265cca0b87SVaclav Hapla options->filename[0] = '\0'; 27ed5e4e85SVaclav Hapla options->meshname[0] = '\0'; 285cca0b87SVaclav Hapla options->format = PETSC_VIEWER_DEFAULT; 295cca0b87SVaclav Hapla options->second_write_read = PETSC_FALSE; 30806c51deSksagiyam options->use_low_level_functions = PETSC_FALSE; 315cca0b87SVaclav Hapla 32d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL)); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex55.c", options->interpolate, &options->interpolate, NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-filename", "The mesh file", "ex55.c", options->filename, options->filename, sizeof(options->filename), NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-meshname", "The mesh file", "ex55.c", options->meshname, options->meshname, sizeof(options->meshname), NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex55.c", options->second_write_read, &options->second_write_read, NULL)); 419566063dSJacob Faibussowitsch PetscCall(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)); 42d0609cedSBarry Smith PetscOptionsEnd(); 435cca0b87SVaclav Hapla PetscFunctionReturn(0); 445cca0b87SVaclav Hapla }; 455cca0b87SVaclav Hapla 46806c51deSksagiyam static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], const char prefix[], AppCtx user, DM *dm_new) 475cca0b87SVaclav Hapla { 485cca0b87SVaclav Hapla DM dmnew; 49ed5e4e85SVaclav Hapla const char savedName[] = "Mesh"; 50ed5e4e85SVaclav Hapla const char loadedName[] = "Mesh_new"; 515cca0b87SVaclav Hapla PetscViewer v; 525cca0b87SVaclav Hapla 535cca0b87SVaclav Hapla PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(v, user.format)); 569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, savedName)); 57806c51deSksagiyam if (user.use_low_level_functions) { 589566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyView(dm, v)); 599566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesView(dm, v)); 609566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsView(dm, v)); 61806c51deSksagiyam } else { 629566063dSJacob Faibussowitsch PetscCall(DMView(dm, v)); 63806c51deSksagiyam } 645cca0b87SVaclav Hapla 659566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(v, FILE_MODE_READ)); 669566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dmnew)); 679566063dSJacob Faibussowitsch PetscCall(DMSetType(dmnew, DMPLEX)); 689566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dmnew, savedName)); 699566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dmnew, prefix)); 70806c51deSksagiyam if (user.use_low_level_functions) { 71c9ad657eSksagiyam PetscSF sfXC; 72c9ad657eSksagiyam 739566063dSJacob Faibussowitsch PetscCall(DMPlexTopologyLoad(dmnew, v, &sfXC)); 749566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesLoad(dmnew, v, sfXC)); 759566063dSJacob Faibussowitsch PetscCall(DMPlexLabelsLoad(dmnew, v, sfXC)); 769566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfXC)); 77806c51deSksagiyam } else { 789566063dSJacob Faibussowitsch PetscCall(DMLoad(dmnew, v)); 79806c51deSksagiyam } 809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dmnew,loadedName)); 815cca0b87SVaclav Hapla 829566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(v)); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&v)); 845cca0b87SVaclav Hapla *dm_new = dmnew; 855cca0b87SVaclav Hapla PetscFunctionReturn(0); 865cca0b87SVaclav Hapla } 875cca0b87SVaclav Hapla 885cca0b87SVaclav Hapla int main(int argc, char **argv) 895cca0b87SVaclav Hapla { 905cca0b87SVaclav Hapla DM dm, dmnew; 915cca0b87SVaclav Hapla PetscPartitioner part; 925cca0b87SVaclav Hapla AppCtx user; 935cca0b87SVaclav Hapla PetscBool flg; 945cca0b87SVaclav Hapla 95*327415f7SBarry Smith PetscFunctionBeginUser; 969566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 979566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 989566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.meshname, user.interpolate, &dm)); 999566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm,"orig_")); 1009566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 1015cca0b87SVaclav Hapla 1025cca0b87SVaclav Hapla if (user.distribute) { 1035cca0b87SVaclav Hapla DM dmdist; 1045cca0b87SVaclav Hapla 1059566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 1069566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 1079566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist)); 1085cca0b87SVaclav Hapla if (dmdist) { 1099566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1105cca0b87SVaclav Hapla dm = dmdist; 1115cca0b87SVaclav Hapla } 1125cca0b87SVaclav Hapla } 1135cca0b87SVaclav Hapla 1149566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm,NULL)); 1159566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 1169566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 1179566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 1185cca0b87SVaclav Hapla 1199566063dSJacob Faibussowitsch PetscCall(DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew)); 1205cca0b87SVaclav Hapla 1215cca0b87SVaclav Hapla if (user.second_write_read) { 1229566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1235cca0b87SVaclav Hapla dm = dmnew; 1249566063dSJacob Faibussowitsch PetscCall(DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew)); 1255cca0b87SVaclav Hapla } 1265cca0b87SVaclav Hapla 1279566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmnew, NULL, "-dm_view")); 1285cca0b87SVaclav Hapla /* TODO: Is it still true? */ 1295cca0b87SVaclav 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. */ 1305cca0b87SVaclav Hapla 1315cca0b87SVaclav Hapla /* This currently makes sense only for sequential meshes. */ 1325cca0b87SVaclav Hapla if (user.compare) { 1339566063dSJacob Faibussowitsch PetscCall(DMPlexEqual(dmnew, dm, &flg)); 13428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal"); 1359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n")); 13638883f1bSVaclav Hapla } 13738883f1bSVaclav Hapla if (user.compare_labels) { 1389566063dSJacob Faibussowitsch PetscCall(DMCompareLabels(dmnew, dm, NULL, NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DMLabels equal\n")); 1405cca0b87SVaclav Hapla } 1415cca0b87SVaclav Hapla 1429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmnew)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 1465cca0b87SVaclav Hapla } 1475cca0b87SVaclav Hapla 1485cca0b87SVaclav Hapla /*TEST 1495cca0b87SVaclav Hapla build: 1505cca0b87SVaclav Hapla requires: hdf5 1515cca0b87SVaclav Hapla # Idempotence of saving/loading 1525cca0b87SVaclav Hapla # Have to replace Exodus file, which is creating uninterpolated edges 1535cca0b87SVaclav Hapla test: 1545cca0b87SVaclav Hapla suffix: 0 15538883f1bSVaclav Hapla TODO: broken 15638883f1bSVaclav Hapla requires: exodusii 1575cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1585cca0b87SVaclav Hapla args: -format hdf5_petsc -compare 1595cca0b87SVaclav Hapla test: 1605cca0b87SVaclav Hapla suffix: 1 16138883f1bSVaclav Hapla TODO: broken 16238883f1bSVaclav Hapla requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES) 1635cca0b87SVaclav Hapla nsize: 2 1645cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1655cca0b87SVaclav Hapla args: -petscpartitioner_type parmetis 1665cca0b87SVaclav Hapla args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail 1675cca0b87SVaclav Hapla testset: 1685cca0b87SVaclav Hapla requires: exodusii 1695cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 1705cca0b87SVaclav Hapla args: -petscpartitioner_type simple 1715cca0b87SVaclav Hapla args: -dm_view ascii::ascii_info_detail 1725cca0b87SVaclav Hapla args: -new_dm_view ascii::ascii_info_detail 1735cca0b87SVaclav Hapla test: 1745cca0b87SVaclav Hapla suffix: 2 1755cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1765cca0b87SVaclav Hapla args: -format {{default hdf5_petsc}separate output} 1775cca0b87SVaclav Hapla args: -interpolate {{0 1}separate output} 1785cca0b87SVaclav Hapla test: 1795cca0b87SVaclav Hapla suffix: 2a 1805cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1815cca0b87SVaclav Hapla args: -format {{hdf5_xdmf hdf5_viz}separate output} 1825cca0b87SVaclav Hapla test: 1835cca0b87SVaclav Hapla suffix: 3 1845cca0b87SVaclav Hapla requires: exodusii 18538883f1bSVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels 1865cca0b87SVaclav Hapla 1875cca0b87SVaclav Hapla # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2 18838883f1bSVaclav Hapla testset: 1895cca0b87SVaclav Hapla suffix: 4 1905cca0b87SVaclav Hapla requires: !complex 19138883f1bSVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 19238883f1bSVaclav Hapla args: -distribute 0 -second_write_read -compare 19338883f1bSVaclav Hapla test: 19438883f1bSVaclav Hapla suffix: hdf5_petsc 19538883f1bSVaclav Hapla nsize: {{1 2}} 19638883f1bSVaclav Hapla args: -format hdf5_petsc -compare_labels 19738883f1bSVaclav Hapla test: 19838883f1bSVaclav Hapla suffix: hdf5_xdmf 19938883f1bSVaclav Hapla nsize: {{1 3 8}} 20038883f1bSVaclav Hapla args: -format hdf5_xdmf 2015cca0b87SVaclav Hapla 202806c51deSksagiyam # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load() 203806c51deSksagiyam # Output must be the same as ex55_2_nsize-2_format-hdf5_petsc_interpolate-0.out 204806c51deSksagiyam test: 205806c51deSksagiyam suffix: 5 206806c51deSksagiyam requires: exodusii 207806c51deSksagiyam nsize: 2 208806c51deSksagiyam args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 209806c51deSksagiyam args: -petscpartitioner_type simple 210806c51deSksagiyam args: -dm_view ascii::ascii_info_detail 211806c51deSksagiyam args: -new_dm_view ascii::ascii_info_detail 212806c51deSksagiyam args: -format hdf5_petsc -use_low_level_functions 213806c51deSksagiyam 2145cca0b87SVaclav Hapla testset: 21538883f1bSVaclav Hapla suffix: 6 21638883f1bSVaclav Hapla requires: hdf5 !complex datafilespath 21738883f1bSVaclav Hapla nsize: {{1 3}} 21838883f1bSVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 21938883f1bSVaclav 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 22038883f1bSVaclav Hapla args: -format hdf5_petsc -second_write_read -compare -compare_labels 22138883f1bSVaclav Hapla args: -interpolate {{0 1}} -distribute {{0 1}} -petscpartitioner_type simple 22238883f1bSVaclav Hapla 22338883f1bSVaclav Hapla testset: 2245cca0b87SVaclav Hapla # the same data and settings as dm_impls_plex_tests-ex18_9% 22538883f1bSVaclav Hapla suffix: 9 2265cca0b87SVaclav Hapla requires: hdf5 !complex datafilespath 2274d056bb6SVaclav Hapla nsize: {{1 2 4}} 2285cca0b87SVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 2295cca0b87SVaclav 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 2305cca0b87SVaclav Hapla args: -format hdf5_xdmf -second_write_read -compare 2315cca0b87SVaclav Hapla test: 23238883f1bSVaclav Hapla suffix: hdf5_seqload 2335cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type simple 2345cca0b87SVaclav Hapla args: -interpolate {{0 1}} 2355cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 2365cca0b87SVaclav Hapla test: 23738883f1bSVaclav Hapla suffix: hdf5_seqload_metis 2385cca0b87SVaclav Hapla requires: parmetis 2395cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2405cca0b87SVaclav Hapla args: -interpolate 1 2415cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 2425cca0b87SVaclav Hapla test: 24338883f1bSVaclav Hapla suffix: hdf5 2445cca0b87SVaclav Hapla args: -interpolate 1 -petscpartitioner_type simple 2455cca0b87SVaclav Hapla test: 24638883f1bSVaclav Hapla suffix: hdf5_repart 2475cca0b87SVaclav Hapla requires: parmetis 2485cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2495cca0b87SVaclav Hapla args: -interpolate 1 2505cca0b87SVaclav Hapla test: 2515cca0b87SVaclav Hapla TODO: Parallel partitioning of uninterpolated meshes not supported 25238883f1bSVaclav Hapla suffix: hdf5_repart_ppu 2535cca0b87SVaclav Hapla requires: parmetis 2545cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2555cca0b87SVaclav Hapla args: -interpolate 0 2565cca0b87SVaclav Hapla 2575cca0b87SVaclav Hapla # reproduce PetscSFView() crash - fixed, left as regression test 2585cca0b87SVaclav Hapla test: 2595cca0b87SVaclav Hapla suffix: new_dm_view 2605cca0b87SVaclav Hapla requires: exodusii 2615cca0b87SVaclav Hapla nsize: 2 2625cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail 26396b7d01aSVaclav Hapla 26496b7d01aSVaclav Hapla # test backward compatibility of petsc_hdf5 format 26596b7d01aSVaclav Hapla testset: 26696b7d01aSVaclav Hapla suffix: 10-v3.16.0-v1.0.0 26796b7d01aSVaclav Hapla requires: hdf5 !complex datafilespath 268e6ef0f20SVaclav Hapla args: -dm_plex_check_all -compare -compare_labels 26996b7d01aSVaclav Hapla args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}} -use_low_level_functions {{0 1}} 27096b7d01aSVaclav Hapla test: 27196b7d01aSVaclav Hapla suffix: a 27296b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5 27396b7d01aSVaclav Hapla test: 27496b7d01aSVaclav Hapla suffix: b 27596b7d01aSVaclav Hapla TODO: broken 27696b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5 27796b7d01aSVaclav Hapla test: 27896b7d01aSVaclav Hapla suffix: c 27996b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5 28096b7d01aSVaclav Hapla test: 28196b7d01aSVaclav Hapla suffix: d 28296b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5 28396b7d01aSVaclav Hapla test: 28496b7d01aSVaclav Hapla suffix: e 28596b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5 28696b7d01aSVaclav Hapla test: 28796b7d01aSVaclav Hapla suffix: f 28896b7d01aSVaclav Hapla args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5 28996b7d01aSVaclav Hapla 2905cca0b87SVaclav Hapla TEST*/ 291