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