version 3.8
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 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
56 [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]]
57 static bool solve(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
58 const std::string& modelParamGroup = "")
59 {
60 Preconditioner precond(A, s.precondIter(), s.relaxation());
61
62 // make a linear operator from a matrix
63 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
64 MatrixAdapter linearOperator(A);
65
66 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
67
68 Vector bTmp(b);
69
70 Dune::InverseOperatorResult result;
71 solver.apply(x, bTmp, result);
72
73 return result.converged;
74 }
75
76 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
77 [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]]
78 static bool solveWithGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
79 const std::string& modelParamGroup = "")
80 {
81 // get the restart threshold
82 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
83
84 Preconditioner precond(A, s.precondIter(), s.relaxation());
85
86 // make a linear operator from a matrix
87 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
88 MatrixAdapter linearOperator(A);
89
90 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
91
92 Vector bTmp(b);
93
94 Dune::InverseOperatorResult result;
95 solver.apply(x, bTmp, result);
96
97 return result.converged;
98 }
99
100 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
101 [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]]
102 static bool solveWithILU0Prec(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
103 const std::string& modelParamGroup = "")
104 {
105 Preconditioner precond(A, s.relaxation());
106
107 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
108 MatrixAdapter operatorA(A);
109
110 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
111
112 Vector bTmp(b);
113
114 Dune::InverseOperatorResult result;
115 solver.apply(x, bTmp, result);
116
117 return result.converged;
118 }
119
120 // solve with RestartedGMRes (needs restartGMRes as additional argument)
121 template<class Preconditioner, class Solver, class SolverInterface, class Matrix, class Vector>
122 [[deprecated("Removed after 3.8. Use solver from istlsolvers.hh")]]
123 static bool solveWithILU0PrecGMRes(const SolverInterface& s, const Matrix& A, Vector& x, const Vector& b,
124 const std::string& modelParamGroup = "")
125 {
126 // get the restart threshold
127 const int restartGMRes = getParamFromGroup<int>(modelParamGroup, "LinearSolver.GMResRestart", 10);
128
129 Preconditioner precond(A, s.relaxation());
130
131 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
132 MatrixAdapter operatorA(A);
133
134 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
135
136 Vector bTmp(b);
137
138 Dune::InverseOperatorResult result;
139 solver.apply(x, bTmp, result);
140
141 return result.converged;
142 }
143
144 // solve with generic parameter tree
145 template<class Preconditioner, class Solver, class Matrix, class Vector>
146 static bool solveWithParamTree(const Matrix& A, Vector& x, const Vector& b,
147 const Dune::ParameterTree& params)
148 {
149 // make a linear operator from a matrix
150 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
151 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
152
153 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub("preconditioner"));
154 Solver solver(linearOperator, precond, params);
155
156 Vector bTmp(b);
157
158 Dune::InverseOperatorResult result;
159 solver.apply(x, bTmp, result);
160
161 return result.converged;
162 }
163};
164
171template<class M>
172constexpr std::size_t preconditionerBlockLevel() noexcept
173{
175}
176
185{
186public:
188
189 template<class Matrix, class Vector>
190 bool solve(const Matrix& A, Vector& x, const Vector& b)
191 {
192 Vector rhs(b);
193 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
194 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
195 jac.pre(x, rhs);
196 jac.apply(x, rhs);
197 jac.post(x);
198 return true;
199 }
200
201 std::string name() const
202 {
203 return "Explicit diagonal matrix solver";
204 }
205};
206
210// \{
211
216template <class LinearSolverTraits>
218{
219public:
221
222 template<class Matrix, class Vector>
223 bool solve(const Matrix& A, Vector& x, const Vector& b)
224 {
225 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
226 using Solver = Dune::BiCGSTABSolver<Vector>;
227 static const auto solverParams = LinearSolverParameters<LinearSolverTraits>::createParameterTree(this->paramGroup());
228 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
229 }
230
231 std::string name() const
232 {
233 return "Uzawa preconditioned BiCGSTAB solver";
234 }
235};
236
241template<class M, class X, class Y, int blockLevel = 2>
242class BlockDiagILU0Preconditioner : public Dune::Preconditioner<X, Y>
243{
244 template<std::size_t i>
245 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
246
247 template<std::size_t i>
248 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
249
250 template<std::size_t i>
251 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
252
253 using ILUTuple = typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
254
255public:
257 using matrix_type = typename std::decay_t<M>;
259 using domain_type = X;
261 using range_type = Y;
263 using field_type = typename X::field_type;
264
271 BlockDiagILU0Preconditioner(const M& m, double w = 1.0)
272 : BlockDiagILU0Preconditioner(m, w, std::make_index_sequence<M::N()>{})
273 {
274 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
275 }
276
277 void pre (X& v, Y& d) final {}
278
279 void apply (X& v, const Y& d) final
280 {
281 using namespace Dune::Hybrid;
282 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](const auto i)
283 {
284 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
285 });
286 }
287
288 void post (X&) final {}
289
291 Dune::SolverCategory::Category category() const final
292 {
293 return Dune::SolverCategory::sequential;
294 }
295
296private:
297 template<std::size_t... Is>
298 BlockDiagILU0Preconditioner (const M& m, double w, std::index_sequence<Is...> is)
299 : ilu_(std::make_tuple(BlockILU<Is>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
300 {}
301
302 ILUTuple ilu_;
303};
304
305
314{
315
316public:
318
319 template<class Matrix, class Vector>
320 bool solve(const Matrix& M, Vector& x, const Vector& b)
321 {
324 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
325 this->maxIter(), this->verbosity());
326 auto bTmp(b);
327 solver.apply(x, bTmp, result_);
328
329 return result_.converged;
330 }
331
332 const Dune::InverseOperatorResult& result() const
333 {
334 return result_;
335 }
336
337 std::string name() const
338 { return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
339
340private:
341 Dune::InverseOperatorResult result_;
342};
343
352{
353
354public:
356
357 template<int precondBlockLevel = 2, class Matrix, class Vector>
358 bool solve(const Matrix& M, Vector& x, const Vector& b)
359 {
362 static const int restartGMRes = getParamFromGroup<int>(this->paramGroup(), "LinearSolver.GMResRestart");
363 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->residReduction(), restartGMRes,
364 this->maxIter(), this->verbosity());
365 auto bTmp(b);
366 solver.apply(x, bTmp, result_);
367
368 return result_.converged;
369 }
370
371 const Dune::InverseOperatorResult& result() const
372 {
373 return result_;
374 }
375
376 std::string name() const
377 { return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
378
379private:
380 Dune::InverseOperatorResult result_;
381};
382
387template<class M, class X, class Y, int blockLevel = 2>
388class BlockDiagAMGPreconditioner : public Dune::Preconditioner<X, Y>
389{
390 template<std::size_t i>
391 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
392
393 template<std::size_t i>
394 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
395
396 template<std::size_t i>
397 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
398
399 template<std::size_t i>
400 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
401
402 template<std::size_t i>
403 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
404
405 template<std::size_t i>
406 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
407
408 using AMGTuple = typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
409
410public:
412 using matrix_type = typename std::decay_t<M>;
414 using domain_type = X;
416 using range_type = Y;
418 using field_type = typename X::field_type;
419
427 template<class LOP, class Criterion, class SmootherArgs>
428 BlockDiagAMGPreconditioner(const LOP& lop, const Criterion& c, const SmootherArgs& sa)
429 : BlockDiagAMGPreconditioner(lop, c, sa, std::make_index_sequence<M::N()>{})
430 {
431 static_assert(blockLevel >= 2, "Only makes sense for MultiTypeBlockMatrix!");
432 }
433
450 void pre (X& v, Y& d) final
451 {
452 using namespace Dune::Hybrid;
453 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
454 {
455 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
456 });
457 }
458
470 void apply (X& v, const Y& d) final
471 {
472 using namespace Dune::Hybrid;
473 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
474 {
475 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
476 });
477 }
478
488 void post (X& v) final
489 {
490 using namespace Dune::Hybrid;
491 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](const auto i)
492 {
493 std::get<decltype(i)::value>(amg_).post(v[i]);
494 });
495 }
496
498 Dune::SolverCategory::Category category() const final
499 {
500 return Dune::SolverCategory::sequential;
501 }
502
503private:
504 template<class LOP, class Criterion, class SmootherArgs, std::size_t... Is>
505 BlockDiagAMGPreconditioner (const LOP& lop, const Criterion& c, const SmootherArgs& sa, std::index_sequence<Is...> is)
506 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
507 {}
508
509 AMGTuple amg_;
510};
511
520{
521 template<class M, std::size_t i>
522 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
523
524 template<class X, std::size_t i>
525 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
526
527 template<class M, class X, std::size_t i>
528 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
529
530 template<class M, class X, std::size_t i>
531 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
532
533 template<class M, std::size_t i>
534 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
535
536 template<class M, class X, std::size_t i>
537 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
538
539public:
541
542 // Solve saddle-point problem using a Schur complement based preconditioner
543 template<class Matrix, class Vector>
544 bool solve(const Matrix& m, Vector& x, const Vector& b)
545 {
548 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
549 params.setDefaultValuesIsotropic(3);
550 params.setDebugLevel(this->verbosity());
551
552 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
553 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
554
555 using namespace Dune::Hybrid;
556 forEach(std::make_index_sequence<Matrix::N()>{}, [&](const auto i)
557 {
558 auto& args = std::get<decltype(i)::value>(smootherArgs);
559 args->iterations = 1;
560 args->relaxationFactor = 1;
561 });
562
563 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
564
565 BlockDiagAMGPreconditioner<Matrix, Vector, Vector> preconditioner(linearOperator, criterion, smootherArgs);
566
568 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->residReduction(),
569 this->maxIter(), this->verbosity());
570 auto bTmp(b);
571 solver.apply(x, bTmp, result_);
572
573 return result_.converged;
574 }
575
576 const Dune::InverseOperatorResult& result() const
577 {
578 return result_;
579 }
580
581 std::string name() const
582 { return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
583
584private:
585 template<template<class M, std::size_t i> class Criterion, class Matrix, class Params, std::size_t... Is>
586 auto makeCriterion_(const Params& p, std::index_sequence<Is...>)
587 {
588 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
589 }
590
591 template<template<class M, class X, std::size_t i> class SmootherArgs, class Matrix, class Vector, std::size_t... Is>
592 auto makeSmootherArgs_(std::index_sequence<Is...>)
593 {
594 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
595 }
596
597 template<template<class M, class X, std::size_t i> class LinearOperator, class Matrix, class Vector, std::size_t... Is>
598 auto makeLinearOperator_(const Matrix& m, std::index_sequence<Is...>)
599 {
600 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
601 }
602
603 Dune::InverseOperatorResult result_;
604};
605
606// \}
607
608} // end namespace Dumux
609
610#endif
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:520
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:576
std::string name() const
Definition: seqsolverbackend.hh:581
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:544
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:389
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:488
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:412
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:470
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:498
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:416
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:450
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:418
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:428
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:414
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:314
std::string name() const
Definition: seqsolverbackend.hh:337
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:332
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:320
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:243
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:263
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:257
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:261
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:277
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:259
void post(X &) final
Definition: seqsolverbackend.hh:288
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:291
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:279
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:271
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:352
std::string name() const
Definition: seqsolverbackend.hh:376
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:371
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:358
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:185
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:190
std::string name() const
Definition: seqsolverbackend.hh:201
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:52
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:123
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:78
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:102
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:57
static bool solveWithParamTree(const Matrix &A, Vector &x, const Vector &b, const Dune::ParameterTree &params)
Definition: seqsolverbackend.hh:146
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:218
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:223
std::string name() const
Definition: seqsolverbackend.hh:231
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:172
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.