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