ExaDG
Loading...
Searching...
No Matches
hypercube_resolution_parameters.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_OPERATORS_HYPERCUBE_RESOLUTION_PARAMETERS_H_
23#define INCLUDE_EXADG_OPERATORS_HYPERCUBE_RESOLUTION_PARAMETERS_H_
24
25// deal.II
26#include <deal.II/base/parameter_handler.h>
27
28// ExaDG
29#include <exadg/grid/grid_data.h>
30#include <exadg/utilities/enum_patterns.h>
31
32namespace ExaDG
33{
34/*
35 * study throughput as a function of polynomial degree or problem size
36 */
37enum class RunType
38{
39 RefineHAndP, // run simulation for a specified range of mesh refinements and polynomial degrees
40 FixedProblemSize, // increase polynomial degree and keep problem size approximately constant
41 IncreasingProblemSize // run at fixed polynomial degree
42};
43
44/*
45 * Determines mesh refinement level l and number of subdivisions in 1d of hyper_cube mesh for a
46 * given number of minimal and maximal number of unknowns and a given RunType.
47 */
48inline void
49fill_resolutions_vector(
50 std::vector<
51 std::tuple<unsigned int /*k*/, unsigned int /*l*/, unsigned int /*subdivisions hypercube*/>> &
52 resolutions,
53 unsigned int const dim,
54 unsigned int const degree,
55 double const dofs_per_element,
56 dealii::types::global_dof_index const n_dofs_min,
57 dealii::types::global_dof_index const n_dofs_max,
58 RunType const & run_type,
59 ElementType const & element_type)
60{
61 unsigned int const resolutions_initial_size = resolutions.size();
62
63 double dofs_per_cube = dofs_per_element;
64
65 // in case of simplex elements, we create a mesh of hypercube elements that gets later subdivided
66 // into simplex elements using dealii::GridTools::subdivided_hyper_cube_with_simplices()
67 if(element_type == ElementType::Simplex)
68 {
69 unsigned int n_cells_per_cube = 1;
70
71 // Determine the number of simplex cells per hypercube cell
72 if(dim == 2)
73 n_cells_per_cube = 2;
74 else if(dim == 3)
75 n_cells_per_cube = 5;
76 else
77 AssertThrow(false, dealii::ExcMessage("not implemented."));
78
79 // From now on, we think of a hyerpcube mesh that gets later subdivided into simplex elements
80 // with n_cells_per_cube simplices per hypercube cell
81 dofs_per_cube *= n_cells_per_cube;
82 }
83
84 dealii::types::global_dof_index const n_cells_min =
85 (n_dofs_min + dofs_per_cube - 1) / dofs_per_cube;
86 dealii::types::global_dof_index const n_cells_max = n_dofs_max / dofs_per_cube;
87
88 // From the maximum number of cells, we derive a maximum refinement level for a uniformly refined
89 // mesh with one coarse-grid cells
90 unsigned int const refine_level_max =
91 int(std::log((double)n_cells_max) / std::log((double)dealii::Utilities::pow(2ULL, dim))) + 1;
92
93 // we start with the coarsest possible mesh and then increase the refine level by 1
94 // until we hit the maximum refinement level
95 unsigned int refine_level = 0;
96
97 // This loop is for a uniformly refined hypercube mesh. Inside the
98 // loop, we test whether additional combinations of coarse-grid cells and mesh refinement levels
99 // are suitable in terms of n_cells_min/n_cells_max. The reason behind is that we want to have
100 // more data points than just data points differing by a factor of 2^dim in the number of cells.
101 while(refine_level <= refine_level_max)
102 {
103 // To obtain additional data points, we test coarse grids with {3,4}^dim elements in 2D and
104 // {3,4,5}^dim elements in 3D. The cases with 3 and 5 subdivisions per coordinate direction are
105 // only relevant for refine_level >= 2.
106
107 // coarse grid with 3^dim cells, and refine_level-2 uniform refinements
108 if(refine_level >= 2)
109 {
110 unsigned int const n_subdivisions_1d = 3;
111 unsigned int const l = refine_level - 2;
112
113 unsigned int const n_cells =
114 dealii::Utilities::pow(n_subdivisions_1d, dim) * dealii::Utilities::pow(2ULL, l * dim);
115
116 if(n_cells >= n_cells_min and n_cells <= n_cells_max)
117 {
118 resolutions.push_back(
119 std::tuple<unsigned int, unsigned int, unsigned int>(degree, l, n_subdivisions_1d));
120
121 if(run_type == RunType::FixedProblemSize)
122 {
123 break;
124 }
125 }
126 }
127
128 // coarse grid with one cell and refine_level uniform refinements
129 {
130 unsigned int const n_subdivisions_1d = 1;
131 unsigned int const l = refine_level;
132
133 unsigned int const n_cells = dealii::Utilities::pow(2ULL, l * dim);
134
135 if(n_cells >= n_cells_min and n_cells <= n_cells_max)
136 {
137 resolutions.push_back(
138 std::tuple<unsigned int, unsigned int, unsigned int>(degree, l, n_subdivisions_1d));
139
140 if(run_type == RunType::FixedProblemSize)
141 {
142 break;
143 }
144 }
145 }
146
147 // coarse grid with 5^dim cells, and refine_level-2 uniform refinements
148 if(dim == 3 and refine_level >= 2)
149 {
150 unsigned int const n_subdivisions_1d = 5;
151 unsigned int const l = refine_level - 2;
152
153 unsigned int const n_cells =
154 dealii::Utilities::pow(n_subdivisions_1d, dim) * dealii::Utilities::pow(2ULL, l * dim);
155
156 if(n_cells >= n_cells_min and n_cells <= n_cells_max)
157 {
158 resolutions.push_back(
159 std::tuple<unsigned int, unsigned int, unsigned int>(degree, l, n_subdivisions_1d));
160
161 if(run_type == RunType::FixedProblemSize)
162 {
163 break;
164 }
165 }
166 }
167
168 // perform one global refinement
169 ++refine_level;
170 }
171
172 if(run_type == RunType::FixedProblemSize)
173 {
174 AssertThrow(resolutions.size() > resolutions_initial_size,
175 dealii::ExcMessage(
176 "No mesh found that meets the requirements regarding problem size. "
177 "Make sure that maximum number of dofs is sufficiently larger than "
178 "minimum number of dofs."));
179 }
180}
181
182struct HypercubeResolutionParameters
183{
184 HypercubeResolutionParameters()
185 {
186 }
187
188 HypercubeResolutionParameters(std::string const & input_file, unsigned int const dim) : dim(dim)
189 {
190 dealii::ParameterHandler prm;
191 add_parameters(prm);
192 prm.parse_input(input_file, "", true, true);
193
194 verify_parameters();
195 }
196
197 void
198 add_parameters(dealii::ParameterHandler & prm)
199 {
200 prm.enter_subsection("Resolution");
201 {
202 prm.add_parameter(
203 "RunType", run_type, "Type of throughput study.", Patterns::Enum<RunType>(), true);
204 prm.add_parameter(
205 "ElementType", element_type, "Type of elements.", Patterns::Enum<ElementType>(), true);
206 prm.add_parameter("DegreeMin",
207 degree_min,
208 "Minimal polynomial degree of shape functions.",
209 dealii::Patterns::Integer(1),
210 true);
211 prm.add_parameter("DegreeMax",
212 degree_max,
213 "Maximal polynomial degree of shape functions.",
214 dealii::Patterns::Integer(1),
215 true);
216 prm.add_parameter("RefineSpaceMin",
217 refine_space_min,
218 "Minimal number of mesh refinements.",
219 dealii::Patterns::Integer(0, 20),
220 true);
221 prm.add_parameter("RefineSpaceMax",
222 refine_space_max,
223 "Maximal number of mesh refinements.",
224 dealii::Patterns::Integer(0, 20),
225 true);
226 prm.add_parameter("DofsMin",
227 n_dofs_min,
228 "Minimal number of degrees of freedom.",
229 dealii::Patterns::Integer(1),
230 true);
231 prm.add_parameter("DofsMax",
232 n_dofs_max,
233 "Maximal number of degrees of freedom.",
234 dealii::Patterns::Integer(1),
235 true);
236 }
237 prm.leave_subsection();
238 }
239
240 void
241 verify_parameters()
242 {
243 if(run_type == RunType::RefineHAndP)
244 {
245 AssertThrow(degree_max >= degree_min, dealii::ExcMessage("Invalid parameters."));
246 AssertThrow(refine_space_max >= refine_space_min, dealii::ExcMessage("Invalid parameters."));
247 }
248 else if(run_type == RunType::FixedProblemSize)
249 {
250 AssertThrow(degree_max >= degree_min, dealii::ExcMessage("Invalid parameters."));
251 AssertThrow(n_dofs_max >= n_dofs_min, dealii::ExcMessage("Invalid parameters."));
252 }
253 else if(run_type == RunType::IncreasingProblemSize)
254 {
255 AssertThrow(
256 degree_min == degree_max,
257 dealii::ExcMessage(
258 "Only a single polynomial degree can be considered for RunType::IncreasingProblemSize"));
259
260 AssertThrow(n_dofs_max >= n_dofs_min, dealii::ExcMessage("Invalid parameters."));
261 }
262 else
263 {
264 AssertThrow(false, dealii::ExcMessage("not implemented."));
265 }
266 }
267
268 void
269 fill_resolution_vector(
270 std::function<double(unsigned int, unsigned int, ElementType)> const & get_dofs_per_element)
271 {
272 if(run_type == RunType::RefineHAndP)
273 {
274 unsigned int const n_cells_1d = 1;
275
276 // k-refinement
277 for(unsigned int degree = degree_min; degree <= degree_max; ++degree)
278 {
279 // h-refinement
280 for(unsigned int refine_space = refine_space_min; refine_space <= refine_space_max;
281 ++refine_space)
282 {
283 resolutions.push_back(
284 std::tuple<unsigned int, unsigned int, unsigned int>(degree, refine_space, n_cells_1d));
285 }
286 }
287 }
288 else if(run_type == RunType::FixedProblemSize or run_type == RunType::IncreasingProblemSize)
289 {
290 for(unsigned int degree = degree_min; degree <= degree_max; ++degree)
291 {
292 double dofs_per_element = get_dofs_per_element(dim, degree, element_type);
293 fill_resolutions_vector(resolutions,
294 dim,
295 degree,
296 dofs_per_element,
297 n_dofs_min,
298 n_dofs_max,
299 run_type,
300 element_type);
301 }
302 }
303 else
304 {
305 AssertThrow(false, dealii::ExcMessage("Not implemented."));
306 }
307 }
308
309 unsigned int dim = 2; // number of space dimensions
310
311 RunType run_type = RunType::RefineHAndP;
312
313 ElementType element_type = ElementType::Hypercube;
314
315 unsigned int degree_min = 3; // minimal polynomial degree
316 unsigned int degree_max = 3; // maximal polynomial degree
317
318 unsigned int refine_space_min = 0; // minimal number of global refinements
319 unsigned int refine_space_max = 0; // maximal number of global refinements
320
321 dealii::types::global_dof_index n_dofs_min = 1e4; // minimal number of unknowns
322 dealii::types::global_dof_index n_dofs_max = 3e4; // maximal number of unknowns
323
324 // a vector storing tuples of the form (degree k, refine level l, n_subdivisions_1d)
325 std::vector<std::tuple<unsigned int, unsigned int, unsigned int>> resolutions;
326};
327} // namespace ExaDG
328
329
330#endif /* INCLUDE_EXADG_OPERATORS_HYPERCUBE_RESOLUTION_PARAMETERS_H_ */
Definition driver.cpp:33