version 3.9
seqsolverbackend.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_SEQ_SOLVER_BACKEND_HH
13#define DUMUX_SEQ_SOLVER_BACKEND_HH
14
15#include <type_traits>
16#include <tuple>
17#include <utility>
18
19#include <dune/istl/preconditioners.hh>
20#include <dune/istl/solvers.hh>
21#include <dune/istl/io.hh>
22#include <dune/common/indices.hh>
23#include <dune/common/hybridutilities.hh>
24
32
33namespace Dumux {
34
52{
53public:
54
55 // solve with generic parameter tree
56 template<class Preconditioner, class Solver, class Matrix, class Vector>
57 static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b,
58 const Dune::ParameterTree& params)
59 {
60 // make a linear operator from a matrix
61 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
62 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
63
64 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner"));
65 Solver solver(linearOperator, precond, params);
66
67 Vector bTmp(b);
68
69 Dune::InverseOperatorResult result;
70 solver.apply(x, bTmp, result);
71
72 return result.converged;
73 }
74};
75
82template<class M>
83constexpr std::size_t preconditionerBlockLevel() noexcept
84{
86}
87
96{
97public:
99
100 template<class Matrix, class Vector>
101 bool solve(const Matrix& A, Vector& x, const Vector& b)
102 {
103 Vector rhs(b);
104 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
105 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
106 jac.pre(x, rhs);
107 jac.apply(x, rhs);
108 jac.post(x);
109 return true;
110 }
111
112 std::string name() const
113 {
114 return "Explicit diagonal matrix solver";
115 }
116};
117
121// \{
122
127template <class LinearSolverTraits>
129{
130public:
132
133 template<class Matrix, class Vector>
134 bool solve(const Matrix& A, Vector& x, const Vector& b)
135 {
136 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
137 using Solver = Dune::BiCGSTABSolver<Vector>;
138 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
139 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
140 }
141
142 std::string name() const
143 {
144 return "Uzawa preconditioned BiCGSTAB solver";
145 }
146};
147
152template<class M, class X, class Y, int blockLevel = 2>
153class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
154{
155 template<std::size_t i>
156 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
157
158 template<std::size_t i>
159 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
160
161 template<std::size_t i>
162 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
163
164 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
165
166public:
168 using matrix_type = typename std::decay_t<M>;
170 using domain_type = X;
172 using range_type = Y;
174 using field_type = typename X::field_type;
175
182 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
183 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
184 {
185 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
186 }
187
188 void pre (X& v, Y& d) final {}
189
190 void apply (X& v, const Y& d) final
191 {
192 using namespace Dune::Hybrid;
193 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
194 {
195 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
196 });
197 }
198
199 void post (X&) final {}
200
202 Dune::SolverCategory::Category category() const final
203 {
204 return Dune::SolverCategory::sequential;
205 }
206
207private:
208 template<std::size_t... Is>
209 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
210 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
211 {}
212
213 ILUTuple ilu_;
214};
215
216
225{
226
227public:
229
230 template<class Matrix, class Vector>
231 bool solve(const Matrix& M, Vector& x, const Vector& b)
232 {
235 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
236 this->maxIter(), this->verbosity());
237 auto bTmp(b);
238 solver.apply(x, bTmp, result_);
239
240 return result_.converged;
241 }
242
243 const Dune::InverseOperatorResult& result() const
244 {
245 return result_;
246 }
247
248 std::string name() const
249 { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
250
251private:
252 Dune::InverseOperatorResult result_;
253};
254
263{
264
265public:
267
268 template<int precondBlockLevel = 2, class Matrix, class Vector>
269 bool solve(const Matrix& M, Vector& x, const Vector& b)
270 {
273 static const int restartGMRes = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.GMResRestart");
274 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes,
275 this->maxIter(), this->verbosity());
276 auto bTmp(b);
277 solver.apply(x, bTmp, result_);
278
279 return result_.converged;
280 }
281
282 const Dune::InverseOperatorResult& result() const
283 {
284 return result_;
285 }
286
287 std::string name() const
288 { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
289
290private:
291 Dune::InverseOperatorResult result_;
292};
293
298template<class M, class X, class Y, int blockLevel = 2>
299class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
300{
301 template<std::size_t i>
302 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
303
304 template<std::size_t i>
305 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
306
307 template<std::size_t i>
308 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
309
310 template<std::size_t i>
311 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
312
313 template<std::size_t i>
314 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
315
316 template<std::size_t i>
317 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
318
319 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
320
321public:
323 using matrix_type = typename std::decay_t<M>;
325 using domain_type = X;
327 using range_type = Y;
329 using field_type = typename X::field_type;
330
338 template<class LOP, class Criterion, class SmootherArgs>
339 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
340 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
341 {
342 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
343 }
344
361 void pre (X& v, Y& d) final
362 {
363 using namespace Dune::Hybrid;
364 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
365 {
366 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
367 });
368 }
369
381 void apply (X& v, const Y& d) final
382 {
383 using namespace Dune::Hybrid;
384 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
385 {
386 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
387 });
388 }
389
399 void post (X& v) final
400 {
401 using namespace Dune::Hybrid;
402 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
403 {
404 std::get<decltype(i)::value>(amg_).post(v[i]);
405 });
406 }
407
409 Dune::SolverCategory::Category category() const final
410 {
411 return Dune::SolverCategory::sequential;
412 }
413
414private:
415 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
416 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
417 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
418 {}
419
420 AMGTuple amg_;
421};
422
431{
432 template<class M, std::size_t i>
433 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
434
435 template<class X, std::size_t i>
436 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
437
438 template<class M, class X, std::size_t i>
439 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
440
441 template<class M, class X, std::size_t i>
442 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
443
444 template<class M, std::size_t i>
445 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
446
447 template<class M, class X, std::size_t i>
448 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
449
450public:
452
453 // Solve saddle-point problem using a Schur complement based preconditioner
454 template<class Matrix, class Vector>
455 bool solve(const Matrix& m, Vector& x, const Vector& b)
456 {
459 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
460 params.setDefaultValuesIsotropic(3);
461 params.setDebugLevel(this->verbosity());
462
463 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
464 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
465
466 using namespace Dune::Hybrid;
467 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
468 {
469 auto& args = std::get<decltype(i)::value>(smootherArgs);
470 args->iterations = 1;
471 args->relaxationFactor = 1;
472 });
473
474 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
475
476 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
477
479 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
480 this->maxIter(), this->verbosity());
481 auto bTmp(b);
482 solver.apply(x, bTmp, result_);
483
484 return result_.converged;
485 }
486
487 const Dune::InverseOperatorResult& result() const
488 {
489 return result_;
490 }
491
492 std::string name() const
493 { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
494
495private:
496 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
497 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
498 {
499 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
500 }
501
502 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
503 auto makeSmootherArgs_(std::index_sequence<Is...>)
504 {
505 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
506 }
507
508 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
509 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
510 {
511 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
512 }
513
514 Dune::InverseOperatorResult result_;
515};
516
517// \}
518
519} // end namespace Dumux
520
521#endif
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:431
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:487
std::string name() const
Definition: seqsolverbackend.hh:492
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:455
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:300
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:399
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:323
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:381
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:409
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:327
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:361
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:329
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:339
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:325
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:225
std::string name() const
Definition: seqsolverbackend.hh:248
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:243
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:231
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:154
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:174
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:168
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:172
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:188
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:170
void post(X &) final
Definition: seqsolverbackend.hh:199
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:202
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:190
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:182
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:263
std::string name() const
Definition: seqsolverbackend.hh:287
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:282
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:269
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:96
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:101
std::string name() const
Definition: seqsolverbackend.hh:112
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:52
static bool solveWithParamTree(const Matrix &A, Vector &x, const Vector &b, const Dune::ParameterTree &params)
Definition: seqsolverbackend.hh:57
Base class for linear solvers.
Definition: solver.hh:27
Scalar residReduction() const
the linear solver residual reduction
Definition: solver.hh:98
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:90
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:78
LinearSolver(const std::string &paramGroup="")
Construct the solver.
Definition: solver.hh:43
int verbosity() const
the verbosity level
Definition: solver.hh:82
static Dune::ParameterTree createParameterTree(const std::string &paramGroup="")
Create a tree containing parameters required for the linear solvers and precondioners of the Dune IST...
Definition: linearsolverparameters.hh:47
Adapter to turn a multi-type matrix into a thread-parallel linear operator. Adapts a matrix to the as...
Definition: parallelmatrixadapter.hh:28
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:70
A Uzawa preconditioned BiCGSTAB solver for saddle-point problems.
Definition: seqsolverbackend.hh:129
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:134
std::string name() const
Definition: seqsolverbackend.hh:142
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:83
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Type traits to be used with matrix types.
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
A parallel version of a linear operator.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:37
Definition: utility.hh:28
Utilities for template meta programming.