ExaDG
Loading...
Searching...
No Matches
inverse_mass_operator.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_OPERATORS_INVERSE_MASS_OPERATOR_H_
23#define EXADG_OPERATORS_INVERSE_MASS_OPERATOR_H_
24
25// deal.II
26#include <deal.II/fe/fe_data.h>
27#include <deal.II/lac/la_parallel_vector.h>
28#include <deal.II/matrix_free/operators.h>
29
30// ExaDG
31#include <exadg/grid/grid_data.h>
32#include <exadg/matrix_free/integrators.h>
33#include <exadg/operators/inverse_mass_parameters.h>
34#include <exadg/operators/mass_operator.h>
35#include <exadg/operators/variable_coefficients.h>
36#include <exadg/solvers_and_preconditioners/preconditioners/block_jacobi_preconditioner.h>
37#include <exadg/solvers_and_preconditioners/preconditioners/jacobi_preconditioner.h>
38#include <exadg/solvers_and_preconditioners/preconditioners/preconditioner_amg.h>
39
40namespace ExaDG
41{
42template<typename Number>
43struct InverseMassOperatorData
44{
45 InverseMassOperatorData()
46 : dof_index(0),
47 quad_index(0),
48 coefficient_is_variable(false),
49 consider_inverse_coefficient(false),
50 variable_coefficients(nullptr)
51 {
52 }
53
54 // Get optimal in the sense of (most likely) fastest implementation type of the inverse mass
55 // operator depending on the approximation space.
56 template<int dim>
57 static InverseMassType
58 get_optimal_inverse_mass_type(dealii::FiniteElement<dim> const & fe)
59 {
60 if(fe.conforms(dealii::FiniteElementData<dim>::L2))
61 {
62 if(fe.reference_cell().is_hyper_cube())
63 {
64 return InverseMassType::MatrixfreeOperator;
65 }
66 else
67 {
68 return InverseMassType::ElementwiseKrylovSolver;
69 }
70 }
71 else
72 {
73 return InverseMassType::GlobalKrylovSolver;
74 }
75 }
76
77 // Parameters referring to dealii::MatrixFree
78 unsigned int dof_index;
79 unsigned int quad_index;
80
81 InverseMassParameters parameters;
82
83 // Enable variable coefficients.
84 bool coefficient_is_variable;
85
86 // Consider the regular form of the coefficient (1) or its inverse (2):
87 // (1) : (u_h , v_h * c)_Omega
88 // (2) : (u_h , v_h / c)_Omega
89 bool consider_inverse_coefficient;
90
91 VariableCoefficients<dealii::VectorizedArray<Number>> const * variable_coefficients;
92};
93
94template<int dim, int n_components, typename Number>
95class InverseMassOperator
96{
97private:
98 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
99
100 typedef InverseMassOperator<dim, n_components, Number> This;
101
102 typedef CellIntegrator<dim, n_components, Number> Integrator;
103
104 // use a template parameter of -1 to select the precompiled version of this operator
105 typedef dealii::MatrixFreeOperators::CellwiseInverseMassMatrix<dim, -1, n_components, Number>
106 InverseMassAsMatrixFreeOperator;
107
108 typedef std::pair<unsigned int, unsigned int> Range;
109
110public:
111 InverseMassOperator();
112
113 unsigned int
114 get_n_iter_global_last() const;
115
116 void
117 initialize(dealii::MatrixFree<dim, Number> const & matrix_free_in,
118 InverseMassOperatorData<Number> const inverse_mass_operator_data,
119 dealii::AffineConstraints<Number> const * constraints = nullptr);
120
125 void
126 update();
127
128 // dst = M^-1 * src
129 void
130 apply(VectorType & dst, VectorType const & src) const;
131
132 // dst = scaling_factor * (M^-1 * src)
133 void
134 apply_scale(VectorType & dst, double const scaling_factor, VectorType const & src) const;
135
136
137private:
138 void
139 cell_loop_matrix_free_operator(dealii::MatrixFree<dim, Number> const &,
140 VectorType & dst,
141 VectorType const & src,
142 Range const & cell_range) const;
143
144 dealii::MatrixFree<dim, Number> const * matrix_free;
145
146 unsigned int dof_index, quad_index;
147
149
150 // Variable coefficients not managed by this class.
151 bool coefficient_is_variable;
152 bool consider_inverse_coefficient;
153
154 VariableCoefficients<dealii::VectorizedArray<Number>> const * variable_coefficients;
155
156 // Solver and preconditioner for solving a global linear system of equations for all degrees of
157 // freedom.
158 std::shared_ptr<PreconditionerBase<Number>> global_preconditioner;
159 std::shared_ptr<Krylov::SolverBase<VectorType>> global_solver;
160
161 // Iterations needed in global Krylov solver at last inverse application.
162 mutable unsigned int n_iter_global_last = dealii::numbers::invalid_unsigned_int;
163
164 // Block-Jacobi preconditioner used as cell-wise inverse for L2-conforming spaces. In this case,
165 // the mass matrix is block-diagonal and a block-Jacobi preconditioner inverts the mass operator
166 // exactly (up to solver tolerances). The implementation of the block-Jacobi preconditioner can be
167 // matrix-based or matrix-free, depending on the parameters specified.
168 std::shared_ptr<BlockJacobiPreconditioner<MassOperator<dim, n_components, Number>>>
169 block_jacobi_preconditioner;
170
171 // MassOperator as underlying operator for the cell-wise or global iterative solves.
173};
174
175} // namespace ExaDG
176
177#endif /* EXADG_OPERATORS_INVERSE_MASS_OPERATOR_H_ */
void update()
Definition inverse_mass_operator.cpp:191
Definition mass_operator.h:53
Definition variable_coefficients.h:40
Definition driver.cpp:33
Definition inverse_mass_operator.h:44
Definition inverse_mass_parameters.h:56