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