xref: /petsc/src/dm/impls/plex/tests/ex51.c (revision aee58fc643cde3a74e29e886bde44205d937e3a7)
1d2522c19Sksagiyam static char help[] = "Tests save/load plex with distribution in HDF5.\n\n";
2d2522c19Sksagiyam 
3d2522c19Sksagiyam #include <petscdmplex.h>
4d2522c19Sksagiyam #include <petscsf.h>
5d2522c19Sksagiyam #include <petsclayouthdf5.h>
6d2522c19Sksagiyam 
7d2522c19Sksagiyam typedef struct {
8d2522c19Sksagiyam   char fname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */
9d2522c19Sksagiyam } AppCtx;
10d2522c19Sksagiyam 
ProcessOptions(MPI_Comm comm,AppCtx * options)11d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12d71ae5a4SJacob Faibussowitsch {
13d2522c19Sksagiyam   PetscBool flg;
14d2522c19Sksagiyam 
15d2522c19Sksagiyam   PetscFunctionBegin;
16d2522c19Sksagiyam   options->fname[0] = '\0';
17d2522c19Sksagiyam   PetscOptionsBegin(comm, "", "DMPlex View/Load Test Options", "DMPLEX");
18d2522c19Sksagiyam   PetscCall(PetscOptionsString("-fname", "The output mesh file", "ex51.c", options->fname, options->fname, sizeof(options->fname), &flg));
19d2522c19Sksagiyam   PetscOptionsEnd();
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21d2522c19Sksagiyam }
22d2522c19Sksagiyam 
main(int argc,char ** argv)23d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
24d71ae5a4SJacob Faibussowitsch {
25d2522c19Sksagiyam   const char        exampleDMPlexName[]       = "exampleDMPlex";
26d2522c19Sksagiyam   const char        exampleDistributionName[] = "exampleDistribution";
27d2522c19Sksagiyam   PetscViewerFormat format                    = PETSC_VIEWER_HDF5_PETSC;
28d2522c19Sksagiyam   AppCtx            user;
29d2522c19Sksagiyam 
30d2522c19Sksagiyam   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
31d2522c19Sksagiyam   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
32d2522c19Sksagiyam   /* Save */
33d2522c19Sksagiyam   {
34d2522c19Sksagiyam     DM          dm;
35d2522c19Sksagiyam     PetscViewer viewer;
36d2522c19Sksagiyam 
37d2522c19Sksagiyam     PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, user.fname, FILE_MODE_WRITE, &viewer));
38d2522c19Sksagiyam     /* Save exampleDMPlex */
39d2522c19Sksagiyam     {
40d2522c19Sksagiyam       DM               pdm;
41d2522c19Sksagiyam       PetscPartitioner part;
42d2522c19Sksagiyam       const PetscInt   faces[2] = {6, 1};
43d2522c19Sksagiyam       PetscSF          sf;
44d2522c19Sksagiyam       PetscInt         overlap = 1;
45d2522c19Sksagiyam 
46*42108689Sksagiyam       PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, 2, PETSC_FALSE, faces, NULL, NULL, NULL, PETSC_TRUE, 0, PETSC_TRUE, &dm));
47d2522c19Sksagiyam       PetscCall(DMPlexGetPartitioner(dm, &part));
48d2522c19Sksagiyam       PetscCall(PetscPartitionerSetFromOptions(part));
49d2522c19Sksagiyam       PetscCall(DMPlexDistribute(dm, overlap, &sf, &pdm));
50d2522c19Sksagiyam       if (pdm) {
51d2522c19Sksagiyam         PetscCall(DMDestroy(&dm));
52d2522c19Sksagiyam         dm = pdm;
53d2522c19Sksagiyam       }
54d2522c19Sksagiyam       PetscCall(PetscSFDestroy(&sf));
55d2522c19Sksagiyam       PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName));
56d2522c19Sksagiyam       PetscCall(DMPlexDistributionSetName(dm, exampleDistributionName));
57d2522c19Sksagiyam       PetscCall(PetscViewerPushFormat(viewer, format));
58d2522c19Sksagiyam       PetscCall(DMPlexTopologyView(dm, viewer));
59d2522c19Sksagiyam       PetscCall(DMPlexLabelsView(dm, viewer));
60d2522c19Sksagiyam       PetscCall(PetscViewerPopFormat(viewer));
61d2522c19Sksagiyam     }
62d2522c19Sksagiyam     /* Save coordinates */
63d2522c19Sksagiyam     PetscCall(PetscViewerPushFormat(viewer, format));
64d2522c19Sksagiyam     PetscCall(DMPlexCoordinatesView(dm, viewer));
65d2522c19Sksagiyam     PetscCall(PetscViewerPopFormat(viewer));
66d2522c19Sksagiyam     PetscCall(PetscObjectSetName((PetscObject)dm, "Save: DM"));
67d2522c19Sksagiyam     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
68d2522c19Sksagiyam     PetscCall(DMDestroy(&dm));
69d2522c19Sksagiyam     PetscCall(PetscViewerDestroy(&viewer));
70d2522c19Sksagiyam   }
71d2522c19Sksagiyam   /* Load */
72d2522c19Sksagiyam   {
73d2522c19Sksagiyam     DM          dm;
74d2522c19Sksagiyam     PetscSF     sfXC;
75d2522c19Sksagiyam     PetscViewer viewer;
76d2522c19Sksagiyam 
77d2522c19Sksagiyam     PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, user.fname, FILE_MODE_READ, &viewer));
78d2522c19Sksagiyam     /* Load exampleDMPlex */
79d2522c19Sksagiyam     PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
80d2522c19Sksagiyam     PetscCall(DMSetType(dm, DMPLEX));
81d2522c19Sksagiyam     PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName));
82d2522c19Sksagiyam     PetscCall(DMPlexDistributionSetName(dm, exampleDistributionName));
83d2522c19Sksagiyam     /* sfXC: X -> C                         */
84d2522c19Sksagiyam     /* X: set of globalPointNumbers, [0, N) */
85d2522c19Sksagiyam     /* C: loaded in-memory plex             */
86d2522c19Sksagiyam     PetscCall(PetscViewerPushFormat(viewer, format));
87d2522c19Sksagiyam     PetscCall(DMPlexTopologyLoad(dm, viewer, &sfXC));
88d2522c19Sksagiyam     PetscCall(PetscViewerPopFormat(viewer));
89d2522c19Sksagiyam     /* Do not distribute (Already distributed just like the saved plex) */
90d2522c19Sksagiyam     /* Load labels */
91d2522c19Sksagiyam     PetscCall(DMPlexLabelsLoad(dm, viewer, sfXC));
92d2522c19Sksagiyam     /* Load coordinates */
93d2522c19Sksagiyam     PetscCall(PetscViewerPushFormat(viewer, format));
94d2522c19Sksagiyam     PetscCall(DMPlexCoordinatesLoad(dm, viewer, sfXC));
95d2522c19Sksagiyam     PetscCall(PetscViewerPopFormat(viewer));
96d2522c19Sksagiyam     PetscCall(DMSetFromOptions(dm));
97d2522c19Sksagiyam     /* Print the exact same plex as the saved one */
98d2522c19Sksagiyam     PetscCall(PetscObjectSetName((PetscObject)dm, "Load: DM"));
99d2522c19Sksagiyam     PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
100d2522c19Sksagiyam     PetscCall(PetscSFDestroy(&sfXC));
101d2522c19Sksagiyam     PetscCall(DMDestroy(&dm));
102d2522c19Sksagiyam     PetscCall(PetscViewerDestroy(&viewer));
103d2522c19Sksagiyam   }
104d2522c19Sksagiyam   /* Finalize */
105d2522c19Sksagiyam   PetscCall(PetscFinalize());
106d2522c19Sksagiyam   return 0;
107d2522c19Sksagiyam }
108d2522c19Sksagiyam 
109d2522c19Sksagiyam /*TEST
110d2522c19Sksagiyam 
111d2522c19Sksagiyam   build:
112d2522c19Sksagiyam     requires: hdf5
113d2522c19Sksagiyam     requires: !complex
114d2522c19Sksagiyam   test:
115d2522c19Sksagiyam     suffix: 0
116d2522c19Sksagiyam     requires: parmetis
117d2522c19Sksagiyam     nsize: {{1 2 4}separate output}
118d2522c19Sksagiyam     args: -fname ex51_dump.h5 -dm_view ascii::ascii_info_detail
119d2522c19Sksagiyam     args: -petscpartitioner_type parmetis
120f1dec35aSVaclav Hapla     args: -dm_plex_view_hdf5_storage_version 2.1.0
121d2522c19Sksagiyam 
122d2522c19Sksagiyam TEST*/
123