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 */ 135cca0b87SVaclav Hapla PetscViewerFormat format; /* Format to write and read */ 145cca0b87SVaclav Hapla PetscBool second_write_read; /* Write and read for the 2nd time */ 15806c51deSksagiyam PetscBool use_low_level_functions; /* Use low level functions for viewing and loading */ 165cca0b87SVaclav Hapla } AppCtx; 175cca0b87SVaclav Hapla 185cca0b87SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 195cca0b87SVaclav Hapla { 205cca0b87SVaclav Hapla PetscErrorCode ierr; 215cca0b87SVaclav Hapla 225cca0b87SVaclav Hapla PetscFunctionBeginUser; 235cca0b87SVaclav Hapla options->compare = PETSC_FALSE; 2438883f1bSVaclav Hapla options->compare_labels = PETSC_FALSE; 255cca0b87SVaclav Hapla options->distribute = PETSC_TRUE; 265cca0b87SVaclav Hapla options->interpolate = PETSC_FALSE; 275cca0b87SVaclav Hapla options->filename[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 325cca0b87SVaclav Hapla ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 33806c51deSksagiyam ierr = PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL);CHKERRQ(ierr); 3438883f1bSVaclav Hapla ierr = PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL);CHKERRQ(ierr); 35806c51deSksagiyam ierr = PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL);CHKERRQ(ierr); 36806c51deSksagiyam ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex55.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr); 37806c51deSksagiyam ierr = PetscOptionsString("-filename", "The mesh file", "ex55.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 38806c51deSksagiyam ierr = PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL);CHKERRQ(ierr); 39806c51deSksagiyam ierr = PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex55.c", options->second_write_read, &options->second_write_read, NULL);CHKERRQ(ierr); 40806c51deSksagiyam ierr = 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);CHKERRQ(ierr); 411e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 425cca0b87SVaclav Hapla PetscFunctionReturn(0); 435cca0b87SVaclav Hapla }; 445cca0b87SVaclav Hapla 45806c51deSksagiyam static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], const char prefix[], AppCtx user, DM *dm_new) 465cca0b87SVaclav Hapla { 475cca0b87SVaclav Hapla DM dmnew; 485cca0b87SVaclav Hapla PetscViewer v; 495cca0b87SVaclav Hapla PetscErrorCode ierr; 505cca0b87SVaclav Hapla 515cca0b87SVaclav Hapla PetscFunctionBeginUser; 525cca0b87SVaclav Hapla ierr = PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v);CHKERRQ(ierr); 53806c51deSksagiyam ierr = PetscViewerPushFormat(v, user.format);CHKERRQ(ierr); 54806c51deSksagiyam if (user.use_low_level_functions) { 55806c51deSksagiyam ierr = DMPlexTopologyView(dm, v);CHKERRQ(ierr); 56806c51deSksagiyam ierr = DMPlexCoordinatesView(dm, v);CHKERRQ(ierr); 57806c51deSksagiyam ierr = DMPlexLabelsView(dm, v);CHKERRQ(ierr); 58806c51deSksagiyam } else { 595cca0b87SVaclav Hapla ierr = DMView(dm, v);CHKERRQ(ierr); 60806c51deSksagiyam } 615cca0b87SVaclav Hapla 625cca0b87SVaclav Hapla ierr = PetscViewerFileSetMode(v, FILE_MODE_READ);CHKERRQ(ierr); 635cca0b87SVaclav Hapla ierr = DMCreate(PETSC_COMM_WORLD, &dmnew);CHKERRQ(ierr); 645cca0b87SVaclav Hapla ierr = DMSetType(dmnew, DMPLEX);CHKERRQ(ierr); 655cca0b87SVaclav Hapla ierr = DMSetOptionsPrefix(dmnew, prefix);CHKERRQ(ierr); 66806c51deSksagiyam if (user.use_low_level_functions) { 67806c51deSksagiyam ierr = DMPlexTopologyLoad(dmnew, v, NULL);CHKERRQ(ierr); 68806c51deSksagiyam ierr = DMPlexCoordinatesLoad(dmnew, v);CHKERRQ(ierr); 69806c51deSksagiyam ierr = DMPlexLabelsLoad(dmnew, v);CHKERRQ(ierr); 70806c51deSksagiyam } else { 715cca0b87SVaclav Hapla ierr = DMLoad(dmnew, v);CHKERRQ(ierr); 72806c51deSksagiyam } 735cca0b87SVaclav Hapla ierr = PetscObjectSetName((PetscObject)dmnew,"Mesh_new");CHKERRQ(ierr); 745cca0b87SVaclav Hapla 755cca0b87SVaclav Hapla ierr = PetscViewerPopFormat(v);CHKERRQ(ierr); 765cca0b87SVaclav Hapla ierr = PetscViewerDestroy(&v);CHKERRQ(ierr); 775cca0b87SVaclav Hapla *dm_new = dmnew; 785cca0b87SVaclav Hapla PetscFunctionReturn(0); 795cca0b87SVaclav Hapla } 805cca0b87SVaclav Hapla 815cca0b87SVaclav Hapla int main(int argc, char **argv) 825cca0b87SVaclav Hapla { 835cca0b87SVaclav Hapla DM dm, dmnew; 845cca0b87SVaclav Hapla PetscPartitioner part; 855cca0b87SVaclav Hapla AppCtx user; 865cca0b87SVaclav Hapla PetscBool flg; 875cca0b87SVaclav Hapla PetscErrorCode ierr; 885cca0b87SVaclav Hapla 895cca0b87SVaclav Hapla ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 905cca0b87SVaclav Hapla ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 91*cd7e8a5eSksagiyam ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, "ex55_plex", user.interpolate, &dm);CHKERRQ(ierr); 925cca0b87SVaclav Hapla ierr = DMSetOptionsPrefix(dm,"orig_");CHKERRQ(ierr); 935cca0b87SVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 945cca0b87SVaclav Hapla 955cca0b87SVaclav Hapla if (user.distribute) { 965cca0b87SVaclav Hapla DM dmdist; 975cca0b87SVaclav Hapla 985cca0b87SVaclav Hapla ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 995cca0b87SVaclav Hapla ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 1005cca0b87SVaclav Hapla ierr = DMPlexDistribute(dm, 0, NULL, &dmdist);CHKERRQ(ierr); 1015cca0b87SVaclav Hapla if (dmdist) { 1025cca0b87SVaclav Hapla ierr = DMDestroy(&dm);CHKERRQ(ierr); 1035cca0b87SVaclav Hapla dm = dmdist; 1045cca0b87SVaclav Hapla } 1055cca0b87SVaclav Hapla } 1065cca0b87SVaclav Hapla 1075cca0b87SVaclav Hapla ierr = DMSetOptionsPrefix(dm,NULL);CHKERRQ(ierr); 1085cca0b87SVaclav Hapla ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 1095cca0b87SVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 1105cca0b87SVaclav Hapla 111806c51deSksagiyam ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew);CHKERRQ(ierr); 1125cca0b87SVaclav Hapla 1135cca0b87SVaclav Hapla if (user.second_write_read) { 1145cca0b87SVaclav Hapla ierr = DMDestroy(&dm);CHKERRQ(ierr); 1155cca0b87SVaclav Hapla dm = dmnew; 116806c51deSksagiyam ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew);CHKERRQ(ierr); 1175cca0b87SVaclav Hapla } 1185cca0b87SVaclav Hapla 1195cca0b87SVaclav Hapla ierr = DMViewFromOptions(dmnew, NULL, "-dm_view");CHKERRQ(ierr); 1205cca0b87SVaclav Hapla /* TODO: Is it still true? */ 1215cca0b87SVaclav 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. */ 1225cca0b87SVaclav Hapla 1235cca0b87SVaclav Hapla /* This currently makes sense only for sequential meshes. */ 1245cca0b87SVaclav Hapla if (user.compare) { 1255cca0b87SVaclav Hapla ierr = DMPlexEqual(dmnew, dm, &flg);CHKERRQ(ierr); 12638883f1bSVaclav Hapla if (!flg) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal"); 12738883f1bSVaclav Hapla ierr = PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");CHKERRQ(ierr); 12838883f1bSVaclav Hapla } 12938883f1bSVaclav Hapla if (user.compare_labels) { 13038883f1bSVaclav Hapla ierr = DMCompareLabels(dmnew, dm, NULL, NULL);CHKERRQ(ierr); 13138883f1bSVaclav Hapla ierr = PetscPrintf(PETSC_COMM_WORLD,"DMLabels equal\n");CHKERRQ(ierr); 1325cca0b87SVaclav Hapla } 1335cca0b87SVaclav Hapla 1345cca0b87SVaclav Hapla ierr = DMDestroy(&dm);CHKERRQ(ierr); 1355cca0b87SVaclav Hapla ierr = DMDestroy(&dmnew);CHKERRQ(ierr); 1365cca0b87SVaclav Hapla ierr = PetscFinalize(); 1375cca0b87SVaclav Hapla return ierr; 1385cca0b87SVaclav Hapla } 1395cca0b87SVaclav Hapla 1405cca0b87SVaclav Hapla /*TEST 1415cca0b87SVaclav Hapla build: 1425cca0b87SVaclav Hapla requires: hdf5 1435cca0b87SVaclav Hapla # Idempotence of saving/loading 1445cca0b87SVaclav Hapla # Have to replace Exodus file, which is creating uninterpolated edges 1455cca0b87SVaclav Hapla test: 1465cca0b87SVaclav Hapla suffix: 0 14738883f1bSVaclav Hapla TODO: broken 14838883f1bSVaclav Hapla requires: exodusii 1495cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1505cca0b87SVaclav Hapla args: -format hdf5_petsc -compare 1515cca0b87SVaclav Hapla test: 1525cca0b87SVaclav Hapla suffix: 1 15338883f1bSVaclav Hapla TODO: broken 15438883f1bSVaclav Hapla requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES) 1555cca0b87SVaclav Hapla nsize: 2 1565cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1575cca0b87SVaclav Hapla args: -petscpartitioner_type parmetis 1585cca0b87SVaclav Hapla args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail 1595cca0b87SVaclav Hapla testset: 1605cca0b87SVaclav Hapla requires: exodusii 1615cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 1625cca0b87SVaclav Hapla args: -petscpartitioner_type simple 1635cca0b87SVaclav Hapla args: -dm_view ascii::ascii_info_detail 1645cca0b87SVaclav Hapla args: -new_dm_view ascii::ascii_info_detail 1655cca0b87SVaclav Hapla test: 1665cca0b87SVaclav Hapla suffix: 2 1675cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1685cca0b87SVaclav Hapla args: -format {{default hdf5_petsc}separate output} 1695cca0b87SVaclav Hapla args: -interpolate {{0 1}separate output} 1705cca0b87SVaclav Hapla test: 1715cca0b87SVaclav Hapla suffix: 2a 1725cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1735cca0b87SVaclav Hapla args: -format {{hdf5_xdmf hdf5_viz}separate output} 1745cca0b87SVaclav Hapla test: 1755cca0b87SVaclav Hapla suffix: 3 1765cca0b87SVaclav Hapla requires: exodusii 17738883f1bSVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels 1785cca0b87SVaclav Hapla 1795cca0b87SVaclav Hapla # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2 18038883f1bSVaclav Hapla testset: 1815cca0b87SVaclav Hapla suffix: 4 1825cca0b87SVaclav Hapla requires: !complex 18338883f1bSVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf 18438883f1bSVaclav Hapla args: -distribute 0 -second_write_read -compare 18538883f1bSVaclav Hapla test: 18638883f1bSVaclav Hapla suffix: hdf5_petsc 18738883f1bSVaclav Hapla nsize: {{1 2}} 18838883f1bSVaclav Hapla args: -format hdf5_petsc -compare_labels 18938883f1bSVaclav Hapla test: 19038883f1bSVaclav Hapla suffix: hdf5_xdmf 19138883f1bSVaclav Hapla nsize: {{1 3 8}} 19238883f1bSVaclav Hapla args: -format hdf5_xdmf 1935cca0b87SVaclav Hapla 194806c51deSksagiyam # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load() 195806c51deSksagiyam # Output must be the same as ex55_2_nsize-2_format-hdf5_petsc_interpolate-0.out 196806c51deSksagiyam test: 197806c51deSksagiyam suffix: 5 198806c51deSksagiyam requires: exodusii 199806c51deSksagiyam nsize: 2 200806c51deSksagiyam args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 201806c51deSksagiyam args: -petscpartitioner_type simple 202806c51deSksagiyam args: -dm_view ascii::ascii_info_detail 203806c51deSksagiyam args: -new_dm_view ascii::ascii_info_detail 204806c51deSksagiyam args: -format hdf5_petsc -use_low_level_functions 205806c51deSksagiyam 2065cca0b87SVaclav Hapla testset: 20738883f1bSVaclav Hapla suffix: 6 20838883f1bSVaclav Hapla requires: hdf5 !complex datafilespath 20938883f1bSVaclav Hapla nsize: {{1 3}} 21038883f1bSVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 21138883f1bSVaclav 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 21238883f1bSVaclav Hapla args: -format hdf5_petsc -second_write_read -compare -compare_labels 21338883f1bSVaclav Hapla args: -interpolate {{0 1}} -distribute {{0 1}} -petscpartitioner_type simple 21438883f1bSVaclav Hapla 21538883f1bSVaclav Hapla testset: 2165cca0b87SVaclav Hapla # the same data and settings as dm_impls_plex_tests-ex18_9% 21738883f1bSVaclav Hapla suffix: 9 2185cca0b87SVaclav Hapla requires: hdf5 !complex datafilespath 2194d056bb6SVaclav Hapla nsize: {{1 2 4}} 2205cca0b87SVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 2215cca0b87SVaclav 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 2225cca0b87SVaclav Hapla args: -format hdf5_xdmf -second_write_read -compare 2235cca0b87SVaclav Hapla test: 22438883f1bSVaclav Hapla suffix: hdf5_seqload 2255cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type simple 2265cca0b87SVaclav Hapla args: -interpolate {{0 1}} 2275cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 2285cca0b87SVaclav Hapla test: 22938883f1bSVaclav Hapla suffix: hdf5_seqload_metis 2305cca0b87SVaclav Hapla requires: parmetis 2315cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2325cca0b87SVaclav Hapla args: -interpolate 1 2335cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 2345cca0b87SVaclav Hapla test: 23538883f1bSVaclav Hapla suffix: hdf5 2365cca0b87SVaclav Hapla args: -interpolate 1 -petscpartitioner_type simple 2375cca0b87SVaclav Hapla test: 23838883f1bSVaclav Hapla suffix: hdf5_repart 2395cca0b87SVaclav Hapla requires: parmetis 2405cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2415cca0b87SVaclav Hapla args: -interpolate 1 2425cca0b87SVaclav Hapla test: 2435cca0b87SVaclav Hapla TODO: Parallel partitioning of uninterpolated meshes not supported 24438883f1bSVaclav Hapla suffix: hdf5_repart_ppu 2455cca0b87SVaclav Hapla requires: parmetis 2465cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 2475cca0b87SVaclav Hapla args: -interpolate 0 2485cca0b87SVaclav Hapla 2495cca0b87SVaclav Hapla # reproduce PetscSFView() crash - fixed, left as regression test 2505cca0b87SVaclav Hapla test: 2515cca0b87SVaclav Hapla suffix: new_dm_view 2525cca0b87SVaclav Hapla requires: exodusii 2535cca0b87SVaclav Hapla nsize: 2 2545cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail 2555cca0b87SVaclav Hapla TEST*/ 256