Actual source code: ex5.c
petsc-3.11.4 2019-09-28
1: static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n";
3: #include <petscdmplex.h>
4: #include <petscviewerhdf5.h>
5: #include <petscsf.h>
7: typedef struct {
8: PetscBool compare; /* Compare the meshes using DMPlexEqual() */
9: PetscBool distribute; /* Distribute the mesh */
10: PetscBool interpolate; /* Generate intermediate mesh elements */
11: char filename[PETSC_MAX_PATH_LEN]; /* Mesh filename */
12: PetscViewerFormat format; /* Format to write and read */
13: PetscBool second_write_read; /* Write and read for the 2nd time */
14: } AppCtx;
16: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
17: {
21: options->compare = PETSC_FALSE;
22: options->distribute = PETSC_TRUE;
23: options->interpolate = PETSC_FALSE;
24: options->filename[0] = '\0';
25: options->format = PETSC_VIEWER_DEFAULT;
26: options->second_write_read = PETSC_FALSE;
28: PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
29: PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex5.c", options->compare, &options->compare, NULL);
30: PetscOptionsBool("-distribute", "Distribute the mesh", "ex5.c", options->distribute, &options->distribute, NULL);
31: PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex5.c", options->interpolate, &options->interpolate, NULL);
32: PetscOptionsString("-filename", "The mesh file", "ex5.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);
33: PetscOptionsEnum("-format", "Format to write and read", "ex5.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL);
34: PetscOptionsBool("-second_write_read", "Write and read for the 2nd time", "ex5.c", options->second_write_read, &options->second_write_read, NULL);
35: PetscOptionsEnd();
36: return(0);
37: };
39: static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], PetscViewerFormat format, const char prefix[], DM *dm_new)
40: {
41: DM dmnew;
42: PetscViewer v;
46: PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v);
47: PetscViewerPushFormat(v, format);
48: DMView(dm, v);
50: PetscViewerFileSetMode(v, FILE_MODE_READ);
51: DMCreate(PETSC_COMM_WORLD, &dmnew);
52: DMSetType(dmnew, DMPLEX);
53: DMSetOptionsPrefix(dmnew, prefix);
54: DMLoad(dmnew, v);
56: PetscViewerPopFormat(v);
57: PetscViewerDestroy(&v);
58: *dm_new = dmnew;
59: return(0);
60: }
62: int main(int argc, char **argv)
63: {
64: DM dm, dmnew;
65: PetscPartitioner part;
66: AppCtx user;
67: PetscBool flg;
70: PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
71: ProcessOptions(PETSC_COMM_WORLD, &user);
72: DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.interpolate, &dm);
73: DMSetOptionsPrefix(dm,"orig_");
74: DMViewFromOptions(dm, NULL, "-dm_view");
76: if (user.distribute) {
77: DM dmdist;
79: DMPlexGetPartitioner(dm, &part);
80: PetscPartitionerSetFromOptions(part);
81: DMPlexDistribute(dm, 0, NULL, &dmdist);
82: if (dmdist) {
83: DMDestroy(&dm);
84: dm = dmdist;
85: }
86: }
88: DMSetOptionsPrefix(dm,NULL);
89: DMSetFromOptions(dm);
90: DMViewFromOptions(dm, NULL, "-dm_view");
92: DMPlexWriteAndReadHDF5(dm, "dmdist.h5", user.format, "new_", &dmnew);
94: if (user.second_write_read) {
95: DMDestroy(&dm);
96: dm = dmnew;
97: DMPlexWriteAndReadHDF5(dm, "dmdist.h5", user.format, "new_", &dmnew);
98: }
100: DMViewFromOptions(dmnew, NULL, "-dm_view");
101: /* TODO: Is it still true? */
102: /* The NATIVE format for coordiante viewing is killing parallel output, since we have a local vector. Map it to global, and it will work. */
104: /* This currently makes sense only for sequential meshes. */
105: if (user.compare) {
106: DMPlexEqual(dmnew, dm, &flg);
107: if (flg) {PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");}
108: else {PetscPrintf(PETSC_COMM_WORLD,"DMs are not equal\n");}
109: }
111: DMDestroy(&dm);
112: DMDestroy(&dmnew);
113: PetscFinalize();
114: return ierr;
115: }
117: /*TEST
118: build:
119: requires: hdf5
120: # Idempotence of saving/loading
121: # Have to replace Exodus file, which is creating uninterpolated edges
122: test:
123: suffix: 0
124: requires: exodusii broken
125: args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
126: args: -format hdf5_petsc -compare
127: test:
128: suffix: 1
129: requires: exodusii parmetis !define(PETSC_USE_64BIT_INDICES) broken
130: nsize: 2
131: args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
132: args: -petscpartitioner_type parmetis
133: args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail
134: testset:
135: requires: exodusii
136: args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
137: args: -petscpartitioner_type simple
138: args: -dm_view ascii::ascii_info_detail
139: args: -new_dm_view ascii::ascii_info_detail
140: test:
141: suffix: 2
142: nsize: {{1 2 4 8}separate output}
143: args: -format {{default hdf5_petsc}separate output}
144: args: -interpolate {{0 1}separate output}
145: test:
146: suffix: 2a
147: nsize: {{1 2 4 8}separate output}
148: args: -format {{hdf5_xdmf hdf5_viz}separate output}
149: test:
150: suffix: 3
151: requires: exodusii
152: args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare
154: # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2
155: test:
156: suffix: 4
157: requires: !complex
158: nsize: {{1 2 3 4 8}}
159: args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5
160: args: -dm_plex_create_from_hdf5_xdmf -distribute 0 -format hdf5_xdmf -second_write_read -compare
162: # reproduce PetscSFView() crash - fixed, left as regression test
163: test:
164: suffix: new_dm_view
165: requires: exodusii
166: nsize: 2
167: args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail
168: TEST*/