xref: /petsc/src/dm/impls/forest/tests/ex3.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "Tests adaptive refinement using DMForest, and uses HDF5.\n\n";
2 
3 #include <petscdmforest.h>
4 #include <petscdmplex.h>
5 #include <petscviewerhdf5.h>
6 
7 int main (int argc, char **argv)
8 {
9   DM             base, forest, plex;
10   PetscSection   s;
11   PetscViewer    viewer;
12   Vec            g = NULL, g2 = NULL;
13   PetscReal      nrm;
14   PetscBool      adapt = PETSC_FALSE, userSection = PETSC_FALSE;
15   PetscInt       vStart, vEnd, v, i;
16 
17   PetscFunctionBeginUser;
18   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
19   PetscCall(PetscOptionsGetBool(NULL, NULL, "-adapt", &adapt, NULL));
20   PetscCall(PetscOptionsGetBool(NULL, NULL, "-user_section", &userSection, NULL));
21 
22   /* Create a base DMPlex mesh */
23   PetscCall(DMCreate(PETSC_COMM_WORLD, &base));
24   PetscCall(DMSetType(base, DMPLEX));
25   PetscCall(DMSetFromOptions(base));
26   PetscCall(DMViewFromOptions(base, NULL, "-dm_view"));
27 
28   /* Covert Plex mesh to Forest and destroy base */
29   PetscCall(DMCreate(PETSC_COMM_WORLD, &forest));
30   PetscCall(DMSetType(forest, DMP4EST));
31   PetscCall(DMForestSetBaseDM(forest, base));
32   PetscCall(DMSetUp(forest));
33   PetscCall(DMDestroy(&base));
34   PetscCall(DMViewFromOptions(forest, NULL, "-dm_view"));
35 
36   if (adapt) {
37     /* Adaptively refine the cell 0 of the mesh */
38     for (i = 0; i < 3; ++i) {
39       DM      postforest;
40       DMLabel adaptLabel = NULL;
41 
42       PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
43       PetscCall(DMLabelSetValue(adaptLabel, 0, DM_ADAPT_REFINE));
44       PetscCall(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest));
45       PetscCall(DMForestSetAdaptivityLabel(postforest, adaptLabel));
46       PetscCall(DMLabelDestroy(&adaptLabel));
47       PetscCall(DMSetUp(postforest));
48       PetscCall(DMDestroy(&forest));
49       forest = postforest;
50     }
51   } else {
52     /* Adaptively refine all cells of the mesh */
53     PetscInt cStart, cEnd, c;
54 
55     for (i = 0; i < 3; ++i) {
56       DM      postforest;
57       DMLabel adaptLabel = NULL;
58 
59       PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
60 
61       PetscCall(DMForestGetCellChart(forest, &cStart, &cEnd));
62       for (c = cStart; c < cEnd; ++c) {
63         PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));
64       }
65 
66       PetscCall(DMForestTemplate(forest, PETSC_COMM_WORLD, &postforest));
67       PetscCall(DMForestSetAdaptivityLabel(postforest, adaptLabel));
68       PetscCall(DMLabelDestroy(&adaptLabel));
69       PetscCall(DMSetUp(postforest));
70       PetscCall(DMDestroy(&forest));
71       forest = postforest;
72     }
73   }
74   PetscCall(DMViewFromOptions(forest, NULL, "-dm_view"));
75 
76   /*  Setup the section*/
77   if (userSection) {
78     PetscCall(DMConvert(forest, DMPLEX, &plex));
79     PetscCall(DMViewFromOptions(plex, NULL, "-dm_view"));
80     PetscCall(DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd));
81     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "Vertices [%" PetscInt_FMT ", %" PetscInt_FMT ")\n", vStart, vEnd));
82     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
83     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) forest), &s));
84     PetscCall(PetscSectionSetNumFields(s, 1));
85     PetscCall(PetscSectionSetChart(s, vStart, vEnd));
86     for (v = vStart; v < vEnd; ++v) {
87       PetscCall(PetscSectionSetDof(s, v, 1));
88       PetscCall(PetscSectionSetFieldDof(s, v, 0, 1));
89     }
90     PetscCall(PetscSectionSetUp(s));
91     PetscCall(DMSetLocalSection(forest, s));
92     PetscCall(PetscObjectViewFromOptions((PetscObject) s, NULL, "-my_section_view"));
93     PetscCall(PetscSectionDestroy(&s));
94     PetscCall(DMDestroy(&plex));
95   } else {
96     PetscFE  fe;
97     PetscInt dim;
98 
99     PetscCall(DMGetDimension(forest, &dim));
100     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe));
101     PetscCall(DMAddField(forest, NULL, (PetscObject) fe));
102     PetscCall(PetscFEDestroy(&fe));
103     PetscCall(DMCreateDS(forest));
104   }
105 
106   /* Create the global vector*/
107   PetscCall(DMCreateGlobalVector(forest, &g));
108   PetscCall(PetscObjectSetName((PetscObject) g, "g"));
109   PetscCall(VecSet(g, 1.0));
110 
111   /* Test global to local*/
112   Vec l;
113   PetscCall(DMCreateLocalVector(forest, &l));
114   PetscCall(VecZeroEntries(l));
115   PetscCall(DMGlobalToLocal(forest, g, INSERT_VALUES, l));
116   PetscCall(VecZeroEntries(g));
117   PetscCall(DMLocalToGlobal(forest, l, INSERT_VALUES, g));
118   PetscCall(VecDestroy(&l));
119 
120   /*  Save a vector*/
121   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_WRITE, &viewer));
122   PetscCall(VecView(g, viewer));
123   PetscCall(PetscViewerDestroy(&viewer));
124 
125   /* Load another vector to load into*/
126   PetscCall(DMCreateGlobalVector(forest, &g2));
127   PetscCall(PetscObjectSetName((PetscObject) g2, "g"));
128   PetscCall(VecZeroEntries(g2));
129 
130   /*  Load a vector*/
131   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "forestHDF.h5", FILE_MODE_READ, &viewer));
132   PetscCall(VecLoad(g2, viewer));
133   PetscCall(PetscViewerDestroy(&viewer));
134 
135   /*  Check if the data is the same*/
136   PetscCall(VecAXPY(g2, -1.0, g));
137   PetscCall(VecNorm(g2, NORM_INFINITY, &nrm));
138   PetscCheck(PetscAbsReal(nrm) <= PETSC_SMALL,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid difference norm %g", (double) nrm);
139 
140   PetscCall(VecDestroy(&g));
141   PetscCall(VecDestroy(&g2));
142   PetscCall(DMDestroy(&forest));
143   PetscCall(PetscFinalize());
144   return 0;
145 }
146 
147 /*TEST
148 
149   build:
150     requires: hdf5 p4est
151 
152   test:
153     suffix: 0
154     nsize: {{1 2 5}}
155     args: -adapt -dm_plex_simplex 0 -dm_plex_box_faces 2,2
156 
157 TEST*/
158