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 unsigned int 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 dealii::types::global_dof_index const n_cells_min =
64 (n_dofs_min + dofs_per_element - 1) / dofs_per_element;
65 dealii::types::global_dof_index const n_cells_max = n_dofs_max / dofs_per_element;
66
67 // Determine the number of cells per cube
68
69 // refers to dealii::GridTools::subdivided_hyper_cube() for Hypercube elements
70 unsigned int n_cells_per_cube = 1;
71
72 // refers to dealii::GridTools::subdivided_hyper_cube_with_simplices()
73 if(element_type == ElementType::Simplex)
74 {
75 if(dim == 2)
76 n_cells_per_cube = 2;
77 else if(dim == 3)
78 n_cells_per_cube = 5;
79 else
80 AssertThrow(false, dealii::ExcMessage("not implemented."));
81 }
82
83 // From the maximum number of cells, we derive a maximum refinement level for a uniformly refined
84 // mesh with n_cells_per_cube coarse-grid cells.
85 unsigned int const refine_level_max =
86 int(std::log((double)n_cells_max / (double)n_cells_per_cube) /
87 std::log((double)dealii::Utilities::pow(2ULL, dim))) +
88 1;
89
90 // we start with the coarsest possible mesh and then increase the refine level by 1
91 // until we hit the maximum refinement level
92 unsigned int refine_level = 0;
93
94 // This loop is for a uniformly refined mesh with n_cells_per_cube coarse-grid cells. Inside the
95 // loop, we test whether additional combinations of coarse-grid cells and mesh refinement levels
96 // are suitable in terms of n_cells_min/n_cells_max. The reason behind is that we want to have
97 // more data points than just data points differing by a factor of 2^dim in the number of cells.
98 while(refine_level <= refine_level_max)
99 {
100 // To obtain additional data points, we test coarse grids with {3,4}^dim elements in 2D and
101 // {3,4,5}^dim elements in 3D. The cases with 3 and 5 subdivisions per coordinate direction are
102 // only relevant for refine_level >= 2.
103
104 // coarse grid with 3^dim cells, and refine_level-2 uniform refinements
105 if(refine_level >= 2)
106 {
107 unsigned int const n_subdivisions_1d = 3;
108 unsigned int const l = refine_level - 2;
109
110 unsigned int const n_cells = n_cells_per_cube *
111 dealii::Utilities::pow(n_subdivisions_1d, dim) *
112 dealii::Utilities::pow(2ULL, l * dim);
113
114 if(n_cells >= n_cells_min and n_cells <= n_cells_max)
115 {
116 resolutions.push_back(
117 std::tuple<unsigned int, unsigned int, unsigned int>(degree, l, n_subdivisions_1d));
118
119 if(run_type == RunType::FixedProblemSize)
120 {
121 break;
122 }
123 }
124 }
125
126 // coarse grid with n_cells_per_cube cells, and refine_level uniform refinements
127 {
128 unsigned int const n_subdivisions_1d = 1;
129 unsigned int const l = refine_level;
130
131 unsigned int const n_cells = n_cells_per_cube * dealii::Utilities::pow(2ULL, l * dim);
132
133 if(n_cells >= n_cells_min and n_cells <= n_cells_max)
134 {
135 resolutions.push_back(
136 std::tuple<unsigned int, unsigned int, unsigned int>(degree, l, n_subdivisions_1d));
137
138 if(run_type == RunType::FixedProblemSize)
139 {
140 break;
141 }
142 }
143 }
144
145 // coarse grid with 5^dim cells, and refine_level-2 uniform refinements
146 if(dim == 3 and refine_level >= 2)
147 {
148 unsigned int const n_subdivisions_1d = 5;
149 unsigned int const l = refine_level - 2;
150
151 unsigned int const n_cells = n_cells_per_cube *
152 dealii::Utilities::pow(n_subdivisions_1d, dim) *
153 dealii::Utilities::pow(2ULL, l * dim);
154
155 if(n_cells >= n_cells_min and n_cells <= n_cells_max)
156 {
157 resolutions.push_back(
158 std::tuple<unsigned int, unsigned int, unsigned int>(degree, l, n_subdivisions_1d));
159
160 if(run_type == RunType::FixedProblemSize)
161 {
162 break;
163 }
164 }
165 }
166
167 // perform one global refinement
168 ++refine_level;
169 }
170
171 if(run_type == RunType::FixedProblemSize)
172 {
173 AssertThrow(resolutions.size() > resolutions_initial_size,
174 dealii::ExcMessage(
175 "No mesh found that meets the requirements regarding problem size. "
176 "Make sure that maximum number of dofs is sufficiently larger than "
177 "minimum number of dofs."));
178 }
179}
180
182{
184 {
185 }
186
187 HypercubeResolutionParameters(std::string const & input_file, unsigned int const dim) : dim(dim)
188 {
189 dealii::ParameterHandler prm;
190 add_parameters(prm);
191 prm.parse_input(input_file, "", true, true);
192
193 verify_parameters();
194 }
195
196 void
197 add_parameters(dealii::ParameterHandler & prm)
198 {
199 prm.enter_subsection("Resolution");
200 {
201 prm.add_parameter(
202 "RunType", run_type, "Type of throughput study.", Patterns::Enum<RunType>(), true);
203 prm.add_parameter(
204 "ElementType", element_type, "Type of elements.", Patterns::Enum<ElementType>(), true);
205 prm.add_parameter("DegreeMin",
206 degree_min,
207 "Minimal polynomial degree of shape functions.",
208 dealii::Patterns::Integer(1),
209 true);
210 prm.add_parameter("DegreeMax",
211 degree_max,
212 "Maximal polynomial degree of shape functions.",
213 dealii::Patterns::Integer(1),
214 true);
215 prm.add_parameter("RefineSpaceMin",
216 refine_space_min,
217 "Minimal number of mesh refinements.",
218 dealii::Patterns::Integer(0, 20),
219 true);
220 prm.add_parameter("RefineSpaceMax",
221 refine_space_max,
222 "Maximal number of mesh refinements.",
223 dealii::Patterns::Integer(0, 20),
224 true);
225 prm.add_parameter("DofsMin",
226 n_dofs_min,
227 "Minimal number of degrees of freedom.",
228 dealii::Patterns::Integer(1),
229 true);
230 prm.add_parameter("DofsMax",
231 n_dofs_max,
232 "Maximal number of degrees of freedom.",
233 dealii::Patterns::Integer(1),
234 true);
235 }
236 prm.leave_subsection();
237 }
238
239 void
240 verify_parameters()
241 {
242 if(run_type == RunType::RefineHAndP)
243 {
244 AssertThrow(degree_max >= degree_min, dealii::ExcMessage("Invalid parameters."));
245 AssertThrow(refine_space_max >= refine_space_min, dealii::ExcMessage("Invalid parameters."));
246 }
247 else if(run_type == RunType::FixedProblemSize)
248 {
249 AssertThrow(degree_max >= degree_min, dealii::ExcMessage("Invalid parameters."));
250 AssertThrow(n_dofs_max >= n_dofs_min, dealii::ExcMessage("Invalid parameters."));
251 }
252 else if(run_type == RunType::IncreasingProblemSize)
253 {
254 AssertThrow(
255 degree_min == degree_max,
256 dealii::ExcMessage(
257 "Only a single polynomial degree can be considered for RunType::IncreasingProblemSize"));
258
259 AssertThrow(n_dofs_max >= n_dofs_min, dealii::ExcMessage("Invalid parameters."));
260 }
261 else
262 {
263 AssertThrow(false, dealii::ExcMessage("not implemented."));
264 }
265 }
266
267 void
268 fill_resolution_vector(
269 std::function<unsigned int(unsigned int, unsigned int, ElementType)> const &
270 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 unsigned int 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
Definition hypercube_resolution_parameters.h:182