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/common/version.hh>
36#include <dune/common/hybridutilities.hh>
66 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
67 static bool solve(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
68 const std::string& modelParamGroup =
"")
70 Preconditioner precond(A, s.precondIter(), s.relaxation());
73 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
74 MatrixAdapter linearOperator(A);
76 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
80 Dune::InverseOperatorResult result;
81 solver.apply(x, bTmp, result);
83 return result.converged;
86 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
87 static bool solveWithGMRes(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
88 const std::string& modelParamGroup =
"")
91 const int restartGMRes = getParamFromGroup<double>(modelParamGroup,
"LinearSolver.GMResRestart");
93 Preconditioner precond(A, s.precondIter(), s.relaxation());
96 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
97 MatrixAdapter linearOperator(A);
99 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
103 Dune::InverseOperatorResult result;
104 solver.apply(x, bTmp, result);
106 return result.converged;
109 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
110 static bool solveWithILU0Prec(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
111 const std::string& modelParamGroup =
"")
113 Preconditioner precond(A, s.relaxation());
115 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
116 MatrixAdapter operatorA(A);
118 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
122 Dune::InverseOperatorResult result;
123 solver.apply(x, bTmp, result);
125 return result.converged;
129 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
131 const std::string& modelParamGroup =
"")
134 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart");
136 Preconditioner precond(A, s.relaxation());
138 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
139 MatrixAdapter operatorA(A);
141 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
145 Dune::InverseOperatorResult result;
146 solver.apply(x, bTmp, result);
148 return result.converged;
175 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
176 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
178 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
179 using Solver = Dune::BiCGSTABSolver<Vector>;
181 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
186 return "ILUn preconditioned BiCGSTAB solver";
213 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
214 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
216 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
217 using Solver = Dune::BiCGSTABSolver<Vector>;
219 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
224 return "SOR preconditioned BiCGSTAB solver";
251 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
252 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
254 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
255 using Solver = Dune::BiCGSTABSolver<Vector>;
257 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
262 return "SSOR preconditioned BiCGSTAB solver";
289 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
290 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
292 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
293 using Solver = Dune::BiCGSTABSolver<Vector>;
295 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
300 return "SSOR preconditioned BiCGSTAB solver";
326 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
327 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
329 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
330 using Solver = Dune::BiCGSTABSolver<Vector>;
332 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
337 return "Jac preconditioned BiCGSTAB solver";
363 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
364 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
366 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
367 using Solver = Dune::CGSolver<Vector>;
369 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
374 return "ILUn preconditioned CG solver";
400 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
401 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
403 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
404 using Solver = Dune::CGSolver<Vector>;
406 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
411 return "SOR preconditioned CG solver";
437 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
438 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
440 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
441 using Solver = Dune::CGSolver<Vector>;
443 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
448 return "SSOR preconditioned CG solver";
474 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
475 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
477 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
478 using Solver = Dune::CGSolver<Vector>;
480 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
485 return "GS preconditioned CG solver";
510 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
511 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
513 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
514 using Solver = Dune::CGSolver<Vector>;
516 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
521 return "GS preconditioned CG solver";
548 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
549 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
551 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
552 using Solver = Dune::RestartedGMResSolver<Vector>;
554 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
559 return "SSOR preconditioned GMRes solver";
585 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
586 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
588 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
589 using Solver = Dune::BiCGSTABSolver<Vector>;
591 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
596 return "ILU0 preconditioned BiCGSTAB solver";
621 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
622 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
624 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
625 using Solver = Dune::CGSolver<Vector>;
627 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
632 return "ILU0 preconditioned BiCGSTAB solver";
658 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
659 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
661 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
662 using Solver = Dune::RestartedGMResSolver<Vector>;
664 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
669 return "ILU0 preconditioned BiCGSTAB solver";
696 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
697 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
699 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
700 using Solver = Dune::RestartedGMResSolver<Vector>;
702 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
707 return "ILUn preconditioned GMRes solver";
724 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
725 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
728 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
737 return "Explicit diagonal matrix solver";
750class SuperLUBackend :
public LinearSolver
756 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
757 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
759 static_assert(isBCRSMatrix<Matrix>::value,
"SuperLU only works with BCRS matrices!");
760 using BlockType =
typename Matrix::block_type;
761 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
762 constexpr auto blockSize = BlockType::rows;
764 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
767 solver.apply(x, bTmp, result_);
770 for (
int i = 0; i < size; i++)
772 for (
int j = 0; j < blockSize; j++)
776 if (isnan(x[i][j]) || isinf(x[i][j]))
778 result_.converged =
false;
784 return result_.converged;
787 std::string name()
const
789 return "SuperLU solver";
792 const Dune::InverseOperatorResult& result()
const
798 Dune::InverseOperatorResult result_;
811class UMFPackBackend :
public LinearSolver
817 template<
int precondBlockLevel = 1,
class Matrix,
class Vector>
818 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
820 static_assert(isBCRSMatrix<Matrix>::value,
"UMFPack only works with BCRS matrices!");
821 using BlockType =
typename Matrix::block_type;
822 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
823 constexpr auto blockSize = BlockType::rows;
825 Dune::UMFPack<Matrix> solver(A, this->verbosity() > 0);
828 solver.apply(x, bTmp, result_);
831 for (
int i = 0; i < size; i++)
833 for (
int j = 0; j < blockSize; j++)
837 if (isnan(x[i][j]) || isinf(x[i][j]))
839 result_.converged =
false;
845 return result_.converged;
848 std::string name()
const
850 return "UMFPack solver";
853 const Dune::InverseOperatorResult& result()
const
859 Dune::InverseOperatorResult result_;
873template<
class M,
class X,
class Y,
int blockLevel = 2>
876 template<std::
size_t i>
877 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
879 template<std::
size_t i>
880 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
882#if DUNE_VERSION_NEWER(DUNE_ISTL,2,6)
883 template<std::
size_t i>
884 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
886 template<std::
size_t i>
887 using BlockILU = Dune::SeqILU0<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
890 using ILUTuple =
typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::size()> >::type;
911 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
914 void pre (X& v, Y& d)
final {}
918 using namespace Dune::Hybrid;
919 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](
const auto i)
921 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
928 Dune::SolverCategory::Category
category() const final
930 return Dune::SolverCategory::sequential;
934 template<std::size_t... Is>
936 : ilu_(std::make_tuple(BlockILU<Is>(m[
Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
957 template<
int precondBlockLevel = 2,
class Matrix,
class Vector>
958 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
961 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
962 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
965 solver.apply(x, bTmp, result_);
967 return result_.converged;
970 const Dune::InverseOperatorResult&
result()
const
976 {
return "block-diagonal ILU0 preconditioned BiCGSTAB solver"; }
979 Dune::InverseOperatorResult result_;
986template<
class M,
class X,
class Y,
int blockLevel = 2>
989 template<std::
size_t i>
990 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
992 template<std::
size_t i>
993 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
995 template<std::
size_t i>
996 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
998 template<std::
size_t i>
999 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1001 template<std::
size_t i>
1002 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1004 template<std::
size_t i>
1005 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1007 using AMGTuple =
typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::size()> >::type;
1026 template<
class LOP,
class Criterion,
class SmootherArgs>
1030 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
1051 using namespace Dune::Hybrid;
1052 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1054 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1071 using namespace Dune::Hybrid;
1072 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1074 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1089 using namespace Dune::Hybrid;
1090 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1092 std::get<decltype(i)::value>(amg_).post(v[i]);
1099 return Dune::SolverCategory::sequential;
1103 template<
class LOP,
class Criterion,
class SmootherArgs, std::size_t... Is>
1105 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1120 template<
class M, std::
size_t i>
1121 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1123 template<
class X, std::
size_t i>
1124 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1126 template<
class M,
class X, std::
size_t i>
1127 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1129 template<
class M,
class X, std::
size_t i>
1130 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1132 template<
class M, std::
size_t i>
1133 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1135 template<
class M,
class X, std::
size_t i>
1136 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1142 template<
int precondBlockLevel = 2,
class Matrix,
class Vector>
1143 bool solve(
const Matrix& m, Vector& x,
const Vector& b)
1147 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1148 params.setDefaultValuesIsotropic(3);
1149 params.setDebugLevel(this->
verbosity());
1151 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::size()>{});
1152 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::size()>{});
1154 using namespace Dune::Hybrid;
1155 forEach(integralRange(Dune::Hybrid::size(m)), [&](
const auto i)
1157 auto& args = std::get<decltype(i)::value>(smootherArgs);
1158 args->iterations = 1;
1159 args->relaxationFactor = 1;
1162 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::size()>{});
1166 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1167 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1170 solver.apply(x, bTmp, result_);
1172 return result_.converged;
1175 const Dune::InverseOperatorResult&
result()
const
1181 {
return "block-diagonal AMG preconditioned BiCGSTAB solver"; }
1184 template<
template<
class M, std::
size_t i>
class Criterion,
class Matrix,
class Params, std::size_t... Is>
1185 auto makeCriterion_(
const Params& p, std::index_sequence<Is...>)
1187 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1190 template<
template<
class M,
class X, std::
size_t i>
class SmootherArgs,
class Matrix,
class Vector, std::size_t... Is>
1191 auto makeSmootherArgs_(std::index_sequence<Is...>)
1193 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1196 template<
template<
class M,
class X, std::
size_t i>
class LinearOperator,
class Matrix,
class Vector, std::size_t... Is>
1197 auto makeLinearOperator_(
const Matrix& m, std::index_sequence<Is...>)
1199 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1202 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.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: common/properties/model.hh:34
Definition: utility.hh:40
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:63
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:130
static bool solveWithGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:87
static bool solveWithILU0Prec(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:110
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:67
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:170
std::string name() const
Definition: seqsolverbackend.hh:184
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:176
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:208
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:214
std::string name() const
Definition: seqsolverbackend.hh:222
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:246
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:252
std::string name() const
Definition: seqsolverbackend.hh:260
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:284
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:290
std::string name() const
Definition: seqsolverbackend.hh:298
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:321
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:327
std::string name() const
Definition: seqsolverbackend.hh:335
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:358
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:364
std::string name() const
Definition: seqsolverbackend.hh:372
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:395
std::string name() const
Definition: seqsolverbackend.hh:409
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:401
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:432
std::string name() const
Definition: seqsolverbackend.hh:446
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:438
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:469
std::string name() const
Definition: seqsolverbackend.hh:483
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:475
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:505
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:511
std::string name() const
Definition: seqsolverbackend.hh:519
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:543
std::string name() const
Definition: seqsolverbackend.hh:557
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:549
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:580
std::string name() const
Definition: seqsolverbackend.hh:594
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:586
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:616
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:622
std::string name() const
Definition: seqsolverbackend.hh:630
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:653
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:659
std::string name() const
Definition: seqsolverbackend.hh:667
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:691
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:697
std::string name() const
Definition: seqsolverbackend.hh:705
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:719
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:725
std::string name() const
Definition: seqsolverbackend.hh:735
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:875
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:900
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:894
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:898
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:914
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:896
void post(X &) final
Definition: seqsolverbackend.hh:925
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:928
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:916
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:908
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:951
std::string name() const
Definition: seqsolverbackend.hh:975
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:970
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:958
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:988
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1087
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1011
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1069
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1097
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1015
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1049
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1017
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1027
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1013
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1119
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1143
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1175
std::string name() const
Definition: seqsolverbackend.hh:1180
Base class for linear solvers.
Definition: dumux/linear/solver.hh:37
double residReduction() const
the linear solver residual reduction
Definition: dumux/linear/solver.hh:97
int maxIter() const
the maximum number of linear solver iterations
Definition: dumux/linear/solver.hh:89
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: dumux/linear/solver.hh:77
LinearSolver(const std::string ¶mGroup="")
Contruct the solver.
Definition: dumux/linear/solver.hh:52
int verbosity() const
the verbosity level
Definition: dumux/linear/solver.hh:81
Base class for linear solvers.