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/version.hh>
38#include <dune/common/hybridutilities.hh>
70 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
71 static bool solve(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
72 const std::string& modelParamGroup =
"")
74 Preconditioner precond(A, s.precondIter(), s.relaxation());
77 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
78 MatrixAdapter linearOperator(A);
80 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
84 Dune::InverseOperatorResult result;
85 solver.apply(x, bTmp, result);
87 return result.converged;
90 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
91 static bool solveWithGMRes(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
92 const std::string& modelParamGroup =
"")
95 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart", 10);
97 Preconditioner precond(A, s.precondIter(), s.relaxation());
100 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
101 MatrixAdapter linearOperator(A);
103 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
107 Dune::InverseOperatorResult result;
108 solver.apply(x, bTmp, result);
110 return result.converged;
113 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
114 static bool solveWithILU0Prec(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
115 const std::string& modelParamGroup =
"")
117 Preconditioner precond(A, s.relaxation());
119 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
120 MatrixAdapter operatorA(A);
122 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
126 Dune::InverseOperatorResult result;
127 solver.apply(x, bTmp, result);
129 return result.converged;
133 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
135 const std::string& modelParamGroup =
"")
138 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart", 10);
140 Preconditioner precond(A, s.relaxation());
142 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
143 MatrixAdapter operatorA(A);
145 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
149 Dune::InverseOperatorResult result;
150 solver.apply(x, bTmp, result);
152 return result.converged;
155#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
157 template<
class Preconditioner,
class Solver,
class Matrix,
class Vector>
158 static bool solveWithParamTree(
const Matrix& A, Vector& x,
const Vector& b,
159 const Dune::ParameterTree& params)
162 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
163 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
165#if DUNE_VERSION_GT(DUNE_ISTL,2,7)
166 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub(
"preconditioner"));
168 auto precond = std::make_shared<Preconditioner>(A, params.sub(
"preconditioner"));
170 Solver solver(linearOperator, precond, params);
174 Dune::InverseOperatorResult result;
175 solver.apply(x, bTmp, result);
177 return result.converged;
216 template<
class Matrix,
class Vector>
217 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
219 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
220 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
221 using Solver = Dune::BiCGSTABSolver<Vector>;
223 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
228 return "ILUn preconditioned BiCGSTAB solver";
254 template<
class Matrix,
class Vector>
255 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
257 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
258 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
259 using Solver = Dune::BiCGSTABSolver<Vector>;
261 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
266 return "SOR preconditioned BiCGSTAB solver";
292 template<
class Matrix,
class Vector>
293 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
295 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
296 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
297 using Solver = Dune::BiCGSTABSolver<Vector>;
299 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
304 return "SSOR preconditioned BiCGSTAB solver";
330 template<
class Matrix,
class Vector>
331 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
333 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
334 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
335 using Solver = Dune::BiCGSTABSolver<Vector>;
337 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
342 return "SSOR preconditioned BiCGSTAB solver";
367 template<
class Matrix,
class Vector>
368 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
370 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
371 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
372 using Solver = Dune::BiCGSTABSolver<Vector>;
374 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
379 return "Jac preconditioned BiCGSTAB solver";
404 template<
class Matrix,
class Vector>
405 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
407 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
408 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
409 using Solver = Dune::CGSolver<Vector>;
411 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
416 return "ILUn preconditioned CG solver";
441 template<
class Matrix,
class Vector>
442 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
444 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
445 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
446 using Solver = Dune::CGSolver<Vector>;
448 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
453 return "SOR preconditioned CG solver";
478 template<
class Matrix,
class Vector>
479 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
481 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
482 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
483 using Solver = Dune::CGSolver<Vector>;
485 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
490 return "SSOR preconditioned CG solver";
515 template<
class Matrix,
class Vector>
516 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
518 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
519 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
520 using Solver = Dune::CGSolver<Vector>;
522 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
527 return "GS preconditioned CG solver";
551 template<
class Matrix,
class Vector>
552 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
554 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
555 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
556 using Solver = Dune::CGSolver<Vector>;
558 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
563 return "GS preconditioned CG solver";
589 template<
class Matrix,
class Vector>
590 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
592 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
593 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
594 using Solver = Dune::RestartedGMResSolver<Vector>;
596 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
601 return "SSOR preconditioned GMRes solver";
626 template<
class Matrix,
class Vector>
627 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
629 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
630 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
631 using Solver = Dune::BiCGSTABSolver<Vector>;
633 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
638 return "ILU0 preconditioned BiCGSTAB solver";
662 template<
class Matrix,
class Vector>
663 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
665 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
666 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
667 using Solver = Dune::CGSolver<Vector>;
669 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
674 return "ILU0 preconditioned BiCGSTAB solver";
699 template<
class Matrix,
class Vector>
700 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
702 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
703 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
704 using Solver = Dune::RestartedGMResSolver<Vector>;
706 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
711 return "ILU0 preconditioned BiCGSTAB solver";
737 template<
class Matrix,
class Vector>
738 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
740 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
741 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
742 using Solver = Dune::RestartedGMResSolver<Vector>;
744 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
749 return "ILUn preconditioned GMRes solver";
765 template<
class Matrix,
class Vector>
766 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
769 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
770 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
779 return "Explicit diagonal matrix solver";
792class SuperLUBackend :
public LinearSolver
797 template<
class Matrix,
class Vector>
798 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
800 static_assert(isBCRSMatrix<Matrix>::value,
"SuperLU only works with BCRS matrices!");
801 using BlockType =
typename Matrix::block_type;
802 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
803 constexpr auto blockSize = BlockType::rows;
805 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
808 solver.apply(x, bTmp, result_);
811 for (
int i = 0; i < size; i++)
813 for (
int j = 0; j < blockSize; j++)
817 if (isnan(x[i][j]) || isinf(x[i][j]))
819 result_.converged =
false;
825 return result_.converged;
828 std::string name()
const
830 return "SuperLU solver";
833 const Dune::InverseOperatorResult& result()
const
839 Dune::InverseOperatorResult result_;
852class UMFPackBackend :
public LinearSolver
857 template<
class Matrix,
class Vector>
858 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
860 static_assert(isBCRSMatrix<Matrix>::value,
"UMFPack only works with BCRS matrices!");
861 using BlockType =
typename Matrix::block_type;
862 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
863 constexpr auto blockSize = BlockType::rows;
865 Dune::UMFPack<Matrix> solver(A, this->verbosity() > 0);
868 solver.apply(x, bTmp, result_);
871 for (
int i = 0; i < size; i++)
873 for (
int j = 0; j < blockSize; j++)
877 if (isnan(x[i][j]) || isinf(x[i][j]))
879 result_.converged =
false;
885 return result_.converged;
888 std::string name()
const
890 return "UMFPack solver";
893 const Dune::InverseOperatorResult& result()
const
899 Dune::InverseOperatorResult result_;
909#if DUNE_VERSION_GTE(DUNE_ISTL,2,7)
914template <
class LinearSolverTraits>
915class UzawaBiCGSTABBackend :
public LinearSolver
920 template<
class Matrix,
class Vector>
921 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
923 using Preconditioner = SeqUzawa<Matrix, Vector, Vector>;
924 using Solver = Dune::BiCGSTABSolver<Vector>;
926 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
929 std::string name()
const
931 return "Uzawa preconditioned BiCGSTAB solver";
940template<
class M,
class X,
class Y,
int blockLevel = 2>
943 template<std::
size_t i>
944 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
946 template<std::
size_t i>
947 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
949#if DUNE_VERSION_NEWER(DUNE_ISTL,2,6)
950 template<std::
size_t i>
951 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
953 template<std::
size_t i>
954 using BlockILU = Dune::SeqILU0<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
957 using ILUTuple =
typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
978 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
981 void pre (X& v, Y& d)
final {}
985 using namespace Dune::Hybrid;
986 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](
const auto i)
988 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
995 Dune::SolverCategory::Category
category() const final
997 return Dune::SolverCategory::sequential;
1001 template<std::size_t... Is>
1003 : ilu_(std::make_tuple(BlockILU<Is>(m[
Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1024 template<
class Matrix,
class Vector>
1025 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1028 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1029 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1032 solver.apply(x, bTmp, result_);
1034 return result_.converged;
1037 const Dune::InverseOperatorResult&
result()
const
1043 {
return "block-diagonal ILU0 preconditioned BiCGSTAB solver"; }
1046 Dune::InverseOperatorResult result_;
1053template<
class M,
class X,
class Y,
int blockLevel = 2>
1056 template<std::
size_t i>
1057 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1059 template<std::
size_t i>
1060 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1062 template<std::
size_t i>
1063 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1065 template<std::
size_t i>
1066 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1068 template<std::
size_t i>
1069 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1071 template<std::
size_t i>
1072 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1074 using AMGTuple =
typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1093 template<
class LOP,
class Criterion,
class SmootherArgs>
1097 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
1118 using namespace Dune::Hybrid;
1119 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1121 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1138 using namespace Dune::Hybrid;
1139 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1141 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1156 using namespace Dune::Hybrid;
1157 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1159 std::get<decltype(i)::value>(amg_).post(v[i]);
1166 return Dune::SolverCategory::sequential;
1170 template<
class LOP,
class Criterion,
class SmootherArgs, std::size_t... Is>
1172 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1187 template<
class M, std::
size_t i>
1188 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1190 template<
class X, std::
size_t i>
1191 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1193 template<
class M,
class X, std::
size_t i>
1194 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1196 template<
class M,
class X, std::
size_t i>
1197 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1199 template<
class M, std::
size_t i>
1200 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1202 template<
class M,
class X, std::
size_t i>
1203 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1209 template<
class Matrix,
class Vector>
1210 bool solve(
const Matrix& m, Vector& x,
const Vector& b)
1214 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1215 params.setDefaultValuesIsotropic(3);
1216 params.setDebugLevel(this->
verbosity());
1218 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1219 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1221 using namespace Dune::Hybrid;
1222 forEach(std::make_index_sequence<Matrix::N()>{}, [&](
const auto i)
1224 auto& args = std::get<decltype(i)::value>(smootherArgs);
1225 args->iterations = 1;
1226 args->relaxationFactor = 1;
1229 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1233 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1234 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1237 solver.apply(x, bTmp, result_);
1239 return result_.converged;
1242 const Dune::InverseOperatorResult&
result()
const
1248 {
return "block-diagonal AMG preconditioned BiCGSTAB solver"; }
1251 template<
template<
class M, std::
size_t i>
class Criterion,
class Matrix,
class Params, std::size_t... Is>
1252 auto makeCriterion_(
const Params& p, std::index_sequence<Is...>)
1254 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1257 template<
template<
class M,
class X, std::
size_t i>
class SmootherArgs,
class Matrix,
class Vector, std::size_t... Is>
1258 auto makeSmootherArgs_(std::index_sequence<Is...>)
1260 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1263 template<
template<
class M,
class X, std::
size_t i>
class LinearOperator,
class Matrix,
class Vector, std::size_t... Is>
1264 auto makeLinearOperator_(
const Matrix& m, std::index_sequence<Is...>)
1266 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1269 Dune::InverseOperatorResult result_;
Type traits to be used with matrix types.
Utilities for template meta programming.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Provides a parallel linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Dumux preconditioners for iterative solvers.
Base class for linear solvers.
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:189
Definition: common/pdesolver.hh:35
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 general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:67
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:134
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:91
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:114
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:71
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:212
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:217
std::string name() const
Definition: seqsolverbackend.hh:226
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:250
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:255
std::string name() const
Definition: seqsolverbackend.hh:264
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:288
std::string name() const
Definition: seqsolverbackend.hh:302
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:293
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:326
std::string name() const
Definition: seqsolverbackend.hh:340
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:331
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:363
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:368
std::string name() const
Definition: seqsolverbackend.hh:377
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:400
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:405
std::string name() const
Definition: seqsolverbackend.hh:414
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:437
std::string name() const
Definition: seqsolverbackend.hh:451
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:442
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:474
std::string name() const
Definition: seqsolverbackend.hh:488
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:479
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:511
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:516
std::string name() const
Definition: seqsolverbackend.hh:525
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:547
std::string name() const
Definition: seqsolverbackend.hh:561
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:552
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:585
std::string name() const
Definition: seqsolverbackend.hh:599
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:590
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:622
std::string name() const
Definition: seqsolverbackend.hh:636
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:627
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:658
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:663
std::string name() const
Definition: seqsolverbackend.hh:672
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:695
std::string name() const
Definition: seqsolverbackend.hh:709
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:700
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:733
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:738
std::string name() const
Definition: seqsolverbackend.hh:747
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:761
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:766
std::string name() const
Definition: seqsolverbackend.hh:777
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:942
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:967
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:961
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:965
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:981
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:963
void post(X &) final
Definition: seqsolverbackend.hh:992
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:995
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:983
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:975
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1018
std::string name() const
Definition: seqsolverbackend.hh:1042
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1037
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1025
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1055
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1154
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1078
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1136
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1164
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1082
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1116
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1084
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1094
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1080
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1186
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1242
std::string name() const
Definition: seqsolverbackend.hh:1247
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1210
Base class for linear solvers.
Definition: solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: solver.hh:131
int maxIter() const
the maximum number of linear solver iterations
Definition: solver.hh:123
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:111
LinearSolver(const std::string ¶mGroup="")
Contruct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:115