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() */ 95cca0b87SVaclav Hapla PetscBool distribute; /* Distribute the mesh */ 105cca0b87SVaclav Hapla PetscBool interpolate; /* Generate intermediate mesh elements */ 115cca0b87SVaclav Hapla char filename[PETSC_MAX_PATH_LEN]; /* Mesh filename */ 125cca0b87SVaclav Hapla PetscViewerFormat format; /* Format to write and read */ 135cca0b87SVaclav Hapla PetscBool second_write_read; /* Write and read for the 2nd time */ 145cca0b87SVaclav Hapla } AppCtx; 155cca0b87SVaclav Hapla 165cca0b87SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 175cca0b87SVaclav Hapla { 185cca0b87SVaclav Hapla PetscErrorCode ierr; 195cca0b87SVaclav Hapla 205cca0b87SVaclav Hapla PetscFunctionBeginUser; 215cca0b87SVaclav Hapla options->compare = PETSC_FALSE; 225cca0b87SVaclav Hapla options->distribute = PETSC_TRUE; 235cca0b87SVaclav Hapla options->interpolate = PETSC_FALSE; 245cca0b87SVaclav Hapla options->filename[0] = '\0'; 255cca0b87SVaclav Hapla options->format = PETSC_VIEWER_DEFAULT; 265cca0b87SVaclav Hapla options->second_write_read = PETSC_FALSE; 275cca0b87SVaclav Hapla 285cca0b87SVaclav Hapla ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 295cca0b87SVaclav Hapla ierr = PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex5.c", options->compare, &options->compare, NULL);CHKERRQ(ierr); 305cca0b87SVaclav Hapla ierr = PetscOptionsBool("-distribute", "Distribute the mesh", "ex5.c", options->distribute, &options->distribute, NULL);CHKERRQ(ierr); 315cca0b87SVaclav Hapla ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex5.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr); 325cca0b87SVaclav Hapla ierr = PetscOptionsString("-filename", "The mesh file", "ex5.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 335cca0b87SVaclav Hapla ierr = PetscOptionsEnum("-format", "Format to write and read", "ex5.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL);CHKERRQ(ierr); 345cca0b87SVaclav Hapla ierr = PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex5.c", options->second_write_read, &options->second_write_read, NULL);CHKERRQ(ierr); 355cca0b87SVaclav Hapla ierr = PetscOptionsEnd(); 365cca0b87SVaclav Hapla PetscFunctionReturn(0); 375cca0b87SVaclav Hapla }; 385cca0b87SVaclav Hapla 395cca0b87SVaclav Hapla static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], PetscViewerFormat format, const char prefix[], DM *dm_new) 405cca0b87SVaclav Hapla { 415cca0b87SVaclav Hapla DM dmnew; 425cca0b87SVaclav Hapla PetscViewer v; 435cca0b87SVaclav Hapla PetscErrorCode ierr; 445cca0b87SVaclav Hapla 455cca0b87SVaclav Hapla PetscFunctionBeginUser; 465cca0b87SVaclav Hapla ierr = PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v);CHKERRQ(ierr); 475cca0b87SVaclav Hapla ierr = PetscViewerPushFormat(v, format);CHKERRQ(ierr); 485cca0b87SVaclav Hapla ierr = DMView(dm, v);CHKERRQ(ierr); 495cca0b87SVaclav Hapla 505cca0b87SVaclav Hapla ierr = PetscViewerFileSetMode(v, FILE_MODE_READ);CHKERRQ(ierr); 515cca0b87SVaclav Hapla ierr = DMCreate(PETSC_COMM_WORLD, &dmnew);CHKERRQ(ierr); 525cca0b87SVaclav Hapla ierr = DMSetType(dmnew, DMPLEX);CHKERRQ(ierr); 535cca0b87SVaclav Hapla ierr = DMSetOptionsPrefix(dmnew, prefix);CHKERRQ(ierr); 545cca0b87SVaclav Hapla ierr = DMLoad(dmnew, v);CHKERRQ(ierr); 555cca0b87SVaclav Hapla ierr = PetscObjectSetName((PetscObject)dmnew,"Mesh_new");CHKERRQ(ierr); 565cca0b87SVaclav Hapla 575cca0b87SVaclav Hapla ierr = PetscViewerPopFormat(v);CHKERRQ(ierr); 585cca0b87SVaclav Hapla ierr = PetscViewerDestroy(&v);CHKERRQ(ierr); 595cca0b87SVaclav Hapla *dm_new = dmnew; 605cca0b87SVaclav Hapla PetscFunctionReturn(0); 615cca0b87SVaclav Hapla } 625cca0b87SVaclav Hapla 635cca0b87SVaclav Hapla int main(int argc, char **argv) 645cca0b87SVaclav Hapla { 655cca0b87SVaclav Hapla DM dm, dmnew; 665cca0b87SVaclav Hapla PetscPartitioner part; 675cca0b87SVaclav Hapla AppCtx user; 685cca0b87SVaclav Hapla PetscBool flg; 695cca0b87SVaclav Hapla PetscErrorCode ierr; 705cca0b87SVaclav Hapla 715cca0b87SVaclav Hapla ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 725cca0b87SVaclav Hapla ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 735cca0b87SVaclav Hapla ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.interpolate, &dm);CHKERRQ(ierr); 745cca0b87SVaclav Hapla ierr = DMSetOptionsPrefix(dm,"orig_");CHKERRQ(ierr); 755cca0b87SVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 765cca0b87SVaclav Hapla 775cca0b87SVaclav Hapla if (user.distribute) { 785cca0b87SVaclav Hapla DM dmdist; 795cca0b87SVaclav Hapla 805cca0b87SVaclav Hapla ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 815cca0b87SVaclav Hapla ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 825cca0b87SVaclav Hapla ierr = DMPlexDistribute(dm, 0, NULL, &dmdist);CHKERRQ(ierr); 835cca0b87SVaclav Hapla if (dmdist) { 845cca0b87SVaclav Hapla ierr = DMDestroy(&dm);CHKERRQ(ierr); 855cca0b87SVaclav Hapla dm = dmdist; 865cca0b87SVaclav Hapla } 875cca0b87SVaclav Hapla } 885cca0b87SVaclav Hapla 895cca0b87SVaclav Hapla ierr = DMSetOptionsPrefix(dm,NULL);CHKERRQ(ierr); 905cca0b87SVaclav Hapla ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 915cca0b87SVaclav Hapla ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 925cca0b87SVaclav Hapla 935cca0b87SVaclav Hapla ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", user.format, "new_", &dmnew);CHKERRQ(ierr); 945cca0b87SVaclav Hapla 955cca0b87SVaclav Hapla if (user.second_write_read) { 965cca0b87SVaclav Hapla ierr = DMDestroy(&dm);CHKERRQ(ierr); 975cca0b87SVaclav Hapla dm = dmnew; 985cca0b87SVaclav Hapla ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", user.format, "new_", &dmnew);CHKERRQ(ierr); 995cca0b87SVaclav Hapla } 1005cca0b87SVaclav Hapla 1015cca0b87SVaclav Hapla ierr = DMViewFromOptions(dmnew, NULL, "-dm_view");CHKERRQ(ierr); 1025cca0b87SVaclav Hapla /* TODO: Is it still true? */ 1035cca0b87SVaclav 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. */ 1045cca0b87SVaclav Hapla 1055cca0b87SVaclav Hapla /* This currently makes sense only for sequential meshes. */ 1065cca0b87SVaclav Hapla if (user.compare) { 1075cca0b87SVaclav Hapla ierr = DMPlexEqual(dmnew, dm, &flg);CHKERRQ(ierr); 1085cca0b87SVaclav Hapla if (flg) {ierr = PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");CHKERRQ(ierr);} 1095cca0b87SVaclav Hapla else {ierr = PetscPrintf(PETSC_COMM_WORLD,"DMs are not equal\n");CHKERRQ(ierr);} 1105cca0b87SVaclav Hapla } 1115cca0b87SVaclav Hapla 1125cca0b87SVaclav Hapla ierr = DMDestroy(&dm);CHKERRQ(ierr); 1135cca0b87SVaclav Hapla ierr = DMDestroy(&dmnew);CHKERRQ(ierr); 1145cca0b87SVaclav Hapla ierr = PetscFinalize(); 1155cca0b87SVaclav Hapla return ierr; 1165cca0b87SVaclav Hapla } 1175cca0b87SVaclav Hapla 1185cca0b87SVaclav Hapla /*TEST 1195cca0b87SVaclav Hapla build: 1205cca0b87SVaclav Hapla requires: hdf5 1215cca0b87SVaclav Hapla # Idempotence of saving/loading 1225cca0b87SVaclav Hapla # Have to replace Exodus file, which is creating uninterpolated edges 1235cca0b87SVaclav Hapla test: 1245cca0b87SVaclav Hapla suffix: 0 1255cca0b87SVaclav Hapla requires: exodusii broken 1265cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1275cca0b87SVaclav Hapla args: -format hdf5_petsc -compare 1285cca0b87SVaclav Hapla test: 1295cca0b87SVaclav Hapla suffix: 1 1305cca0b87SVaclav Hapla requires: exodusii parmetis !define(PETSC_USE_64BIT_INDICES) broken 1315cca0b87SVaclav Hapla nsize: 2 1325cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail 1335cca0b87SVaclav Hapla args: -petscpartitioner_type parmetis 1345cca0b87SVaclav Hapla args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail 1355cca0b87SVaclav Hapla testset: 1365cca0b87SVaclav Hapla requires: exodusii 1375cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 1385cca0b87SVaclav Hapla args: -petscpartitioner_type simple 1395cca0b87SVaclav Hapla args: -dm_view ascii::ascii_info_detail 1405cca0b87SVaclav Hapla args: -new_dm_view ascii::ascii_info_detail 1415cca0b87SVaclav Hapla test: 1425cca0b87SVaclav Hapla suffix: 2 1435cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1445cca0b87SVaclav Hapla args: -format {{default hdf5_petsc}separate output} 1455cca0b87SVaclav Hapla args: -interpolate {{0 1}separate output} 1465cca0b87SVaclav Hapla test: 1475cca0b87SVaclav Hapla suffix: 2a 1485cca0b87SVaclav Hapla nsize: {{1 2 4 8}separate output} 1495cca0b87SVaclav Hapla args: -format {{hdf5_xdmf hdf5_viz}separate output} 1505cca0b87SVaclav Hapla test: 1515cca0b87SVaclav Hapla suffix: 3 1525cca0b87SVaclav Hapla requires: exodusii 1535cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare 1545cca0b87SVaclav Hapla 1555cca0b87SVaclav Hapla # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2 1565cca0b87SVaclav Hapla test: 1575cca0b87SVaclav Hapla suffix: 4 1585cca0b87SVaclav Hapla requires: !complex 1595cca0b87SVaclav Hapla nsize: {{1 2 3 4 8}} 1605cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 1615cca0b87SVaclav Hapla args: -dm_plex_create_from_hdf5_xdmf -distribute 0 -format hdf5_xdmf -second_write_read -compare 1625cca0b87SVaclav Hapla 1635cca0b87SVaclav Hapla testset: 1645cca0b87SVaclav Hapla # the same data and settings as dm_impls_plex_tests-ex18_9% 1655cca0b87SVaclav Hapla requires: hdf5 !complex datafilespath 166*4d056bb6SVaclav Hapla nsize: {{1 2 4}} 1675cca0b87SVaclav Hapla args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry 1685cca0b87SVaclav 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 1695cca0b87SVaclav Hapla args: -format hdf5_xdmf -second_write_read -compare 1705cca0b87SVaclav Hapla test: 1715cca0b87SVaclav Hapla suffix: 9_hdf5_seqload 1725cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type simple 1735cca0b87SVaclav Hapla args: -interpolate {{0 1}} 1745cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 1755cca0b87SVaclav Hapla test: 1765cca0b87SVaclav Hapla suffix: 9_hdf5_seqload_metis 1775cca0b87SVaclav Hapla requires: parmetis 1785cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 1795cca0b87SVaclav Hapla args: -interpolate 1 1805cca0b87SVaclav Hapla args: -dm_plex_hdf5_force_sequential 1815cca0b87SVaclav Hapla test: 1825cca0b87SVaclav Hapla suffix: 9_hdf5 1835cca0b87SVaclav Hapla args: -interpolate 1 -petscpartitioner_type simple 1845cca0b87SVaclav Hapla test: 1855cca0b87SVaclav Hapla suffix: 9_hdf5_repart 1865cca0b87SVaclav Hapla requires: parmetis 1875cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 1885cca0b87SVaclav Hapla args: -interpolate 1 1895cca0b87SVaclav Hapla test: 1905cca0b87SVaclav Hapla TODO: Parallel partitioning of uninterpolated meshes not supported 1915cca0b87SVaclav Hapla suffix: 9_hdf5_repart_ppu 1925cca0b87SVaclav Hapla requires: parmetis 1935cca0b87SVaclav Hapla args: -distribute -petscpartitioner_type parmetis 1945cca0b87SVaclav Hapla args: -interpolate 0 1955cca0b87SVaclav Hapla 1965cca0b87SVaclav Hapla # reproduce PetscSFView() crash - fixed, left as regression test 1975cca0b87SVaclav Hapla test: 1985cca0b87SVaclav Hapla suffix: new_dm_view 1995cca0b87SVaclav Hapla requires: exodusii 2005cca0b87SVaclav Hapla nsize: 2 2015cca0b87SVaclav Hapla args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail 2025cca0b87SVaclav Hapla TEST*/ 203