ExaDG
Loading...
Searching...
No Matches
elementwise_preconditioners.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_SOLVERS_AND_PRECONDITIONERS_PRECONDITIONER_ELEMENTWISE_PRECONDITIONERS_H_
23#define INCLUDE_EXADG_SOLVERS_AND_PRECONDITIONERS_PRECONDITIONER_ELEMENTWISE_PRECONDITIONERS_H_
24
25// deal.II
26#include <deal.II/matrix_free/operators.h>
27
28// ExaDG
29#include <exadg/matrix_free/integrators.h>
30#include <exadg/solvers_and_preconditioners/solvers/elementwise_krylov_solvers.h>
31
32namespace ExaDG
33{
34namespace Elementwise
35{
36/*
37 * Preconditioners
38 */
39template<typename Number>
41{
42public:
43 PreconditionerBase() : update_needed(true)
44 {
45 }
46
47 virtual ~PreconditionerBase()
48 {
49 }
50
51 virtual void
52 setup(unsigned int const cell) = 0;
53
54 virtual void
55 update() = 0;
56
57 bool
58 needs_update()
59 {
60 return update_needed;
61 }
62
63 virtual void
64 vmult(Number * dst, Number const * src) const = 0;
65
66protected:
67 bool update_needed;
68};
69
73template<typename Number>
75{
76public:
77 PreconditionerIdentity(unsigned int const size) : M(size)
78 {
79 this->update_needed = false;
80 }
81
83 {
84 }
85
86 void
87 setup(unsigned int const /* cell */) final
88 {
89 // nothing to do
90 }
91
92 void
93 update() final
94 {
95 // nothing to do
96 }
97
98 virtual void
99 vmult(Number * dst, Number const * src) const final
100 {
101 Number one;
102 one = 1.0;
103 equ(dst, one, src, M);
104 }
105
106private:
107 unsigned int const M;
108};
109
113template<int dim, int n_components, typename Number, typename Operator>
114class JacobiPreconditioner : public Elementwise::PreconditionerBase<dealii::VectorizedArray<Number>>
115{
116 typedef CellIntegrator<dim, n_components, Number> Integrator;
117
118public:
119 JacobiPreconditioner(dealii::MatrixFree<dim, Number> const & matrix_free,
120 unsigned int const dof_index,
121 unsigned int const quad_index,
122 Operator const & underlying_operator_in,
123 bool const initialize)
124 : underlying_operator(underlying_operator_in)
125 {
126 integrator = std::make_shared<Integrator>(matrix_free, dof_index, quad_index);
127
128 underlying_operator.initialize_dof_vector(global_inverse_diagonal);
129
130 if(initialize)
131 {
132 this->update();
133 }
134 }
135
136 void
137 setup(unsigned int cell) final
138 {
139 integrator->reinit(cell);
140 integrator->read_dof_values(global_inverse_diagonal, 0);
141 }
142
143 void
144 update() final
145 {
146 underlying_operator.calculate_inverse_diagonal(global_inverse_diagonal);
147
148 this->update_needed = false;
149 }
150
154 void
155 vmult(dealii::VectorizedArray<Number> * dst,
156 dealii::VectorizedArray<Number> const * src) const final
157 {
158 for(unsigned int i = 0; i < integrator->dofs_per_cell; ++i)
159 {
160 dst[i] = integrator->begin_dof_values()[i] * src[i];
161 }
162 }
163
164private:
165 std::shared_ptr<Integrator> integrator;
166
167 Operator const & underlying_operator;
168
169 dealii::LinearAlgebra::distributed::Vector<Number> global_inverse_diagonal;
170};
171
177template<int dim, int n_components, typename Number>
179 : public Elementwise::PreconditionerBase<dealii::VectorizedArray<Number>>
180{
181public:
182 typedef CellIntegrator<dim, n_components, Number> Integrator;
183
184 // use a template parameter of -1 to select the precompiled version of this operator
185 typedef dealii::MatrixFreeOperators::CellwiseInverseMassMatrix<dim, -1, n_components, Number>
186 CellwiseInverseMass;
187
188 InverseMassPreconditioner(dealii::MatrixFree<dim, Number> const & matrix_free,
189 unsigned int const dof_index,
190 unsigned int const quad_index)
191 {
192 integrator = std::make_shared<Integrator>(matrix_free, dof_index, quad_index);
193 inverse = std::make_shared<CellwiseInverseMass>(*integrator);
194
195 dealii::FiniteElement<dim> const & fe = matrix_free.get_dof_handler(dof_index).get_fe();
196
197 // The inverse mass preconditioner is only available for discontinuous Galerkin discretizations.
198 AssertThrow(
199 fe.conforms(dealii::FiniteElementData<dim>::L2),
200 dealii::ExcMessage(
201 "The elementwise inverse mass preconditioner is only implemented for DG (L2-conforming) elements."));
202
203 // Currently, the inverse mass realized as matrix-free operator evaluation is only available
204 // in deal.II for tensor-product elements.
205 AssertThrow(
206 fe.base_element(0).dofs_per_cell == dealii::Utilities::pow(fe.degree + 1, dim),
207 dealii::ExcMessage(
208 "The elementwise inverse mass preconditioner is only implemented for tensor-product DG elements."));
209
210 // Currently, the inverse mass realized as matrix-free operator evaluation is only available
211 // in deal.II if n_q_points_1d = n_nodes_1d.
212 AssertThrow(
213 matrix_free.get_shape_info(0, quad_index).data[0].n_q_points_1d == fe.degree + 1,
214 dealii::ExcMessage(
215 "The elementwise inverse mass preconditioner is only available if n_q_points_1d = n_nodes_1d."));
216
217 this->update_needed = false;
218 }
219
220 void
221 setup(unsigned int const cell) final
222 {
223 integrator->reinit(cell);
224 }
225
226 void
227 update() final
228 {
229 // no updates needed as long as the MatrixFree/Integrator object is up-to-date (which is not the
230 // responsibility of the present class).
231 }
232
236 void
237 vmult(dealii::VectorizedArray<Number> * dst,
238 dealii::VectorizedArray<Number> const * src) const final
239 {
240 inverse->apply(src, dst);
241 }
242
243private:
244 std::shared_ptr<Integrator> integrator;
245
246 std::shared_ptr<CellwiseInverseMass> inverse;
247};
248
249} // namespace Elementwise
250} // namespace ExaDG
251
252
253#endif /* INCLUDE_EXADG_SOLVERS_AND_PRECONDITIONERS_PRECONDITIONER_ELEMENTWISE_PRECONDITIONERS_H_ \
254 */
Definition elementwise_preconditioners.h:180
void vmult(dealii::VectorizedArray< Number > *dst, dealii::VectorizedArray< Number > const *src) const final
Definition elementwise_preconditioners.h:237
Definition elementwise_preconditioners.h:115
void vmult(dealii::VectorizedArray< Number > *dst, dealii::VectorizedArray< Number > const *src) const final
Definition elementwise_preconditioners.h:155
Definition elementwise_preconditioners.h:41
Definition elementwise_preconditioners.h:75
Definition driver.cpp:33