xref: /petsc/src/dm/impls/plex/tests/ex55.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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 */
14806c51deSksagiyam   PetscBool use_low_level_functions;      /* Use low level functions for viewing and loading */
155cca0b87SVaclav Hapla } AppCtx;
165cca0b87SVaclav Hapla 
175cca0b87SVaclav Hapla static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
185cca0b87SVaclav Hapla {
195cca0b87SVaclav Hapla   PetscErrorCode ierr;
205cca0b87SVaclav Hapla 
215cca0b87SVaclav Hapla   PetscFunctionBeginUser;
225cca0b87SVaclav Hapla   options->compare = PETSC_FALSE;
235cca0b87SVaclav Hapla   options->distribute = PETSC_TRUE;
245cca0b87SVaclav Hapla   options->interpolate = PETSC_FALSE;
255cca0b87SVaclav Hapla   options->filename[0] = '\0';
265cca0b87SVaclav Hapla   options->format = PETSC_VIEWER_DEFAULT;
275cca0b87SVaclav Hapla   options->second_write_read = PETSC_FALSE;
28806c51deSksagiyam   options->use_low_level_functions = PETSC_FALSE;
295cca0b87SVaclav Hapla 
305cca0b87SVaclav Hapla   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
31806c51deSksagiyam   ierr = PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL);CHKERRQ(ierr);
32806c51deSksagiyam   ierr = PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL);CHKERRQ(ierr);
33806c51deSksagiyam   ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex55.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
34806c51deSksagiyam   ierr = PetscOptionsString("-filename", "The mesh file", "ex55.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
35806c51deSksagiyam   ierr = PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL);CHKERRQ(ierr);
36806c51deSksagiyam   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);
37806c51deSksagiyam   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);
38*1e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
395cca0b87SVaclav Hapla   PetscFunctionReturn(0);
405cca0b87SVaclav Hapla };
415cca0b87SVaclav Hapla 
42806c51deSksagiyam static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], const char prefix[], AppCtx user, DM *dm_new)
435cca0b87SVaclav Hapla {
445cca0b87SVaclav Hapla   DM             dmnew;
455cca0b87SVaclav Hapla   PetscViewer    v;
465cca0b87SVaclav Hapla   PetscErrorCode ierr;
475cca0b87SVaclav Hapla 
485cca0b87SVaclav Hapla   PetscFunctionBeginUser;
495cca0b87SVaclav Hapla   ierr = PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v);CHKERRQ(ierr);
50806c51deSksagiyam   ierr = PetscViewerPushFormat(v, user.format);CHKERRQ(ierr);
51806c51deSksagiyam   if (user.use_low_level_functions) {
52806c51deSksagiyam     ierr = DMPlexTopologyView(dm, v);CHKERRQ(ierr);
53806c51deSksagiyam     ierr = DMPlexCoordinatesView(dm, v);CHKERRQ(ierr);
54806c51deSksagiyam     ierr = DMPlexLabelsView(dm, v);CHKERRQ(ierr);
55806c51deSksagiyam   } else {
565cca0b87SVaclav Hapla     ierr = DMView(dm, v);CHKERRQ(ierr);
57806c51deSksagiyam   }
585cca0b87SVaclav Hapla 
595cca0b87SVaclav Hapla   ierr = PetscViewerFileSetMode(v, FILE_MODE_READ);CHKERRQ(ierr);
605cca0b87SVaclav Hapla   ierr = DMCreate(PETSC_COMM_WORLD, &dmnew);CHKERRQ(ierr);
615cca0b87SVaclav Hapla   ierr = DMSetType(dmnew, DMPLEX);CHKERRQ(ierr);
625cca0b87SVaclav Hapla   ierr = DMSetOptionsPrefix(dmnew, prefix);CHKERRQ(ierr);
63806c51deSksagiyam   if (user.use_low_level_functions) {
64806c51deSksagiyam     ierr = DMPlexTopologyLoad(dmnew, v, NULL);CHKERRQ(ierr);
65806c51deSksagiyam     ierr = DMPlexCoordinatesLoad(dmnew, v);CHKERRQ(ierr);
66806c51deSksagiyam     ierr = DMPlexLabelsLoad(dmnew, v);CHKERRQ(ierr);
67806c51deSksagiyam   } else {
685cca0b87SVaclav Hapla     ierr = DMLoad(dmnew, v);CHKERRQ(ierr);
69806c51deSksagiyam   }
705cca0b87SVaclav Hapla   ierr = PetscObjectSetName((PetscObject)dmnew,"Mesh_new");CHKERRQ(ierr);
715cca0b87SVaclav Hapla 
725cca0b87SVaclav Hapla   ierr = PetscViewerPopFormat(v);CHKERRQ(ierr);
735cca0b87SVaclav Hapla   ierr = PetscViewerDestroy(&v);CHKERRQ(ierr);
745cca0b87SVaclav Hapla   *dm_new = dmnew;
755cca0b87SVaclav Hapla   PetscFunctionReturn(0);
765cca0b87SVaclav Hapla }
775cca0b87SVaclav Hapla 
785cca0b87SVaclav Hapla int main(int argc, char **argv)
795cca0b87SVaclav Hapla {
805cca0b87SVaclav Hapla   DM             dm, dmnew;
815cca0b87SVaclav Hapla   PetscPartitioner part;
825cca0b87SVaclav Hapla   AppCtx         user;
835cca0b87SVaclav Hapla   PetscBool      flg;
845cca0b87SVaclav Hapla   PetscErrorCode ierr;
855cca0b87SVaclav Hapla 
865cca0b87SVaclav Hapla   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
875cca0b87SVaclav Hapla   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
885cca0b87SVaclav Hapla   ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.interpolate, &dm);CHKERRQ(ierr);
895cca0b87SVaclav Hapla   ierr = DMSetOptionsPrefix(dm,"orig_");CHKERRQ(ierr);
905cca0b87SVaclav Hapla   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
915cca0b87SVaclav Hapla 
925cca0b87SVaclav Hapla   if (user.distribute) {
935cca0b87SVaclav Hapla     DM dmdist;
945cca0b87SVaclav Hapla 
955cca0b87SVaclav Hapla     ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr);
965cca0b87SVaclav Hapla     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
975cca0b87SVaclav Hapla     ierr = DMPlexDistribute(dm, 0, NULL, &dmdist);CHKERRQ(ierr);
985cca0b87SVaclav Hapla     if (dmdist) {
995cca0b87SVaclav Hapla       ierr = DMDestroy(&dm);CHKERRQ(ierr);
1005cca0b87SVaclav Hapla       dm   = dmdist;
1015cca0b87SVaclav Hapla     }
1025cca0b87SVaclav Hapla   }
1035cca0b87SVaclav Hapla 
1045cca0b87SVaclav Hapla   ierr = DMSetOptionsPrefix(dm,NULL);CHKERRQ(ierr);
1055cca0b87SVaclav Hapla   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
1065cca0b87SVaclav Hapla   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
1075cca0b87SVaclav Hapla 
108806c51deSksagiyam   ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew);CHKERRQ(ierr);
1095cca0b87SVaclav Hapla 
1105cca0b87SVaclav Hapla   if (user.second_write_read) {
1115cca0b87SVaclav Hapla     ierr = DMDestroy(&dm);CHKERRQ(ierr);
1125cca0b87SVaclav Hapla     dm = dmnew;
113806c51deSksagiyam     ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew);CHKERRQ(ierr);
1145cca0b87SVaclav Hapla   }
1155cca0b87SVaclav Hapla 
1165cca0b87SVaclav Hapla   ierr = DMViewFromOptions(dmnew, NULL, "-dm_view");CHKERRQ(ierr);
1175cca0b87SVaclav Hapla   /* TODO: Is it still true? */
1185cca0b87SVaclav 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. */
1195cca0b87SVaclav Hapla 
1205cca0b87SVaclav Hapla   /* This currently makes sense only for sequential meshes. */
1215cca0b87SVaclav Hapla   if (user.compare) {
1225cca0b87SVaclav Hapla     ierr = DMPlexEqual(dmnew, dm, &flg);CHKERRQ(ierr);
1235cca0b87SVaclav Hapla     if (flg) {ierr = PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");CHKERRQ(ierr);}
1245cca0b87SVaclav Hapla     else     {ierr = PetscPrintf(PETSC_COMM_WORLD,"DMs are not equal\n");CHKERRQ(ierr);}
1255cca0b87SVaclav Hapla   }
1265cca0b87SVaclav Hapla 
1275cca0b87SVaclav Hapla   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1285cca0b87SVaclav Hapla   ierr = DMDestroy(&dmnew);CHKERRQ(ierr);
1295cca0b87SVaclav Hapla   ierr = PetscFinalize();
1305cca0b87SVaclav Hapla   return ierr;
1315cca0b87SVaclav Hapla }
1325cca0b87SVaclav Hapla 
1335cca0b87SVaclav Hapla /*TEST
1345cca0b87SVaclav Hapla   build:
1355cca0b87SVaclav Hapla     requires: hdf5
1365cca0b87SVaclav Hapla   # Idempotence of saving/loading
1375cca0b87SVaclav Hapla   #   Have to replace Exodus file, which is creating uninterpolated edges
1385cca0b87SVaclav Hapla   test:
1395cca0b87SVaclav Hapla     suffix: 0
1405cca0b87SVaclav Hapla     requires: exodusii broken
1415cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
1425cca0b87SVaclav Hapla     args: -format hdf5_petsc -compare
1435cca0b87SVaclav Hapla   test:
1445cca0b87SVaclav Hapla     suffix: 1
1455cca0b87SVaclav Hapla     requires: exodusii parmetis !define(PETSC_USE_64BIT_INDICES) broken
1465cca0b87SVaclav Hapla     nsize: 2
1475cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
1485cca0b87SVaclav Hapla     args: -petscpartitioner_type parmetis
1495cca0b87SVaclav Hapla     args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail
1505cca0b87SVaclav Hapla   testset:
1515cca0b87SVaclav Hapla     requires: exodusii
1525cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
1535cca0b87SVaclav Hapla     args: -petscpartitioner_type simple
1545cca0b87SVaclav Hapla     args: -dm_view ascii::ascii_info_detail
1555cca0b87SVaclav Hapla     args: -new_dm_view ascii::ascii_info_detail
1565cca0b87SVaclav Hapla     test:
1575cca0b87SVaclav Hapla       suffix: 2
1585cca0b87SVaclav Hapla       nsize: {{1 2 4 8}separate output}
1595cca0b87SVaclav Hapla       args: -format {{default hdf5_petsc}separate output}
1605cca0b87SVaclav Hapla       args: -interpolate {{0 1}separate output}
1615cca0b87SVaclav Hapla     test:
1625cca0b87SVaclav Hapla       suffix: 2a
1635cca0b87SVaclav Hapla       nsize: {{1 2 4 8}separate output}
1645cca0b87SVaclav Hapla       args: -format {{hdf5_xdmf hdf5_viz}separate output}
1655cca0b87SVaclav Hapla   test:
1665cca0b87SVaclav Hapla     suffix: 3
1675cca0b87SVaclav Hapla     requires: exodusii
1685cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare
1695cca0b87SVaclav Hapla 
1705cca0b87SVaclav Hapla   # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2
1715cca0b87SVaclav Hapla   test:
1725cca0b87SVaclav Hapla     suffix: 4
1735cca0b87SVaclav Hapla     requires: !complex
1745cca0b87SVaclav Hapla     nsize: {{1 2 3 4 8}}
1755cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5
1765cca0b87SVaclav Hapla     args: -dm_plex_create_from_hdf5_xdmf -distribute 0 -format hdf5_xdmf -second_write_read -compare
1775cca0b87SVaclav Hapla 
178806c51deSksagiyam   # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load()
179806c51deSksagiyam   # Output must be the same as ex55_2_nsize-2_format-hdf5_petsc_interpolate-0.out
180806c51deSksagiyam   test:
181806c51deSksagiyam     suffix: 5
182806c51deSksagiyam     requires: exodusii
183806c51deSksagiyam     nsize: 2
184806c51deSksagiyam     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
185806c51deSksagiyam     args: -petscpartitioner_type simple
186806c51deSksagiyam     args: -dm_view ascii::ascii_info_detail
187806c51deSksagiyam     args: -new_dm_view ascii::ascii_info_detail
188806c51deSksagiyam     args: -format hdf5_petsc -use_low_level_functions
189806c51deSksagiyam 
1905cca0b87SVaclav Hapla   testset:
1915cca0b87SVaclav Hapla     # the same data and settings as dm_impls_plex_tests-ex18_9%
1925cca0b87SVaclav Hapla     requires: hdf5 !complex datafilespath
1935cca0b87SVaclav Hapla     #TODO DMPlexCheckPointSF() fails for nsize 4
1945cca0b87SVaclav Hapla     nsize: {{1 2}}
1955cca0b87SVaclav Hapla     args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
1965cca0b87SVaclav 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
1975cca0b87SVaclav Hapla     args: -format hdf5_xdmf -second_write_read -compare
1985cca0b87SVaclav Hapla     test:
1995cca0b87SVaclav Hapla       suffix: 9_hdf5_seqload
2005cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type simple
2015cca0b87SVaclav Hapla       args: -interpolate {{0 1}}
2025cca0b87SVaclav Hapla       args: -dm_plex_hdf5_force_sequential
2035cca0b87SVaclav Hapla     test:
2045cca0b87SVaclav Hapla       suffix: 9_hdf5_seqload_metis
2055cca0b87SVaclav Hapla       requires: parmetis
2065cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type parmetis
2075cca0b87SVaclav Hapla       args: -interpolate 1
2085cca0b87SVaclav Hapla       args: -dm_plex_hdf5_force_sequential
2095cca0b87SVaclav Hapla     test:
2105cca0b87SVaclav Hapla       suffix: 9_hdf5
2115cca0b87SVaclav Hapla       args: -interpolate 1 -petscpartitioner_type simple
2125cca0b87SVaclav Hapla     test:
2135cca0b87SVaclav Hapla       suffix: 9_hdf5_repart
2145cca0b87SVaclav Hapla       requires: parmetis
2155cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type parmetis
2165cca0b87SVaclav Hapla       args: -interpolate 1
2175cca0b87SVaclav Hapla     test:
2185cca0b87SVaclav Hapla       TODO: Parallel partitioning of uninterpolated meshes not supported
2195cca0b87SVaclav Hapla       suffix: 9_hdf5_repart_ppu
2205cca0b87SVaclav Hapla       requires: parmetis
2215cca0b87SVaclav Hapla       args: -distribute -petscpartitioner_type parmetis
2225cca0b87SVaclav Hapla       args: -interpolate 0
2235cca0b87SVaclav Hapla 
2245cca0b87SVaclav Hapla   # reproduce PetscSFView() crash - fixed, left as regression test
2255cca0b87SVaclav Hapla   test:
2265cca0b87SVaclav Hapla     suffix: new_dm_view
2275cca0b87SVaclav Hapla     requires: exodusii
2285cca0b87SVaclav Hapla     nsize: 2
2295cca0b87SVaclav Hapla     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail
2305cca0b87SVaclav Hapla TEST*/
231