24#ifndef DUMUX_SEQ_SOLVER_BACKEND_HH
25#define DUMUX_SEQ_SOLVER_BACKEND_HH
31#include <dune/istl/preconditioners.hh>
32#include <dune/istl/solvers.hh>
33#include <dune/istl/superlu.hh>
34#include <dune/istl/umfpack.hh>
35#include <dune/istl/io.hh>
36#include <dune/common/indices.hh>
37#include <dune/common/hybridutilities.hh>
69 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
70 static bool solve(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
71 const std::string& modelParamGroup =
"")
73 Preconditioner precond(A, s.precondIter(), s.relaxation());
76 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
77 MatrixAdapter linearOperator(A);
79 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
83 Dune::InverseOperatorResult result;
84 solver.apply(x, bTmp, result);
86 return result.converged;
89 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
90 static bool solveWithGMRes(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
91 const std::string& modelParamGroup =
"")
94 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart", 10);
96 Preconditioner precond(A, s.precondIter(), s.relaxation());
99 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
100 MatrixAdapter linearOperator(A);
102 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
106 Dune::InverseOperatorResult result;
107 solver.apply(x, bTmp, result);
109 return result.converged;
112 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
113 static bool solveWithILU0Prec(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
114 const std::string& modelParamGroup =
"")
116 Preconditioner precond(A, s.relaxation());
118 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
119 MatrixAdapter operatorA(A);
121 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
125 Dune::InverseOperatorResult result;
126 solver.apply(x, bTmp, result);
128 return result.converged;
132 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
134 const std::string& modelParamGroup =
"")
137 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart", 10);
139 Preconditioner precond(A, s.relaxation());
141 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
142 MatrixAdapter operatorA(A);
144 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
148 Dune::InverseOperatorResult result;
149 solver.apply(x, bTmp, result);
151 return result.converged;
154#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
156 template<
class Preconditioner,
class Solver,
class Matrix,
class Vector>
157 static bool solveWithParamTree(
const Matrix& A, Vector& x,
const Vector& b,
158 const Dune::ParameterTree& params)
161 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
162 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
164#if DUNE_VERSION_GTE(DUNE_ISTL,2,8)
165 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub(
"preconditioner"));
167 auto precond = std::make_shared<Preconditioner>(A, params.sub(
"preconditioner"));
169 Solver solver(linearOperator, precond, params);
173 Dune::InverseOperatorResult result;
174 solver.apply(x, bTmp, result);
176 return result.converged;
215 template<
class Matrix,
class Vector>
216 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
218 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
219 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
220 using Solver = Dune::BiCGSTABSolver<Vector>;
222 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
227 return "ILUn preconditioned BiCGSTAB solver";
253 template<
class Matrix,
class Vector>
254 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
256 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
257 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
258 using Solver = Dune::BiCGSTABSolver<Vector>;
260 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
265 return "SOR preconditioned BiCGSTAB solver";
291 template<
class Matrix,
class Vector>
292 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
294 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
295 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
296 using Solver = Dune::BiCGSTABSolver<Vector>;
298 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
303 return "SSOR preconditioned BiCGSTAB solver";
329 template<
class Matrix,
class Vector>
330 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
332 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
333 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
334 using Solver = Dune::BiCGSTABSolver<Vector>;
336 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
341 return "SSOR preconditioned BiCGSTAB solver";
366 template<
class Matrix,
class Vector>
367 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
369 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
370 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
371 using Solver = Dune::BiCGSTABSolver<Vector>;
373 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
378 return "Jac preconditioned BiCGSTAB solver";
403 template<
class Matrix,
class Vector>
404 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
406 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
407 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
408 using Solver = Dune::CGSolver<Vector>;
410 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
415 return "ILUn preconditioned CG solver";
440 template<
class Matrix,
class Vector>
441 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
443 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
444 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
445 using Solver = Dune::CGSolver<Vector>;
447 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
452 return "SOR preconditioned CG solver";
477 template<
class Matrix,
class Vector>
478 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
480 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
481 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
482 using Solver = Dune::CGSolver<Vector>;
484 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
489 return "SSOR preconditioned CG solver";
514 template<
class Matrix,
class Vector>
515 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
517 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
518 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
519 using Solver = Dune::CGSolver<Vector>;
521 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
526 return "GS preconditioned CG solver";
550 template<
class Matrix,
class Vector>
551 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
553 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
554 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
555 using Solver = Dune::CGSolver<Vector>;
557 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
562 return "GS preconditioned CG solver";
588 template<
class Matrix,
class Vector>
589 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
591 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
592 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
593 using Solver = Dune::RestartedGMResSolver<Vector>;
595 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
600 return "SSOR preconditioned GMRes solver";
625 template<
class Matrix,
class Vector>
626 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
628 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
629 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
630 using Solver = Dune::BiCGSTABSolver<Vector>;
632 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
637 return "ILU0 preconditioned BiCGSTAB solver";
661 template<
class Matrix,
class Vector>
662 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
664 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
665 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
666 using Solver = Dune::CGSolver<Vector>;
668 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
673 return "ILU0 preconditioned BiCGSTAB solver";
698 template<
class Matrix,
class Vector>
699 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
701 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
702 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
703 using Solver = Dune::RestartedGMResSolver<Vector>;
705 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
710 return "ILU0 preconditioned BiCGSTAB solver";
736 template<
class Matrix,
class Vector>
737 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
739 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
740 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
741 using Solver = Dune::RestartedGMResSolver<Vector>;
743 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
748 return "ILUn preconditioned GMRes solver";
764 template<
class Matrix,
class Vector>
765 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
768 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
769 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
778 return "Explicit diagonal matrix solver";
791class SuperLUBackend :
public LinearSolver
796 template<
class Matrix,
class Vector>
797 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
799 static_assert(isBCRSMatrix<Matrix>::value,
"SuperLU only works with BCRS matrices!");
800 using BlockType =
typename Matrix::block_type;
801 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
802 constexpr auto blockSize = BlockType::rows;
804 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
807 solver.apply(x, bTmp, result_);
810 for (
int i = 0; i < size; i++)
816 if (isnan(x[i][j]) || isinf(x[i][j]))
818 result_.converged =
false;
824 return result_.converged;
827 std::string name()
const
829 return "SuperLU solver";
832 const Dune::InverseOperatorResult& result()
const
838 Dune::InverseOperatorResult result_;
865class UMFPackBackend :
public LinearSolver
869 UMFPackBackend(
const std::string& paramGroup =
"")
870 : LinearSolver(paramGroup)
872 ordering_ = getParamFromGroup<int>(this->paramGroup(),
"LinearSolver.UMFPackOrdering", 1);
876 void setOrdering(
int i)
881 {
return ordering_; }
883 template<
class Matrix,
class Vector>
884 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
886 static_assert(isBCRSMatrix<Matrix>::value,
"UMFPack only works with BCRS matrices!");
887 using BlockType =
typename Matrix::block_type;
888 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
889 constexpr auto blockSize = BlockType::rows;
891 Dune::UMFPack<Matrix> solver;
892 solver.setVerbosity(this->verbosity() > 0);
893 solver.setOption(UMFPACK_ORDERING, ordering_);
897 solver.apply(x, bTmp, result_);
900 for (
int i = 0; i < size; i++)
906 if (isnan(x[i][j]) || isinf(x[i][j]))
908 result_.converged =
false;
914 return result_.converged;
917 std::string name()
const
919 return "UMFPack solver";
922 const Dune::InverseOperatorResult& result()
const
928 Dune::InverseOperatorResult result_;
943template <
class LinearSolverTraits>
949 template<
class Matrix,
class Vector>
950 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
953 using Solver = Dune::BiCGSTABSolver<Vector>;
955 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
960 return "Uzawa preconditioned BiCGSTAB solver";
968template<
class M,
class X,
class Y,
int blockLevel = 2>
971 template<std::
size_t i>
972 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
974 template<std::
size_t i>
975 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
977 template<std::
size_t i>
978 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
980 using ILUTuple =
typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
1001 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
1004 void pre (X& v, Y& d)
final {}
1008 using namespace Dune::Hybrid;
1009 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](
const auto i)
1011 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
1020 return Dune::SolverCategory::sequential;
1024 template<std::size_t... Is>
1026 : ilu_(std::make_tuple(BlockILU<Is>(m[
Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1046 template<
class Matrix,
class Vector>
1047 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1050 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1051 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1054 solver.apply(x, bTmp, result_);
1056 return result_.converged;
1059 const Dune::InverseOperatorResult&
result()
const
1065 {
return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1068 Dune::InverseOperatorResult result_;
1084 template<
int precondBlockLevel = 2,
class Matrix,
class Vector>
1085 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1088 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1089 static const int restartGMRes = getParamFromGroup<int>(this->
paramGroup(),
"LinearSolver.GMResRestart");
1090 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->
residReduction(), restartGMRes,
1093 solver.apply(x, bTmp, result_);
1095 return result_.converged;
1098 const Dune::InverseOperatorResult&
result()
const
1104 {
return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1107 Dune::InverseOperatorResult result_;
1114template<
class M,
class X,
class Y,
int blockLevel = 2>
1117 template<std::
size_t i>
1118 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1120 template<std::
size_t i>
1121 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1123 template<std::
size_t i>
1124 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1126 template<std::
size_t i>
1127 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1129 template<std::
size_t i>
1130 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1132 template<std::
size_t i>
1133 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1135 using AMGTuple =
typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1154 template<
class LOP,
class Criterion,
class SmootherArgs>
1158 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
1179 using namespace Dune::Hybrid;
1180 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1182 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1199 using namespace Dune::Hybrid;
1200 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1202 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1217 using namespace Dune::Hybrid;
1218 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1220 std::get<decltype(i)::value>(amg_).post(v[i]);
1227 return Dune::SolverCategory::sequential;
1231 template<
class LOP,
class Criterion,
class SmootherArgs, std::size_t... Is>
1233 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1248 template<
class M, std::
size_t i>
1249 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1251 template<
class X, std::
size_t i>
1252 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1254 template<
class M,
class X, std::
size_t i>
1255 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1257 template<
class M,
class X, std::
size_t i>
1258 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1260 template<
class M, std::
size_t i>
1261 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1263 template<
class M,
class X, std::
size_t i>
1264 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1270 template<
class Matrix,
class Vector>
1271 bool solve(
const Matrix& m, Vector& x,
const Vector& b)
1275 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1276 params.setDefaultValuesIsotropic(3);
1277 params.setDebugLevel(this->
verbosity());
1279 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1280 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1282 using namespace Dune::Hybrid;
1283 forEach(std::make_index_sequence<Matrix::N()>{}, [&](
const auto i)
1285 auto& args = std::get<decltype(i)::value>(smootherArgs);
1286 args->iterations = 1;
1287 args->relaxationFactor = 1;
1290 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1294 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1295 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1298 solver.apply(x, bTmp, result_);
1300 return result_.converged;
1303 const Dune::InverseOperatorResult&
result()
const
1309 {
return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1312 template<
template<
class M, std::
size_t i>
class Criterion,
class Matrix,
class Params, std::size_t... Is>
1313 auto makeCriterion_(
const Params& p, std::index_sequence<Is...>)
1315 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1318 template<
template<
class M,
class X, std::
size_t i>
class SmootherArgs,
class Matrix,
class Vector, std::size_t... Is>
1319 auto makeSmootherArgs_(std::index_sequence<Is...>)
1321 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1324 template<
template<
class M,
class X, std::
size_t i>
class LinearOperator,
class Matrix,
class Vector, std::size_t... Is>
1325 auto makeLinearOperator_(
const Matrix& m, std::index_sequence<Is...>)
1327 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1330 Dune::InverseOperatorResult result_;
Dumux preconditioners for iterative solvers.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Base class for linear solvers.
Provides a parallel linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Utilities for template meta programming.
Type traits to be used with matrix types.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:188
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:185
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
Definition: common/pdesolver.hh:36
Helper type to determine whether a given type is a Dune::MultiTypeBlockMatrix.
Definition: matrix.hh:49
Definition: utility.hh:40
static Dune::ParameterTree createParameterTree(const std::string ¶mGroup="")
Create a tree containing parameters required for the linear solvers and precondioners of the Dune IST...
Definition: linearsolverparameters.hh:46
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:80
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:66
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:133
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:90
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:113
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:70
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:211
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:216
std::string name() const
Definition: seqsolverbackend.hh:225
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:249
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:254
std::string name() const
Definition: seqsolverbackend.hh:263
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:287
std::string name() const
Definition: seqsolverbackend.hh:301
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:292
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:325
std::string name() const
Definition: seqsolverbackend.hh:339
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:330
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:362
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:367
std::string name() const
Definition: seqsolverbackend.hh:376
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:399
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:404
std::string name() const
Definition: seqsolverbackend.hh:413
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:436
std::string name() const
Definition: seqsolverbackend.hh:450
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:441
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:473
std::string name() const
Definition: seqsolverbackend.hh:487
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:478
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:510
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:515
std::string name() const
Definition: seqsolverbackend.hh:524
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:546
std::string name() const
Definition: seqsolverbackend.hh:560
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:551
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:584
std::string name() const
Definition: seqsolverbackend.hh:598
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:589
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:621
std::string name() const
Definition: seqsolverbackend.hh:635
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:626
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:657
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:662
std::string name() const
Definition: seqsolverbackend.hh:671
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:694
std::string name() const
Definition: seqsolverbackend.hh:708
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:699
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:732
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:737
std::string name() const
Definition: seqsolverbackend.hh:746
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:760
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:765
std::string name() const
Definition: seqsolverbackend.hh:776
A Uzawa preconditioned BiCGSTAB solver for saddle-point problems.
Definition: seqsolverbackend.hh:945
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:950
std::string name() const
Definition: seqsolverbackend.hh:958
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:970
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:990
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:984
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:988
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:1004
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:986
void post(X &) final
Definition: seqsolverbackend.hh:1015
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1018
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:1006
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:998
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1041
std::string name() const
Definition: seqsolverbackend.hh:1064
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1059
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1047
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1079
std::string name() const
Definition: seqsolverbackend.hh:1103
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1098
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1085
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1116
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1215
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1139
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1197
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1225
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1143
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1177
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1145
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1155
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1141
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1247
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1303
std::string name() const
Definition: seqsolverbackend.hh:1308
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1271
Base class for linear solvers.
Definition: solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: solver.hh:99
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:91
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:79
LinearSolver(const std::string ¶mGroup="")
Contruct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:83