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_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_USER_INTERFACE_INPUT_PARAMETERS_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_USER_INTERFACE_INPUT_PARAMETERS_H_
24
25// deal.II
26#include <deal.II/base/conditional_ostream.h>
27
28// ExaDG
29#include <exadg/grid/grid_data.h>
30#include <exadg/incompressible_navier_stokes/user_interface/enum_types.h>
31#include <exadg/incompressible_navier_stokes/user_interface/viscosity_model_data.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/newton/newton_solver_data.h>
35#include <exadg/solvers_and_preconditioners/preconditioners/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 IncNS
44{
45class Parameters
46{
47public:
48 // standard constructor that initializes parameters
49 Parameters();
50
51 void
52 check(dealii::ConditionalOStream const & pcout) const;
53
54 bool
55 convective_problem() const;
56
57 bool
58 viscous_problem() const;
59
60 bool
61 viscosity_is_variable() const;
62
63 bool
64 non_explicit_convective_problem() const;
65
66 bool
67 implicit_nonlinear_convective_problem() const;
68
69 bool
70 nonlinear_viscous_problem() const;
71
72 bool
73 nonlinear_problem_has_to_be_solved() const;
74
75 bool
76 involves_h_multigrid() const;
77
78 unsigned int
79 get_degree_p(unsigned int const degree_u) const;
80
81 void
82 print(dealii::ConditionalOStream const & pcout, std::string const & name) const;
83
84private:
85 void
86 print_parameters_mathematical_model(dealii::ConditionalOStream const & pcout) const;
87
88 void
89 print_parameters_physical_quantities(dealii::ConditionalOStream const & pcout) const;
90
91 void
92 print_parameters_temporal_discretization(dealii::ConditionalOStream const & pcout) const;
93
94 void
95 print_parameters_spatial_discretization(dealii::ConditionalOStream const & pcout) const;
96
97 void
98 print_parameters_numerical_parameters(dealii::ConditionalOStream const & pcout) const;
99
100 void
101 print_parameters_pressure_poisson(dealii::ConditionalOStream const & pcout) const;
102
103 void
104 print_parameters_projection_step(dealii::ConditionalOStream const & pcout) const;
105
106 void
107 print_parameters_momentum_step(dealii::ConditionalOStream const & pcout) const;
108
109 void
110 print_parameters_dual_splitting(dealii::ConditionalOStream const & pcout) const;
111
112 void
113 print_parameters_pressure_correction(dealii::ConditionalOStream const & pcout) const;
114
115 void
116 print_parameters_coupled_solver(dealii::ConditionalOStream const & pcout) const;
117
118 // coupled solver
119 bool
120 involves_h_multigrid_velocity_block() const;
121
122 bool
123 involves_h_multigrid_pressure_block() const;
124
125
126 // penalty step (coupled solver or projection methods)
127 bool
128 involves_h_multigrid_penalty_step() const;
129
130 // projection methods
131 bool
132 involves_h_multigrid_pressure_step() const;
133
134 bool
135 involves_h_multigrid_momentum_step() const;
136
137public:
138 /**************************************************************************************/
139 /* */
140 /* MATHEMATICAL MODEL */
141 /* */
142 /**************************************************************************************/
143
144 // description: see enum declaration
145 ProblemType problem_type;
146
147 // description: see enum declaration
148 EquationType equation_type;
149
150 // description: see enum declaration
151 FormulationViscousTerm formulation_viscous_term;
152
153 // description: see enum declaration
154 FormulationConvectiveTerm formulation_convective_term;
155
156 // use stable outflow boundary condition for convective term according to
157 // Gravemeier et al. (2012)
158 bool use_outflow_bc_convective_term;
159
160 // if the body force vector on the right-hand side of the momentum equation of the
161 // Navier-Stokes equations is unequal zero, set right_hand_side = true
162 // This parameter also has to be true when considering the Boussinesq term.
163 bool right_hand_side;
164
165 // Boussinesq term (natural convection through buoyancy forces)
166 bool boussinesq_term;
167
168 // If Boussinesq term is activated: solves only for dynamic pressure variations if true,
169 // and includes hydrostatic component if false.
170 bool boussinesq_dynamic_part_only;
171
172 /**************************************************************************************/
173 /* */
174 /* Arbitrary Lagrangian-Eulerian formulation (ALE) */
175 /* */
176 /**************************************************************************************/
177
178 // use true to activate ALE formulation, otherwise the standard Eulerian formulation
179 // with fixed mesh will be used
180 bool ale_formulation;
181
182 MeshMovementType mesh_movement_type;
183
184 bool neumann_with_variable_normal_vector;
185
186 /**************************************************************************************/
187 /* */
188 /* PHYSICAL QUANTITIES */
189 /* */
190 /**************************************************************************************/
191
192 // start time of simulation
193 double start_time;
194
195 // end time of simulation
196 double end_time;
197
198 // kinematic viscosity
199 double viscosity;
200
201 // density (not required by fluid solver which is formulated in terms of the kinematic
202 // viscosity only, but for the calculation of the fluid stress for FSI problems)
203 double density;
204
205 // Boussinesq term
206 double thermal_expansion_coefficient;
207 double reference_temperature;
208
209 /**************************************************************************************/
210 /* */
211 /* TEMPORAL DISCRETIZATION */
212 /* */
213 /**************************************************************************************/
214
215 // description: see enum declaration
216 SolverType solver_type;
217
218 // description: see enum declaration
219 TemporalDiscretization temporal_discretization;
220
221 // description: see enum declaration
222 TreatmentOfConvectiveTerm treatment_of_convective_term;
223
224 // description: see enum declaration
225 TimeStepCalculation calculation_of_time_step_size;
226
227 // use adaptive time stepping?
228 bool adaptive_time_stepping;
229
230 // This parameter defines by which factor the time step size is allowed to increase
231 // or to decrease in case of adaptive time step, e.g., if one wants to avoid large
232 // jumps in the time step size. A factor of 1 implies that the time step size can not
233 // change at all, while a factor towards infinity implies that arbitrary changes in
234 // the time step size are allowed from one time step to the next.
235 double adaptive_time_stepping_limiting_factor;
236
237 // specify a maximum time step size in case of adaptive time stepping since the adaptive
238 // time stepping algorithm would choose arbitrarily large time step sizes if the velocity field
239 // is close to zero. This variable is only used for adaptive time stepping.
240 double time_step_size_max;
241
242 // Different variants are available for calculating the time step size based on a local CFL
243 // criterion.
244 CFLConditionType adaptive_time_stepping_cfl_type;
245
246 // maximum velocity needed when calculating the time step according to cfl-condition
247 double max_velocity;
248
249 // cfl number: note that this cfl number is the first in a series of cfl numbers
250 // when performing temporal convergence tests, i.e., cfl_real = cfl, cfl/2, cfl/4, ...
251 // ("global" CFL number, can be larger than critical CFL in case
252 // of operator-integration-factor splitting)
253 double cfl;
254
255 // dt = CFL/k_u^{exp} * h / || u ||
256 double cfl_exponent_fe_degree_velocity;
257
258 // C_eff: constant that has to be specified for time step calculation method
259 // MaxEfficiency, which means that the time step is selected such that the errors of
260 // the temporal and spatial discretization are comparable
261 double c_eff;
262
263 // user specified time step size: note that this time_step_size is the first
264 // in a series of time_step_size's when performing temporal convergence tests,
265 // i.e., delta_t = time_step_size, time_step_size/2, ...
266 double time_step_size;
267
268 // maximum number of time steps
269 unsigned int max_number_of_time_steps;
270
271 // number of refinements for temporal discretization
272 unsigned int n_refine_time;
273
274 // order of BDF time integration scheme and extrapolation scheme
275 unsigned int order_time_integrator;
276
277 // start time integrator with low order time integrator, i.e., first order Euler method
278 bool start_with_low_order;
279
280 // description: see enum declaration
281 ConvergenceCriterionSteadyProblem convergence_criterion_steady_problem;
282
283 // Pseudo-timestepping for steady-state problems: These tolerances are only relevant
284 // when using an unsteady solver to solve the steady Navier-Stokes equations.
285 //
286 // option ResidualNavierStokes:
287 // - these tolerances refer to the norm of the residual of the steady
288 // Navier-Stokes equations.
289 //
290 // option SolutionIncrement:
291 // - these tolerances refer to the norm of the increment of the solution
292 // vector from one time step to the next.
293 double abs_tol_steady;
294 double rel_tol_steady;
295
296 // show solver performance (wall time, number of iterations) every ... timesteps
297 SolverInfoData solver_info_data;
298
299 // set this variable to true to start the simulation from restart files
300 bool restarted_simulation;
301
302 // restart
303 RestartData restart_data;
304
305
306 /**************************************************************************************/
307 /* */
308 /* SPATIAL DISCRETIZATION */
309 /* */
310 /**************************************************************************************/
311
312 // Grid data
313 GridData grid;
314
315 // Mapping
316 unsigned int mapping_degree;
317
318 // mapping degree for coarser grids in h-multigrid
319 unsigned int mapping_degree_coarse_grids;
320
321 // type of spatial discretization approach
322 SpatialDiscretization spatial_discretization;
323
324 // Polynomial degree of velocity shape functions. In case of H-div Raviart-Thomas 'degree_u' is
325 // the polynomial degree in normal direction and (degree_u - 1) is the degree in tangential
326 // direction.
327 unsigned int degree_u;
328
329 // Polynomial degree of pressure shape functions
330 DegreePressure degree_p;
331
332 // convective term: upwind factor describes the scaling factor in front of the
333 // stabilization term (which is strictly dissipative) of the numerical function
334 // of the convective term. For the divergence formulation of the convective term with
335 // local Lax-Friedrichs flux, a value of upwind_factor = 1.0 corresponds to the
336 // theoretical value (e.g., maximum eigenvalue of the flux Jacobian, lambda = 2 |u*n|)
337 // but a lower value (e.g., upwind_factor = 0.5, lambda = |u*n|) might be much more
338 // advantages in terms of computational costs by allowing significantly larger time
339 // step sizes.
340 double upwind_factor;
341
342 // description: see enum declaration
343 TypeDirichletBCs type_dirichlet_bc_convective;
344
345 // description: see enum declaration
346 InteriorPenaltyFormulation IP_formulation_viscous;
347
348 // description: see enum declaration
349 PenaltyTermDivergenceFormulation penalty_term_div_formulation;
350
351 // interior penalty parameter scaling factor for Helmholtz equation of viscous step
352 double IP_factor_viscous;
353
354 // integration by parts of grad(P)
355 bool gradp_integrated_by_parts;
356
357 // type for formulation
358 FormulationPressureGradientTerm gradp_formulation;
359
360 // use boundary data if integrated by parts
361 bool gradp_use_boundary_data;
362
363 // integration by parts of div(U)
364 bool divu_integrated_by_parts;
365
366 // type of formulation
367 FormulationVelocityDivergenceTerm divu_formulation;
368
369 // use boundary data if integrated by parts
370 bool divu_use_boundary_data;
371
372 // For certain setups and types of boundary conditions, the pressure level is undefined,
373 // e.g., if the velocity is prescribed on the whole boundary or in case of periodic
374 // boundary conditions. This variable defines the method used to adjust the pressure level
375 // in or to obtain a well-defined pressure solution.
376 // If you observe convergence in the velocity error, but no convergence in the pressure
377 // error, this might likely be related to the fact that this parameter has not been set
378 // as intended.
379 AdjustPressureLevel adjust_pressure_level;
380
381 // use div-div penalty term
382 bool use_divergence_penalty;
383
384 // penalty factor of divergence penalty term
385 double divergence_penalty_factor;
386
387 // use continuity penalty term
388 bool use_continuity_penalty;
389
390 // penalty factor of continuity penalty term
391 double continuity_penalty_factor;
392
393 // Divergence and continuity penalty terms are applied in a postprocessing step
394 // if set to true.
395 // Otherwise, penalty terms are added to the monolithic systems of equations in case of
396 // the monolithic solver, or to the projection step in case of the dual splitting scheme.
397 // For the pressure-correction scheme, this parameter is irrelevant since the projection
398 // step is the last step of the splitting scheme anyway.
399 bool apply_penalty_terms_in_postprocessing_step;
400
401 // specify which components of the velocity field to penalize, i.e., normal
402 // components only or all components
403 ContinuityPenaltyComponents continuity_penalty_components;
404
405 // specify whether boundary conditions prescribed for the velocity should be
406 // used in continuity penalty operator. Otherwise, boundary faces are ignored.
407 bool continuity_penalty_use_boundary_data;
408
409 // type of penalty parameter (see enum declaration for more information)
410 TypePenaltyParameter type_penalty_parameter;
411
412 /**************************************************************************************/
413 /* */
414 /* Variable viscosity models */
415 /* */
416 /**************************************************************************************/
417
418 TreatmentOfVariableViscosity treatment_of_variable_viscosity;
419 TurbulenceModelData turbulence_model_data;
420 GeneralizedNewtonianModelData generalized_newtonian_model_data;
421
422 /**************************************************************************************/
423 /* */
424 /* NUMERICAL PARAMETERS */
425 /* */
426 /**************************************************************************************/
427
428 // Implement block diagonal (block Jacobi) preconditioner in a matrix-free way
429 // by solving the block Jacobi problems elementwise using iterative solvers and
430 // matrix-free operator evaluation. By default, this variable should be set to true
431 // because the matrix-based variant (which is used otherwise) is very slow and the
432 // matrix-free variant can be expected to be much faster.
433 // Only in case that convergence problems occur or for reasons of testing/debugging
434 // the matrix-based variant should be used.
435 bool implement_block_diagonal_preconditioner_matrix_free;
436
437 // By default, the matrix-free implementation performs separate loops over all cells,
438 // interior faces, and boundary faces. For a certain type of operations, however, it
439 // is necessary to perform the face-loop as a loop over all faces of a cell with an
440 // outer loop over all cells, e.g., preconditioners operating on the level of
441 // individual cells (for example block Jacobi). With this parameter, the loop structure
442 // can be changed to such an algorithm (cell_based_face_loops).
443 bool use_cell_based_face_loops;
444
445 // Solver data for block Jacobi preconditioner. Accordingly, this parameter is only
446 // relevant if the block diagonal preconditioner is implemented in a matrix-free way
447 // using an elementwise iterative solution procedure for which solver tolerances have to
448 // be provided by the user. It was found that rather coarse relative solver tolerances
449 // of around 1.e-2 are enough and for lower tolerances one would 'over-solve' the local
450 // preconditioning problems without a further benefit in global iteration counts.
451 SolverData solver_data_block_diagonal;
452
453 // Quadrature rule used to integrate the linearized convective term. This parameter is
454 // therefore only relevant if linear systems of equations have to be solved involving
455 // the convective term. The quadrature rule specified by this parameter is also used in
456 // case of a linearly implicit formulation of the convective term.
457 // Note that if the convective term is not involved in the momentum operator, this
458 // parameter does not have any effect and the standard quadrature rule will be used
459 // for the momentum operator.
460 // For reasons of computational efficiency, it might be advantageous to use a standard
461 // quadrature rule for the linear(ized) momentum operator in order to speed up
462 // the computation. However, it was found that choosing a lower order quadrature rule
463 // for the linearized problem only, increases the number of iterations significantly. It
464 // was found that the quadrature rules used for the nonlinear and linear problems should
465 // be the same. Hence, although this parameter speeds up the operator evaluation (i.e.
466 // the wall time per iteration), it is unclear whether a lower order quadrature rule
467 // really allows to achieve a more efficient method overall.
468 QuadratureRuleLinearization quad_rule_linearization;
469
470 /**************************************************************************************/
471 /* */
472 /* Solver parameters for mass matrix problem */
473 /* */
474 /**************************************************************************************/
475 // These parameters are relevant if the inverse mass can not be realized as a
476 // matrix-free operator evaluation. The typical use case is a DG formulation with non
477 // hypercube elements (e.g. simplex).
478
479 InverseMassParameters inverse_mass_operator;
480
481 // This parameter is only relevant if an H(div)-conforming formulation is chosen.
482 InverseMassParametersHdiv inverse_mass_operator_hdiv;
483
484 /**************************************************************************************/
485 /* */
486 /* InverseMassPreconditioner */
487 /* */
488 /**************************************************************************************/
489
490 // This parameter is used for all inverse mass preconditioners used for incompressible
491 // Navier-Stokes solvers. Note that there is a separate parameter above for the inverse
492 // mass operator (which needs to be inverted "exactly" for reasons of accuracy).
493 // This parameter is only relevant if the mass operator is block diagonal and if the
494 // inverse mass operator can not be realized as a matrix-free operator evaluation. The
495 // typical use case is a DG formulation with non hypercube elements (e.g. simplex). In
496 // these cases, however, you should think about replacing the inverse mass preconditioner
497 // by a simple point Jacobi preconditioner as an efficient alternative to the (exact)
498 // inverse mass preconditioner. This type of preconditioner is not available for an
499 // H(div)-conforming formulation.
500
501 InverseMassParameters inverse_mass_preconditioner;
502
503 /**************************************************************************************/
504 /* */
505 /* PROJECTION METHODS */
506 /* */
507 /**************************************************************************************/
508
509 // PRESSURE POISSON EQUATION
510
511 // interior penalty parameter scaling factor for pressure Poisson equation
512 double IP_factor_pressure;
513
514 // description: see enum declaration
515 SolverPressurePoisson solver_pressure_poisson;
516
517 // solver data for pressure Poisson equation
518 SolverData solver_data_pressure_poisson;
519
520 // description: see enum declaration
521 PreconditionerPressurePoisson preconditioner_pressure_poisson;
522
523 // update of preconditioner for this equation is currently not provided and not needed
524
525 // description: see declaration of MultigridData
526 MultigridData multigrid_data_pressure_poisson;
527
528 // Update preconditioner before solving the linear system of equations.
529 bool update_preconditioner_pressure_poisson;
530
531 // Update preconditioner every ... time steps.
532 // This variable is only used if update of preconditioner is true.
533 unsigned int update_preconditioner_pressure_poisson_every_time_steps;
534
535 // PROJECTION STEP
536
537 // description: see enum declaration
538 SolverProjection solver_projection;
539
540 // solver data for projection step
541 SolverData solver_data_projection;
542
543 // description: see enum declaration
544 PreconditionerProjection preconditioner_projection;
545
546 // description: see declaration of MultigridData
547 MultigridData multigrid_data_projection;
548
549 // Update preconditioner before solving the linear system of equations.
550 // Note that this variable is only used when using an iterative method
551 // to solve the global projection equation.
552 bool update_preconditioner_projection;
553
554 // Update preconditioner every ... time steps.
555 // This variable is only used if update of preconditioner is true.
556 unsigned int update_preconditioner_projection_every_time_steps;
557
558 // description: see enum declaration (only relevant if block diagonal is used as
559 // preconditioner)
560 Elementwise::Preconditioner preconditioner_block_diagonal_projection;
561
562 // solver data for block Jacobi preconditioner (only relevant if elementwise
563 // iterative solution procedure is used for block diagonal preconditioner)
564 SolverData solver_data_block_diagonal_projection;
565
566 // MOMENTUM STEP
567
568 // Newton solver data
569 Newton::SolverData newton_solver_data_momentum;
570
571 // description: see enum declaration
572 SolverMomentum solver_momentum;
573
574 // Solver data for (linearized) momentum equation
575 SolverData solver_data_momentum;
576
577 // description: see enum declaration
578 MomentumPreconditioner preconditioner_momentum;
579
580 // update preconditioner before solving the linear system of equations
581 // only necessary if the operator changes during the simulation
582 bool update_preconditioner_momentum;
583
584 // Update preconditioner every ... Newton iterations (only relevant for
585 // nonlinear problems, i.e., if the convective term is formulated implicitly)
586 // This variable is only used if update_preconditioner_coupled = true.
587 unsigned int update_preconditioner_momentum_every_newton_iter;
588
589 // Update preconditioner every ... time steps.
590 // This variable is only used if update_preconditioner_coupled = true.
591 unsigned int update_preconditioner_momentum_every_time_steps;
592
593 // description: see declaration of MultigridData
594 MultigridData multigrid_data_momentum;
595
596 // description: see enum declaration
597 MultigridOperatorType multigrid_operator_type_momentum;
598
599 /**************************************************************************************/
600 /* */
601 /* HIGH-ORDER DUAL SPLITTING SCHEME */
602 /* */
603 /**************************************************************************************/
604
605 // FORMULATIONS
606
607 // order of extrapolation of viscous term and convective term in pressure Neumann BC
608 unsigned int order_extrapolation_pressure_nbc;
609
610 // description: see enum declaration
611 // The formulation of the convective term in the boundary conditions for the dual
612 // splitting scheme (there are two occurrences, the g_u_hat term
613 // arising from the divergence term on the right-hand side of the pressure Poisson
614 // equation as well as the pressure NBC for the dual splitting scheme)can be chosen
615 // independently from the type of formulation used for the discretization of the
616 // convective term in the Navier-Stokes equations.
617 // As a default parameter, FormulationConvectiveTerm::ConvectiveFormulation is
618 // used (exploiting that div(u)=0 holds in the continuous case).
619 FormulationConvectiveTerm formulation_convective_term_bc;
620
621 /**************************************************************************************/
622 /* */
623 /* PRESSURE-CORRECTION SCHEME */
624 /* */
625 /**************************************************************************************/
626
627 // order of pressure extrapolation in case of incremental formulation
628 // a value of 0 corresponds to non-incremental formulation
629 // and a value >=1 to incremental formulation
630 unsigned int order_pressure_extrapolation;
631
632 // rotational formulation
633 bool rotational_formulation;
634
635
636 /**************************************************************************************/
637 /* */
638 /* COUPLED NAVIER-STOKES SOLVER */
639 /* */
640 /**************************************************************************************/
641
642 // use a scaling of the continuity equation
643 bool use_scaling_continuity;
644
645 // scaling factor continuity equation
646 double scaling_factor_continuity;
647
648 // solver tolerances Newton solver
649 Newton::SolverData newton_solver_data_coupled;
650
651 // description: see enum declaration
652 SolverCoupled solver_coupled;
653
654 // Solver data for coupled solver
655 SolverData solver_data_coupled;
656
657 // description: see enum declaration
658 PreconditionerCoupled preconditioner_coupled;
659
660 // Update preconditioner
661 bool update_preconditioner_coupled;
662
663 // Update preconditioner every ... Newton iterations (only relevant for
664 // nonlinear problems, i.e., if the convective term is formulated implicitly)
665 // This variable is only used if update_preconditioner_coupled = true.
666 unsigned int update_preconditioner_coupled_every_newton_iter;
667
668 // Update preconditioner every ... time steps.
669 // This variable is only used if update_preconditioner_coupled = true.
670 unsigned int update_preconditioner_coupled_every_time_steps;
671
672 // description: see enum declaration
673 MomentumPreconditioner preconditioner_velocity_block;
674
675 // description: see enum declaration
676 MultigridOperatorType multigrid_operator_type_velocity_block;
677
678 // description: see declaration
679 MultigridData multigrid_data_velocity_block;
680
681 // The momentum block is inverted "exactly" in block preconditioner
682 // by solving the velocity convection-diffusion problem to a given
683 // relative tolerance
684 bool iterative_solve_of_velocity_block;
685
686 // solver data for velocity block (only relevant if velocity block
687 // is inverted exactly)
688 SolverData solver_data_velocity_block;
689
690 // description: see enum declaration
691 SchurComplementPreconditioner preconditioner_pressure_block;
692
693 // description: see declaration
694 MultigridData multigrid_data_pressure_block;
695
696 // The Laplace operator is inverted "exactly" in block preconditioner
697 // by solving the Laplace problem to a given relative tolerance
698 bool iterative_solve_of_pressure_block;
699
700 // solver data for Schur complement
701 // (only relevant if iterative_solve_of_pressure_block == true)
702 SolverData solver_data_pressure_block;
703};
704
705} // namespace IncNS
706} // namespace ExaDG
707
708#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_USER_INTERFACE_INPUT_PARAMETERS_H_ */
Definition driver.cpp:33
Definition grid_data.h:88
Definition viscosity_model_data.h:92
Definition viscosity_model_data.h:51
Definition inverse_mass_parameters.h:75
Definition inverse_mass_parameters.h:49
Definition multigrid_parameters.h:259
Definition newton_solver_data.h:36
Definition restart_data.h:38
Definition solver_data.h:34
Definition solver_info_data.h:38