ExaDG
Loading...
Searching...
No Matches
coarse_grid_solvers.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_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_COARSE_GRID_SOLVERS_H_
23#define EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_COARSE_GRID_SOLVERS_H_
24
25// deal.II
26#include <deal.II/lac/la_parallel_vector.h>
27#include <deal.II/lac/solver_cg.h>
28#include <deal.II/lac/solver_control.h>
29#include <deal.II/lac/solver_gmres.h>
30#include <deal.II/multigrid/mg_base.h>
31#include <deal.II/numerics/vector_tools_mean_value.h>
32
33// ExaDG
34#include <exadg/solvers_and_preconditioners/multigrid/smoothers/chebyshev_smoother.h>
35#include <exadg/solvers_and_preconditioners/preconditioners/block_jacobi_preconditioner.h>
36#include <exadg/solvers_and_preconditioners/preconditioners/jacobi_preconditioner.h>
37#include <exadg/solvers_and_preconditioners/preconditioners/preconditioner_amg.h>
38#include <exadg/solvers_and_preconditioners/preconditioners/preconditioner_base.h>
39#include <exadg/solvers_and_preconditioners/solvers/iterative_solvers_dealii_wrapper.h>
40#include <exadg/solvers_and_preconditioners/solvers/solver_data.h>
41#include <exadg/solvers_and_preconditioners/utilities/linear_algebra_utilities.h>
42
43namespace ExaDG
44{
49template<typename Operator>
52 dealii::LinearAlgebra::distributed::Vector<typename Operator::value_type>>
53{
54public:
55 virtual ~CoarseGridSolverBase(){};
56
57 virtual void
58 update() = 0;
59};
60
61template<typename Operator>
62class MGCoarseKrylov : public CoarseGridSolverBase<Operator>
63{
64public:
65 typedef typename Operator::value_type Number;
66
67 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
68
70 {
75 : solver_type(MultigridCoarseGridSolver::CG),
76 solver_data(SolverData(1e4, 1.e-12, 1.e-3, LinearSolver::Undefined, 100)),
77 operator_is_singular(false),
78 preconditioner(MultigridCoarseGridPreconditioner::None),
79 amg_data(AMGData())
80 {
81 }
82
83 // Type of Krylov solver
84 MultigridCoarseGridSolver solver_type;
85
86 // Solver data
87 SolverData solver_data;
88
89 // in case of singular operators (with constant vectors forming the nullspace) the rhs vector
90 // has to be projected onto the space of vectors with zero mean prior to solving the coarse
91 // grid problem
92 bool operator_is_singular;
93
94 // Preconditioner
95 MultigridCoarseGridPreconditioner preconditioner;
96
97 // Configuration of AMG settings
98 AMGData amg_data;
99 };
100
101 MGCoarseKrylov(Operator const & pde_operator_in,
102 bool const initialize,
103 AdditionalData const & additional_data,
104 MPI_Comm const & comm)
105 : pde_operator(pde_operator_in), additional_data(additional_data), mpi_comm(comm)
106 {
107 if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::PointJacobi)
108 {
109 preconditioner = std::make_shared<JacobiPreconditioner<Operator>>(pde_operator, initialize);
110
111 std::shared_ptr<JacobiPreconditioner<Operator>> jacobi =
112 std::dynamic_pointer_cast<JacobiPreconditioner<Operator>>(preconditioner);
113 AssertDimension(jacobi->get_size_of_diagonal(), pde_operator.m());
114 }
115 else if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::BlockJacobi)
116 {
117 preconditioner =
118 std::make_shared<BlockJacobiPreconditioner<Operator>>(pde_operator, initialize);
119 }
120 else if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::AMG)
121 {
122 preconditioner =
123 std::make_shared<PreconditionerAMG<Operator, Number>>(pde_operator,
124 initialize,
125 additional_data.amg_data);
126 }
127 else
128 {
129 AssertThrow(
130 additional_data.preconditioner == MultigridCoarseGridPreconditioner::None or
131 additional_data.preconditioner == MultigridCoarseGridPreconditioner::PointJacobi or
132 additional_data.preconditioner == MultigridCoarseGridPreconditioner::BlockJacobi or
133 additional_data.preconditioner == MultigridCoarseGridPreconditioner::AMG,
134 dealii::ExcMessage("Specified preconditioner for PCG coarse grid solver not implemented."));
135 }
136 }
137
138 void
139 update() final
140 {
141 if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::None)
142 {
143 // do nothing
144 }
145 else if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::PointJacobi or
146 additional_data.preconditioner == MultigridCoarseGridPreconditioner::BlockJacobi or
147 additional_data.preconditioner == MultigridCoarseGridPreconditioner::AMG)
148 {
149 preconditioner->update();
150 }
151 else
152 {
153 AssertThrow(false, dealii::ExcMessage("Not implemented."));
154 }
155 }
156
157 void
158 operator()(unsigned int const, VectorType & dst, VectorType const & src) const final
159 {
160 VectorType r(src);
161 if(additional_data.operator_is_singular)
162 dealii::VectorTools::subtract_mean_value(r);
163
164 if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::AMG)
165 {
166 std::shared_ptr<PreconditionerAMG<Operator, Number>> preconditioner_amg =
167 std::dynamic_pointer_cast<PreconditionerAMG<Operator, Number>>(preconditioner);
168
169 AssertThrow(preconditioner_amg != nullptr,
170 dealii::ExcMessage("Constructed preconditioner is not of type "
171 "PreconditionerAMG<Operator, Number>."));
172
173 preconditioner_amg->apply_krylov_solver_with_amg_preconditioner(dst,
174 r,
175 additional_data.solver_type,
176 additional_data.solver_data);
177 }
178 else
179 {
180 // Overwrite `linear_solver` as enum also enables non-Krylov coarse grid solvers.
181 SolverData solver_data = additional_data.solver_data;
182 if(additional_data.solver_type == MultigridCoarseGridSolver::CG)
183 {
184 solver_data.linear_solver = LinearSolver::CG;
185 }
186 else if(additional_data.solver_type == MultigridCoarseGridSolver::GMRES)
187 {
188 solver_data.linear_solver = LinearSolver::GMRES;
189 }
190 else
191 {
192 AssertThrow(false, dealii::ExcMessage("Not implemented."));
193 }
194
195 bool constexpr compute_performance_metrics = false;
196 bool constexpr compute_eigenvalues = false;
197 bool use_preconditioner;
198 if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::None)
199 {
200 use_preconditioner = false;
201 }
202 else if(additional_data.preconditioner == MultigridCoarseGridPreconditioner::PointJacobi or
203 additional_data.preconditioner == MultigridCoarseGridPreconditioner::BlockJacobi)
204 {
205 use_preconditioner = true;
206 }
207 else
208 {
209 AssertThrow(false, dealii::ExcMessage("Not implemented."));
210 }
211
212 typedef Krylov::KrylovSolver<Operator, PreconditionerBase<Number>, VectorType> SolverType;
213
214 SolverType solver(pde_operator,
215 *preconditioner,
216 solver_data,
217 use_preconditioner,
218 compute_performance_metrics,
219 compute_eigenvalues);
220
221 // Note that the preconditioner has already been updated
222 solver.solve(dst, r);
223 }
224 }
225
226private:
227 Operator const & pde_operator;
228
229 std::shared_ptr<PreconditionerBase<Number>> preconditioner;
230
231 AdditionalData additional_data;
232
233 MPI_Comm const mpi_comm;
234};
235
236
237template<typename Operator>
238class MGCoarseChebyshev : public CoarseGridSolverBase<Operator>
239{
240public:
241 typedef typename Operator::value_type MultigridNumber;
242
243 typedef dealii::LinearAlgebra::distributed::Vector<MultigridNumber> VectorType;
244
245 typedef dealii::PreconditionChebyshev<Operator, VectorType, dealii::DiagonalMatrix<VectorType>>
246 DealiiChebyshev;
247
248 MGCoarseChebyshev(Operator const & coarse_operator_in,
249 bool const initialize_preconditioner_in,
250 double const relative_tolerance_in,
251 MultigridCoarseGridPreconditioner const & preconditioner,
252 bool const operator_is_singular_in)
253 : coarse_operator(coarse_operator_in),
254 relative_tolerance(relative_tolerance_in),
255 operator_is_singular(operator_is_singular_in)
256 {
257 AssertThrow(preconditioner == MultigridCoarseGridPreconditioner::PointJacobi,
258 dealii::ExcMessage(
259 "Only PointJacobi preconditioner implemented for Chebyshev coarse-grid solver."));
260
261 if(initialize_preconditioner_in)
262 {
263 update();
264 }
265 }
266
267 void
268 update() final
269 {
270 // use Chebyshev smoother of high degree to solve the coarse grid problem approximately
271 typename DealiiChebyshev::AdditionalData dealii_additional_data;
272
273 std::shared_ptr<dealii::DiagonalMatrix<VectorType>> diagonal_matrix =
274 std::make_shared<dealii::DiagonalMatrix<VectorType>>();
275 VectorType & diagonal_vector = diagonal_matrix->get_vector();
276
277 coarse_operator.initialize_dof_vector(diagonal_vector);
278 coarse_operator.calculate_inverse_diagonal(diagonal_vector);
279
280 std::pair<double, double> eigenvalues =
281 estimate_eigenvalues_cg(coarse_operator, diagonal_vector, operator_is_singular);
282
283 double const factor = 1.1;
284
285 dealii_additional_data.preconditioner = diagonal_matrix;
286 dealii_additional_data.max_eigenvalue = factor * eigenvalues.second;
287 dealii_additional_data.smoothing_range = eigenvalues.second / eigenvalues.first * factor;
288
289 double sigma = (1. - std::sqrt(1. / dealii_additional_data.smoothing_range)) /
290 (1. + std::sqrt(1. / dealii_additional_data.smoothing_range));
291
292 // calculate/estimate the number of Chebyshev iterations needed to reach a specified relative
293 // solver tolerance
294 double const eps = relative_tolerance;
295
296 dealii_additional_data.degree = static_cast<unsigned int>(
297 std::log(1. / eps + std::sqrt(1. / eps / eps - 1.)) / std::log(1. / sigma));
298 dealii_additional_data.eig_cg_n_iterations = 0;
299
300 chebyshev_smoother = std::make_shared<DealiiChebyshev>();
301 chebyshev_smoother->initialize(coarse_operator, dealii_additional_data);
302 }
303
304 void
305 operator()(unsigned int const level, VectorType & dst, const VectorType & src) const final
306 {
307 AssertThrow(chebyshev_smoother.get() != 0,
308 dealii::ExcMessage("MGCoarseChebyshev: chebyshev_smoother is not initialized."));
309
310 AssertThrow(level == 0, dealii::ExcNotImplemented());
311
312 chebyshev_smoother->vmult(dst, src);
313 }
314
315private:
316 Operator const & coarse_operator;
317 double const relative_tolerance;
318 bool const operator_is_singular;
319
320 std::shared_ptr<DealiiChebyshev> chebyshev_smoother;
321};
322
327template<typename Operator>
328class MGCoarseAMG : public CoarseGridSolverBase<Operator>
329{
330private:
331 typedef typename Operator::value_type Number;
332
333 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
334
335public:
336 MGCoarseAMG(Operator const & op, bool const initialize, AMGData data = AMGData())
337 {
338 amg_preconditioner =
339 std::make_shared<PreconditionerAMG<Operator, Number>>(op, initialize, data);
340 }
341
342 void
343 update() final
344 {
345 amg_preconditioner->update();
346 }
347
348 void
349 operator()(unsigned int const /*level*/, VectorType & dst, VectorType const & src) const final
350 {
351 amg_preconditioner->vmult(dst, src);
352 }
353
354private:
355 std::shared_ptr<PreconditionerAMG<Operator, Number>> amg_preconditioner;
356};
357
358} // namespace ExaDG
359
360#endif /* EXADG_SOLVERS_AND_PRECONDITIONERS_MULTIGRID_COARSE_GRID_SOLVERS_H_ */
Definition coarse_grid_solvers.h:53
Definition coarse_grid_solvers.h:63
Definition multigrid_preconditioner_base.h:55
Definition driver.cpp:33
Definition multigrid_parameters.h:98
Definition coarse_grid_solvers.h:70
AdditionalData()
Definition coarse_grid_solvers.h:74
Definition solver_data.h:92