ExaDG
Loading...
Searching...
No Matches
jacobi_smoother.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_SOLVERS_AND_PRECONDITIONERS_JACOBISMOOTHER_H_
23#define INCLUDE_SOLVERS_AND_PRECONDITIONERS_JACOBISMOOTHER_H_
24
25#include <exadg/solvers_and_preconditioners/multigrid/multigrid_parameters.h>
26#include <exadg/solvers_and_preconditioners/multigrid/smoothers/smoother_base.h>
27#include <exadg/solvers_and_preconditioners/preconditioners/additive_schwarz_preconditioner.h>
28#include <exadg/solvers_and_preconditioners/preconditioners/block_jacobi_preconditioner.h>
29#include <exadg/solvers_and_preconditioners/preconditioners/jacobi_preconditioner.h>
30
31namespace ExaDG
32{
33template<typename Operator, typename VectorType>
34class JacobiSmoother : public SmootherBase<VectorType>
35{
36public:
37 JacobiSmoother() : underlying_operator(nullptr), preconditioner(nullptr)
38 {
39 }
40
42 {
43 delete preconditioner;
44 preconditioner = nullptr;
45 }
46
47 JacobiSmoother(JacobiSmoother const &) = delete;
48
50 operator=(JacobiSmoother const &) = delete;
51
53 {
58 : preconditioner(PreconditionerSmoother::PointJacobi),
59 number_of_smoothing_steps(5),
60 damping_factor(1.0)
61 {
62 }
63
64 // preconditioner
65 PreconditionerSmoother preconditioner;
66
67 // number of iterations per smoothing step
68 unsigned int number_of_smoothing_steps;
69
70 // damping factor
71 double damping_factor;
72 };
73
74 void
75 setup(Operator const & operator_in,
76 bool const initialize_preconditioner,
77 AdditionalData const & additional_data_in)
78 {
79 underlying_operator = &operator_in;
80
81 data = additional_data_in;
82
83 if(data.preconditioner == PreconditionerSmoother::PointJacobi)
84 {
85 preconditioner =
86 new JacobiPreconditioner<Operator>(*underlying_operator, initialize_preconditioner);
87 }
88 else if(data.preconditioner == PreconditionerSmoother::BlockJacobi)
89 {
90 preconditioner =
91 new BlockJacobiPreconditioner<Operator>(*underlying_operator, initialize_preconditioner);
92 }
93 else if(data.preconditioner == PreconditionerSmoother::AdditiveSchwarz)
94 {
95 preconditioner = new AdditiveSchwarzPreconditioner<Operator>(*underlying_operator,
96 initialize_preconditioner);
97 }
98 else
99 {
100 AssertThrow(data.preconditioner == PreconditionerSmoother::PointJacobi or
101 data.preconditioner == PreconditionerSmoother::BlockJacobi,
102 dealii::ExcMessage(
103 "Specified type of preconditioner for Jacobi smoother not implemented."));
104 }
105 }
106
107 void
108 update() final
109 {
110 if(preconditioner != nullptr)
111 preconditioner->update();
112 }
113
114 /*
115 * Approximately solve linear system of equations (b=src, x=dst)
116 *
117 * A*x = b (r=b-A*x)
118 *
119 * using the iteration
120 *
121 * x^{k+1} = x^{k} + omega * P^{-1} * r^{k}
122 *
123 * where
124 *
125 * omega: damping factor
126 * P: preconditioner
127 */
128 void
129 vmult(VectorType & dst, VectorType const & src) const final
130 {
131 dst = 0;
132
133 VectorType tmp(src), residual(src);
134
135 for(unsigned int k = 0; k < data.number_of_smoothing_steps; ++k)
136 {
137 if(k > 0)
138 {
139 // calculate residual r^{k} = src - A * x^{k}
140 underlying_operator->vmult(residual, dst);
141 residual.sadd(-1.0, 1.0, src);
142 }
143 else // we do not have to evaluate the residual for k=0 since dst = 0
144 {
145 residual = src;
146 }
147
148 // apply preconditioner: tmp = P^{-1} * residual
149 preconditioner->vmult(tmp, residual);
150
151 // x^{k+1} = x^{k} + damping_factor * tmp
152 dst.add(data.damping_factor, tmp);
153 }
154 }
155
156 void
157 step(VectorType & dst, VectorType const & src) const final
158 {
159 VectorType tmp(src), residual(src);
160
161 for(unsigned int k = 0; k < data.number_of_smoothing_steps; ++k)
162 {
163 // calculate residual r^{k} = src - A * x^{k}
164 underlying_operator->vmult(residual, dst);
165 residual.sadd(-1.0, 1.0, src);
166
167 // apply preconditioner: tmp = P^{-1} * residual
168 preconditioner->vmult(tmp, residual);
169
170 // x^{k+1} = x^{k} + damping_factor * tmp
171 dst.add(data.damping_factor, tmp);
172 }
173 }
174
175private:
176 Operator const * underlying_operator;
177
178 AdditionalData data;
179
180 PreconditionerBase<typename Operator::value_type> * preconditioner;
181};
182} // namespace ExaDG
183
184
185#endif /* INCLUDE_SOLVERS_AND_PRECONDITIONERS_JACOBISMOOTHER_H_ */
Definition jacobi_preconditioner.h:35
Definition jacobi_smoother.h:35
Definition smoother_base.h:29
Definition driver.cpp:33
Definition jacobi_smoother.h:53
AdditionalData()
Definition jacobi_smoother.h:57