xref: /petsc/src/dm/impls/plex/tests/ex55.c (revision d70f29a362ae60d541be4e9a72e9494be00f9e3d)
1 static char help[] = "Load and save the mesh and fields to HDF5 and ExodusII\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscviewerhdf5.h>
5 #include <petscsf.h>
6 
7 typedef struct {
8   PetscBool compare;                      /* Compare the meshes using DMPlexEqual() */
9   PetscBool compare_labels;               /* Compare labels in the meshes using DMCompareLabels() */
10   PetscBool distribute;                   /* Distribute the mesh */
11   PetscBool interpolate;                  /* Generate intermediate mesh elements */
12   char      filename[PETSC_MAX_PATH_LEN]; /* Mesh filename */
13   char      meshname[PETSC_MAX_PATH_LEN]; /* Mesh name */
14   PetscViewerFormat format;               /* Format to write and read */
15   PetscBool second_write_read;            /* Write and read for the 2nd time */
16   PetscBool use_low_level_functions;      /* Use low level functions for viewing and loading */
17 } AppCtx;
18 
19 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
20 {
21   PetscErrorCode ierr;
22 
23   PetscFunctionBeginUser;
24   options->compare = PETSC_FALSE;
25   options->compare_labels = PETSC_FALSE;
26   options->distribute = PETSC_TRUE;
27   options->interpolate = PETSC_FALSE;
28   options->filename[0] = '\0';
29   options->meshname[0] = '\0';
30   options->format = PETSC_VIEWER_DEFAULT;
31   options->second_write_read = PETSC_FALSE;
32   options->use_low_level_functions = PETSC_FALSE;
33 
34   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
35   ierr = PetscOptionsBool("-compare", "Compare the meshes using DMPlexEqual()", "ex55.c", options->compare, &options->compare, NULL);CHKERRQ(ierr);
36   ierr = PetscOptionsBool("-compare_labels", "Compare labels in the meshes using DMCompareLabels()", "ex55.c", options->compare_labels, &options->compare_labels, NULL);CHKERRQ(ierr);
37   ierr = PetscOptionsBool("-distribute", "Distribute the mesh", "ex55.c", options->distribute, &options->distribute, NULL);CHKERRQ(ierr);
38   ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex55.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
39   ierr = PetscOptionsString("-filename", "The mesh file", "ex55.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
40   ierr = PetscOptionsString("-meshname", "The mesh file", "ex55.c", options->meshname, options->meshname, sizeof(options->meshname), NULL);CHKERRQ(ierr);
41   ierr = PetscOptionsEnum("-format", "Format to write and read", "ex55.c", PetscViewerFormats, (PetscEnum)options->format, (PetscEnum*)&options->format, NULL);CHKERRQ(ierr);
42   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);
43   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);
44   ierr = PetscOptionsEnd();CHKERRQ(ierr);
45   PetscFunctionReturn(0);
46 };
47 
48 static PetscErrorCode DMPlexWriteAndReadHDF5(DM dm, const char filename[], const char prefix[], AppCtx user, DM *dm_new)
49 {
50   DM             dmnew;
51   const char     savedName[]  = "Mesh";
52   const char     loadedName[] = "Mesh_new";
53   PetscViewer    v;
54   PetscErrorCode ierr;
55 
56   PetscFunctionBeginUser;
57   ierr = PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), filename, FILE_MODE_WRITE, &v);CHKERRQ(ierr);
58   ierr = PetscViewerPushFormat(v, user.format);CHKERRQ(ierr);
59   ierr = PetscObjectSetName((PetscObject) dm, savedName);CHKERRQ(ierr);
60   if (user.use_low_level_functions) {
61     ierr = DMPlexTopologyView(dm, v);CHKERRQ(ierr);
62     ierr = DMPlexCoordinatesView(dm, v);CHKERRQ(ierr);
63     ierr = DMPlexLabelsView(dm, v);CHKERRQ(ierr);
64   } else {
65     ierr = DMView(dm, v);CHKERRQ(ierr);
66   }
67 
68   ierr = PetscViewerFileSetMode(v, FILE_MODE_READ);CHKERRQ(ierr);
69   ierr = DMCreate(PETSC_COMM_WORLD, &dmnew);CHKERRQ(ierr);
70   ierr = DMSetType(dmnew, DMPLEX);CHKERRQ(ierr);
71   ierr = PetscObjectSetName((PetscObject) dmnew, savedName);CHKERRQ(ierr);
72   ierr = DMSetOptionsPrefix(dmnew, prefix);CHKERRQ(ierr);
73   if (user.use_low_level_functions) {
74     PetscSF  sfXC;
75 
76     ierr = DMPlexTopologyLoad(dmnew, v, &sfXC);CHKERRQ(ierr);
77     ierr = DMPlexCoordinatesLoad(dmnew, v, sfXC);CHKERRQ(ierr);
78     ierr = DMPlexLabelsLoad(dmnew, v, sfXC);CHKERRQ(ierr);
79     ierr = PetscSFDestroy(&sfXC);CHKERRQ(ierr);
80   } else {
81     ierr = DMLoad(dmnew, v);CHKERRQ(ierr);
82   }
83   ierr = PetscObjectSetName((PetscObject)dmnew,loadedName);CHKERRQ(ierr);
84 
85   ierr = PetscViewerPopFormat(v);CHKERRQ(ierr);
86   ierr = PetscViewerDestroy(&v);CHKERRQ(ierr);
87   *dm_new = dmnew;
88   PetscFunctionReturn(0);
89 }
90 
91 int main(int argc, char **argv)
92 {
93   DM             dm, dmnew;
94   PetscPartitioner part;
95   AppCtx         user;
96   PetscBool      flg;
97   PetscErrorCode ierr;
98 
99   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
100   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
101   ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, user.filename, user.meshname, user.interpolate, &dm);CHKERRQ(ierr);
102   ierr = DMSetOptionsPrefix(dm,"orig_");CHKERRQ(ierr);
103   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
104 
105   if (user.distribute) {
106     DM dmdist;
107 
108     ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr);
109     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
110     ierr = DMPlexDistribute(dm, 0, NULL, &dmdist);CHKERRQ(ierr);
111     if (dmdist) {
112       ierr = DMDestroy(&dm);CHKERRQ(ierr);
113       dm   = dmdist;
114     }
115   }
116 
117   ierr = DMSetOptionsPrefix(dm,NULL);CHKERRQ(ierr);
118   ierr = DMPlexDistributeSetDefault(dm, PETSC_FALSE);CHKERRQ(ierr);
119   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
120   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
121 
122   ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew);CHKERRQ(ierr);
123 
124   if (user.second_write_read) {
125     ierr = DMDestroy(&dm);CHKERRQ(ierr);
126     dm = dmnew;
127     ierr = DMPlexWriteAndReadHDF5(dm, "dmdist.h5", "new_", user, &dmnew);CHKERRQ(ierr);
128   }
129 
130   ierr = DMViewFromOptions(dmnew, NULL, "-dm_view");CHKERRQ(ierr);
131   /* TODO: Is it still true? */
132   /* The NATIVE format for coordiante viewing is killing parallel output, since we have a local vector. Map it to global, and it will work. */
133 
134   /* This currently makes sense only for sequential meshes. */
135   if (user.compare) {
136     ierr = DMPlexEqual(dmnew, dm, &flg);CHKERRQ(ierr);
137     PetscCheckFalse(!flg,PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DMs are not equal");
138     ierr = PetscPrintf(PETSC_COMM_WORLD,"DMs equal\n");CHKERRQ(ierr);
139   }
140   if (user.compare_labels) {
141     ierr = DMCompareLabels(dmnew, dm, NULL, NULL);CHKERRQ(ierr);
142     ierr = PetscPrintf(PETSC_COMM_WORLD,"DMLabels equal\n");CHKERRQ(ierr);
143   }
144 
145   ierr = DMDestroy(&dm);CHKERRQ(ierr);
146   ierr = DMDestroy(&dmnew);CHKERRQ(ierr);
147   ierr = PetscFinalize();
148   return ierr;
149 }
150 
151 /*TEST
152   build:
153     requires: hdf5
154   # Idempotence of saving/loading
155   #   Have to replace Exodus file, which is creating uninterpolated edges
156   test:
157     suffix: 0
158     TODO: broken
159     requires: exodusii
160     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
161     args: -format hdf5_petsc -compare
162   test:
163     suffix: 1
164     TODO: broken
165     requires: exodusii parmetis !defined(PETSC_USE_64BIT_INDICES)
166     nsize: 2
167     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/Rect-tri3.exo -dm_view ascii::ascii_info_detail
168     args: -petscpartitioner_type parmetis
169     args: -format hdf5_petsc -new_dm_view ascii::ascii_info_detail
170   testset:
171     requires: exodusii
172     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
173     args: -petscpartitioner_type simple
174     args: -dm_view ascii::ascii_info_detail
175     args: -new_dm_view ascii::ascii_info_detail
176     test:
177       suffix: 2
178       nsize: {{1 2 4 8}separate output}
179       args: -format {{default hdf5_petsc}separate output}
180       args: -interpolate {{0 1}separate output}
181     test:
182       suffix: 2a
183       nsize: {{1 2 4 8}separate output}
184       args: -format {{hdf5_xdmf hdf5_viz}separate output}
185   test:
186     suffix: 3
187     requires: exodusii
188     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -compare -compare_labels
189 
190   # Load HDF5 file in XDMF format in parallel, write, read dm1, write, read dm2, and compare dm1 and dm2
191   testset:
192     suffix: 4
193     requires: !complex
194     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
195     args: -distribute 0 -second_write_read -compare
196     test:
197       suffix: hdf5_petsc
198       nsize: {{1 2}}
199       args: -format hdf5_petsc -compare_labels
200     test:
201       suffix: hdf5_xdmf
202       nsize: {{1 3 8}}
203       args: -format hdf5_xdmf
204 
205   # Use low level functions, DMPlexTopologyView()/Load(), DMPlexCoordinatesView()/Load(), and DMPlexLabelsView()/Load()
206   # Output must be the same as ex55_2_nsize-2_format-hdf5_petsc_interpolate-0.out
207   test:
208     suffix: 5
209     requires: exodusii
210     nsize: 2
211     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
212     args: -petscpartitioner_type simple
213     args: -dm_view ascii::ascii_info_detail
214     args: -new_dm_view ascii::ascii_info_detail
215     args: -format hdf5_petsc -use_low_level_functions
216 
217   testset:
218     suffix: 6
219     requires: hdf5 !complex datafilespath
220     nsize: {{1 3}}
221     args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
222     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
223     args: -format hdf5_petsc -second_write_read -compare -compare_labels
224     args: -interpolate {{0 1}} -distribute {{0 1}} -petscpartitioner_type simple
225 
226   testset:
227     # the same data and settings as dm_impls_plex_tests-ex18_9%
228     suffix: 9
229     requires: hdf5 !complex datafilespath
230     nsize: {{1 2 4}}
231     args: -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_geometry
232     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
233     args: -format hdf5_xdmf -second_write_read -compare
234     test:
235       suffix: hdf5_seqload
236       args: -distribute -petscpartitioner_type simple
237       args: -interpolate {{0 1}}
238       args: -dm_plex_hdf5_force_sequential
239     test:
240       suffix: hdf5_seqload_metis
241       requires: parmetis
242       args: -distribute -petscpartitioner_type parmetis
243       args: -interpolate 1
244       args: -dm_plex_hdf5_force_sequential
245     test:
246       suffix: hdf5
247       args: -interpolate 1 -petscpartitioner_type simple
248     test:
249       suffix: hdf5_repart
250       requires: parmetis
251       args: -distribute -petscpartitioner_type parmetis
252       args: -interpolate 1
253     test:
254       TODO: Parallel partitioning of uninterpolated meshes not supported
255       suffix: hdf5_repart_ppu
256       requires: parmetis
257       args: -distribute -petscpartitioner_type parmetis
258       args: -interpolate 0
259 
260   # reproduce PetscSFView() crash - fixed, left as regression test
261   test:
262     suffix: new_dm_view
263     requires: exodusii
264     nsize: 2
265     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo -new_dm_view ascii:ex5_new_dm_view.log:ascii_info_detail
266 
267   # test backward compatibility of petsc_hdf5 format
268   testset:
269     suffix: 10-v3.16.0-v1.0.0
270     requires: hdf5 !complex datafilespath
271     args: -dm_plex_check_all -compare -compare_labels
272     args: -dm_plex_view_hdf5_storage_version {{1.0.0 2.0.0}} -use_low_level_functions {{0 1}}
273     test:
274       suffix: a
275       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/annulus-20.h5
276     test:
277       suffix: b
278       TODO: broken
279       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/barycentricallyrefinedcube.h5
280     test:
281       suffix: c
282       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/blockcylinder-50.h5
283     test:
284       suffix: d
285       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/cube-hexahedra-refined.h5
286     test:
287       suffix: e
288       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/hybrid_hexwedge.h5
289     test:
290       suffix: f
291       args: -filename ${DATAFILESPATH}/meshes/hdf5-petsc/petsc-v3.16.0/v1.0.0/square.h5
292 
293 TEST*/
294