ExaDG
Loading...
Searching...
No Matches
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_CONVECTION_DIFFUSION_INPUT_PARAMETERS_H_
23#define INCLUDE_CONVECTION_DIFFUSION_INPUT_PARAMETERS_H_
24
25// deal.II
26#include <deal.II/base/exceptions.h>
27
28// ExaDG
29#include <exadg/convection_diffusion/user_interface/enum_types.h>
30#include <exadg/grid/grid_data.h>
31#include <exadg/operators/adaptive_mesh_refinement.h>
32#include <exadg/operators/inverse_mass_parameters.h>
33#include <exadg/solvers_and_preconditioners/multigrid/multigrid_parameters.h>
34#include <exadg/solvers_and_preconditioners/preconditioners/enum_types.h>
35#include <exadg/solvers_and_preconditioners/solvers/enum_types.h>
36#include <exadg/solvers_and_preconditioners/solvers/solver_data.h>
37#include <exadg/time_integration/enum_types.h>
38#include <exadg/time_integration/restart_data.h>
39#include <exadg/time_integration/solver_info_data.h>
40
41namespace ExaDG
42{
43namespace ConvDiff
44{
46{
47public:
48 // standard constructor that initializes parameters with default values
49 Parameters();
50
51 // check correctness of parameters
52 void
53 check() const;
54
55 bool
56 involves_h_multigrid() const;
57
58 bool
59 linear_system_including_convective_term_has_to_be_solved() const;
60
61 bool
62 convective_problem() const;
63
64 bool
65 diffusive_problem() const;
66
67 bool
68 linear_system_has_to_be_solved() const;
69
70 TypeVelocityField
71 get_type_velocity_field() const;
72
73 // print functions
74 void
75 print(dealii::ConditionalOStream const & pcout, std::string const & name) const;
76
77private:
78 void
79 print_parameters_mathematical_model(dealii::ConditionalOStream const & pcout) const;
80
81 void
82 print_parameters_physical_quantities(dealii::ConditionalOStream const & pcout) const;
83
84 void
85 print_parameters_temporal_discretization(dealii::ConditionalOStream const & pcout) const;
86
87 void
88 print_parameters_spatial_discretization(dealii::ConditionalOStream const & pcout) const;
89
90 void
91 print_parameters_solver(dealii::ConditionalOStream const & pcout) const;
92
93 void
94 print_parameters_numerical_parameters(dealii::ConditionalOStream const & pcout) const;
95
96public:
97 /**************************************************************************************/
98 /* */
99 /* MATHEMATICAL MODEL */
100 /* */
101 /**************************************************************************************/
102
103 // description: see enum declaration
104 ProblemType problem_type;
105
106 // description: see enum declaration
107 EquationType equation_type;
108
109 // Use true if an analytical function is used to prescribe the velocity field
110 bool analytical_velocity_field;
111
112 // Use Arbitrary Lagrangian-Eulerian (ALE) formulation
113 bool ale_formulation;
114
115 // set right_hand_side = true if the right-hand side f is unequal zero
116 bool right_hand_side;
117
118 // type of formulation of convective term
119 FormulationConvectiveTerm formulation_convective_term;
120
121 /**************************************************************************************/
122 /* */
123 /* PHYSICAL QUANTITIES */
124 /* */
125 /**************************************************************************************/
126
127 // start time of simulation
128 double start_time;
129
130 // end time of simulation
131 double end_time;
132
133 // kinematic diffusivity
134 double diffusivity;
135
136
137
138 /**************************************************************************************/
139 /* */
140 /* TEMPORAL DISCRETIZATION */
141 /* */
142 /**************************************************************************************/
143
144 // temporal discretization method
145 TemporalDiscretization temporal_discretization;
146
147 // description: see enum declaration (only relevant for explicit time integration)
148 TimeIntegratorRK time_integrator_rk;
149
150 // order of time integration scheme (only relevant for BDF time integration)
151 unsigned int order_time_integrator;
152
153 // start with low order (only relevant for BDF time integration)
154 bool start_with_low_order;
155
156 // description: see enum declaration (this parameter is ignored for steady problems or
157 // unsteady problems with explicit Runge-Kutta time integration scheme). In case of
158 // a purely diffusive problem, one also does not have to specify this parameter.
159 TreatmentOfConvectiveTerm treatment_of_convective_term;
160
161 // calculation of time step size
162 TimeStepCalculation calculation_of_time_step_size;
163
164 // use adaptive time stepping?
165 bool adaptive_time_stepping;
166
167 // This parameter defines by which factor the time step size is allowed to increase
168 // or to decrease in case of adaptive time step, e.g., if one wants to avoid large
169 // jumps in the time step size. A factor of 1 implies that the time step size can not
170 // change at all, while a factor towards infinity implies that arbitrary changes in
171 // the time step size are allowed from one time step to the next.
172 double adaptive_time_stepping_limiting_factor;
173
174 // specify a maximum time step size in case of adaptive time stepping since the adaptive
175 // time stepping algorithm would choose arbitrarily large time step sizes if the velocity field
176 // is close to zero. This variable is only used for adaptive time stepping.
177 double time_step_size_max;
178
179 // Different variants are available for calculating the time step size based on a local CFL
180 // criterion.
181 CFLConditionType adaptive_time_stepping_cfl_type;
182
183 // user specified time step size: note that this time_step_size is the first
184 // in a series of time_step_size's when performing temporal convergence tests,
185 // i.e., delta_t = time_step_size, time_step_size/2, ...
186 double time_step_size;
187
188 // maximum number of time steps
189 unsigned int max_number_of_time_steps;
190
191 // number of refinements for temporal discretization
192 unsigned int n_refine_time;
193
194 // cfl number ("global" CFL number, can be larger than critical CFL in case
195 // of operator-integration-factor splitting)
196 double cfl;
197
198 // estimation of maximum velocity required for CFL condition
199 double max_velocity;
200
201 // diffusion number (relevant number for limitation of time step size
202 // when treating the diffusive term explicitly)
203 double diffusion_number;
204
205 // C_eff: constant that has to be specified for time step calculation method
206 // MaxEfficiency, which means that the time step is selected such that the errors of
207 // the temporal and spatial discretization are comparable
208 double c_eff;
209
210 // exponent of fe_degree used in the calculation of the convective time step size
211 double exponent_fe_degree_convection;
212
213 // exponent of fe_degree used in the calculation of the diffusion time step size
214 double exponent_fe_degree_diffusion;
215
216 // set this variable to true to start the simulation from restart files
217 bool restarted_simulation;
218
219 // Restart
220 RestartData restart_data;
221
222 /**************************************************************************************/
223 /* */
224 /* SPATIAL DISCRETIZATION */
225 /* */
226 /**************************************************************************************/
227
228 // Grid data
229 GridData grid;
230
231 // Mapping
232 unsigned int mapping_degree;
233
234 // mapping degree for coarser grids in h-multigrid
235 unsigned int mapping_degree_coarse_grids;
236
237 // polynomial degree of shape functions
238 unsigned int degree;
239
240 // enable adaptive mesh refinement
241 bool enable_adaptivity;
243
244 // description: see enum declaration
245 NumericalFluxConvectiveOperator numerical_flux_convective_operator;
246
247 // diffusive term: Symmetric interior penalty discretization Galerkin (SIPG)
248 // interior penalty parameter scaling factor: default value is 1.0
249 double IP_factor;
250
251 /**************************************************************************************/
252 /* */
253 /* Solver parameters for mass matrix problem */
254 /* */
255 /**************************************************************************************/
256 // These parameters are only relevant if the inverse mass can not be realized as a
257 // matrix-free operator evaluation. The typical use case is a DG formulation with non
258 // hypercube elements (e.g. simplex elements).
259
260 InverseMassParameters inverse_mass_operator;
261
262 /**************************************************************************************/
263 /* */
264 /* InverseMassPreconditioner */
265 /* */
266 /**************************************************************************************/
267
268 // This parameter is used for the inverse mass preconditioner. Note that there is a
269 // separate parameter above for the inverse mass operator (which needs to be inverted
270 // "exactly" for reasons of accuracy).
271 // This parameter is only relevant if the inversion of the block-diagonal mass operator
272 // can not be realized as a matrix-free operator evaluation. The typical use case is a
273 // DG formulation with non hypercube elements (e.g. simplex). In this case, however,
274 // you should think about replacing the inverse mass preconditioner by a simple point
275 // Jacobi preconditioner as an efficient alternative to the (exact) inverse mass
276 // preconditioner.
277
278 InverseMassParameters inverse_mass_preconditioner;
279
280 /**************************************************************************************/
281 /* */
282 /* SOLVER */
283 /* */
284 /**************************************************************************************/
285
286 // description: see enum declaration
287 Solver solver;
288
289 // solver data
290 SolverData solver_data;
291
292 // description: see enum declaration
293 Preconditioner preconditioner;
294
295 // update preconditioner in case of varying parameters
296 bool update_preconditioner;
297
298 // update preconditioner every ... time step. Only relevant if update preconditioner
299 // is set to true.
300 unsigned int update_preconditioner_every_time_steps;
301
302 // Implement block diagonal (block Jacobi) preconditioner in a matrix-free way
303 // by solving the block Jacobi problems elementwise using iterative solvers and
304 // matrix-free operator evaluation
305 bool implement_block_diagonal_preconditioner_matrix_free;
306
307 // description: see enum declaration
308 Elementwise::Solver solver_block_diagonal;
309
310 // description: see enum declaration
311 Elementwise::Preconditioner preconditioner_block_diagonal;
312
313 // solver data for block Jacobi preconditioner (only relevant for elementwise
314 // iterative solution procedure)
315 SolverData solver_data_block_diagonal;
316
317 // description: see enum declaration
318 MultigridOperatorType mg_operator_type;
319
320 // description: see declaration of MultigridData
321 MultigridData multigrid_data;
322
323 // show solver performance (wall time, number of iterations)
324 SolverInfoData solver_info_data;
325
326 /**************************************************************************************/
327 /* */
328 /* NUMERICAL PARAMETERS */
329 /* */
330 /**************************************************************************************/
331
332 // By default, the matrix-free implementation performs separate loops over all cells,
333 // interior faces, and boundary faces. For a certain type of operations, however, it
334 // is necessary to perform the face-loop as a loop over all faces of a cell with an
335 // outer loop over all cells, e.g., preconditioners operating on the level of
336 // individual cells (for example block Jacobi). With this parameter, the loop structure
337 // can be changed to such an algorithm (cell_based_face_loops).
338 bool use_cell_based_face_loops;
339
340 // Evaluate convective term and diffusive term at once instead of implementing each
341 // operator separately and subsequently looping over all operators. This parameter is
342 // only relevant in case of fully explicit time stepping. In case of semi-implicit or
343 // fully implicit time integration the combined operator will always be used.
344 bool use_combined_operator;
345
346 // In case that the velocity field is prescribed analytically, it might be advantageous
347 // from the point of view of computational costs to store the velocity field in a DoF
348 // vector instead of repeatedly calling dealii::Function<dim>::value() whenever evaluating the
349 // operator in iterative solvers. In other words, depending on the computer hardware it
350 // might be more efficient to load a DoF vector and interpolate into the quadrature points
351 // than calculating the velocity field by means of dealii::Function<dim>::value().
352 // This strategy makes only sense in case of steady state problems or unsteady problems
353 // with an implicit treatment of the convective term, i.e., in cases where the convective
354 // term has to be evaluated more than once at a given time t.
355 bool store_analytical_velocity_in_dof_vector;
356
357 // use 3/2 overintegration rule for convective term
358 bool use_overintegration;
359};
360
361} // namespace ConvDiff
362} // namespace ExaDG
363
364#endif /* INCLUDE_CONVECTION_DIFFUSION_INPUT_PARAMETERS_H_ */
Definition parameters.h:46
Definition driver.cpp:33
Definition adaptive_mesh_refinement.h:44
Definition grid_data.h:88
Definition inverse_mass_parameters.h:49
Definition multigrid_parameters.h:259
Definition restart_data.h:38
Definition solver_data.h:34
Definition solver_info_data.h:38