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