12#ifndef DUMUX_SEQ_SOLVER_BACKEND_HH
13#define DUMUX_SEQ_SOLVER_BACKEND_HH
19#include <dune/istl/preconditioners.hh>
20#include <dune/istl/solvers.hh>
21#include <dune/istl/superlu.hh>
22#include <dune/istl/umfpack.hh>
23#include <dune/istl/io.hh>
24#include <dune/common/indices.hh>
25#include <dune/common/hybridutilities.hh>
57 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
58 static bool solve(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
59 const std::string& modelParamGroup =
"")
61 Preconditioner precond(A, s.precondIter(), s.relaxation());
64 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
65 MatrixAdapter linearOperator(A);
67 Solver solver(linearOperator, precond, s.residReduction(), s.maxIter(), s.verbosity());
71 Dune::InverseOperatorResult result;
72 solver.apply(x, bTmp, result);
74 return result.converged;
77 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
78 static bool solveWithGMRes(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
79 const std::string& modelParamGroup =
"")
82 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart", 10);
84 Preconditioner precond(A, s.precondIter(), s.relaxation());
87 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
88 MatrixAdapter linearOperator(A);
90 Solver solver(linearOperator, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
94 Dune::InverseOperatorResult result;
95 solver.apply(x, bTmp, result);
97 return result.converged;
100 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
101 static bool solveWithILU0Prec(
const SolverInterface& s,
const Matrix& A, Vector& x,
const Vector& b,
102 const std::string& modelParamGroup =
"")
104 Preconditioner precond(A, s.relaxation());
106 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
107 MatrixAdapter operatorA(A);
109 Solver solver(operatorA, precond, s.residReduction(), s.maxIter(), s.verbosity());
113 Dune::InverseOperatorResult result;
114 solver.apply(x, bTmp, result);
116 return result.converged;
120 template<
class Preconditioner,
class Solver,
class SolverInterface,
class Matrix,
class Vector>
122 const std::string& modelParamGroup =
"")
125 const int restartGMRes = getParamFromGroup<int>(modelParamGroup,
"LinearSolver.GMResRestart", 10);
127 Preconditioner precond(A, s.relaxation());
129 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
130 MatrixAdapter operatorA(A);
132 Solver solver(operatorA, precond, s.residReduction(), restartGMRes, s.maxIter(), s.verbosity());
136 Dune::InverseOperatorResult result;
137 solver.apply(x, bTmp, result);
139 return result.converged;
143 template<
class Preconditioner,
class Solver,
class Matrix,
class Vector>
145 const Dune::ParameterTree& params)
148 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
149 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
151 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub(
"preconditioner"));
152 Solver solver(linearOperator, precond, params);
156 Dune::InverseOperatorResult result;
157 solver.apply(x, bTmp, result);
159 return result.converged;
197 template<
class Matrix,
class Vector>
198 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
200 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
201 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
202 using Solver = Dune::BiCGSTABSolver<Vector>;
204 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
209 return "ILUn preconditioned BiCGSTAB solver";
235 template<
class Matrix,
class Vector>
236 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
238 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
239 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
240 using Solver = Dune::BiCGSTABSolver<Vector>;
242 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
247 return "SOR preconditioned BiCGSTAB solver";
273 template<
class Matrix,
class Vector>
274 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
276 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
277 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
278 using Solver = Dune::BiCGSTABSolver<Vector>;
280 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
285 return "SSOR preconditioned BiCGSTAB solver";
311 template<
class Matrix,
class Vector>
312 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
314 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
315 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
316 using Solver = Dune::BiCGSTABSolver<Vector>;
318 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
323 return "SSOR preconditioned BiCGSTAB solver";
348 template<
class Matrix,
class Vector>
349 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
351 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
352 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
353 using Solver = Dune::BiCGSTABSolver<Vector>;
355 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
360 return "Jac preconditioned BiCGSTAB solver";
385 template<
class Matrix,
class Vector>
386 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
388 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
389 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
390 using Solver = Dune::CGSolver<Vector>;
392 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
397 return "ILUn preconditioned CG solver";
422 template<
class Matrix,
class Vector>
423 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
425 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
426 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
427 using Solver = Dune::CGSolver<Vector>;
429 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
434 return "SOR preconditioned CG solver";
459 template<
class Matrix,
class Vector>
460 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
462 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
463 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
464 using Solver = Dune::CGSolver<Vector>;
466 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
471 return "SSOR preconditioned CG solver";
496 template<
class Matrix,
class Vector>
497 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
499 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
500 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
501 using Solver = Dune::CGSolver<Vector>;
503 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
508 return "GS preconditioned CG solver";
532 template<
class Matrix,
class Vector>
533 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
535 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
536 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
537 using Solver = Dune::CGSolver<Vector>;
539 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
544 return "GS preconditioned CG solver";
570 template<
class Matrix,
class Vector>
571 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
573 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
574 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
575 using Solver = Dune::RestartedGMResSolver<Vector>;
577 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
582 return "SSOR preconditioned GMRes solver";
607 template<
class Matrix,
class Vector>
608 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
610 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
611 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
612 using Solver = Dune::BiCGSTABSolver<Vector>;
614 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
619 return "ILU0 preconditioned BiCGSTAB solver";
643 template<
class Matrix,
class Vector>
644 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
646 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
647 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
648 using Solver = Dune::CGSolver<Vector>;
650 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
655 return "ILU0 preconditioned BiCGSTAB solver";
680 template<
class Matrix,
class Vector>
681 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
683 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
684 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
685 using Solver = Dune::RestartedGMResSolver<Vector>;
687 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
692 return "ILU0 preconditioned BiCGSTAB solver";
718 template<
class Matrix,
class Vector>
719 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
721 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
722 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
723 using Solver = Dune::RestartedGMResSolver<Vector>;
725 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->paramGroup());
730 return "ILUn preconditioned GMRes solver";
746 template<
class Matrix,
class Vector>
747 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
750 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
751 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
760 return "Explicit diagonal matrix solver";
773class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] SuperLUBackend :
public LinearSolver
778 template<
class Matrix,
class Vector>
779 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
781 static_assert(isBCRSMatrix<Matrix>::value,
"SuperLU only works with BCRS matrices!");
782 using BlockType =
typename Matrix::block_type;
783 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
784 constexpr auto blockSize = BlockType::rows;
786 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
789 solver.apply(x, bTmp, result_);
792 for (
int i = 0; i < size; i++)
798 if (isnan(x[i][j]) || isinf(x[i][j]))
800 result_.converged =
false;
806 return result_.converged;
809 std::string name()
const
811 return "SuperLU solver";
814 const Dune::InverseOperatorResult& result()
const
820 Dune::InverseOperatorResult result_;
847class [[deprecated("Removed after 3.7. Use solver from istlsolvers.hh")]] UMFPackBackend :
public LinearSolver
851 UMFPackBackend(
const std::string& paramGroup =
"")
852 : LinearSolver(paramGroup)
854 ordering_ = getParamFromGroup<int>(this->paramGroup(),
"LinearSolver.UMFPackOrdering", 1);
858 void setOrdering(
int i)
863 {
return ordering_; }
865 template<
class Matrix,
class Vector>
866 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
868 static_assert(isBCRSMatrix<Matrix>::value,
"UMFPack only works with BCRS matrices!");
869 using BlockType =
typename Matrix::block_type;
870 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
871 constexpr auto blockSize = BlockType::rows;
873 Dune::UMFPack<Matrix> solver;
874 solver.setVerbosity(this->verbosity() > 0);
875 solver.setOption(UMFPACK_ORDERING, ordering_);
879 solver.apply(x, bTmp, result_);
882 for (
int i = 0; i < size; i++)
888 if (isnan(x[i][j]) || isinf(x[i][j]))
890 result_.converged =
false;
896 return result_.converged;
899 std::string name()
const
901 return "UMFPack solver";
904 const Dune::InverseOperatorResult& result()
const
910 Dune::InverseOperatorResult result_;
925template <
class LinearSolverTraits>
931 template<
class Matrix,
class Vector>
932 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
935 using Solver = Dune::BiCGSTABSolver<Vector>;
937 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
942 return "Uzawa preconditioned BiCGSTAB solver";
950template<
class M,
class X,
class Y,
int blockLevel = 2>
953 template<std::
size_t i>
954 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
956 template<std::
size_t i>
957 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
959 template<std::
size_t i>
960 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
962 using ILUTuple =
typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
983 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
986 void pre (X& v, Y& d)
final {}
990 using namespace Dune::Hybrid;
991 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](
const auto i)
993 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
1002 return Dune::SolverCategory::sequential;
1006 template<std::size_t... Is>
1008 : ilu_(std::make_tuple(BlockILU<Is>(m[
Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1028 template<
class Matrix,
class Vector>
1029 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1033 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1036 solver.apply(x, bTmp, result_);
1038 return result_.converged;
1041 const Dune::InverseOperatorResult&
result()
const
1047 {
return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1050 Dune::InverseOperatorResult result_;
1066 template<
int precondBlockLevel = 2,
class Matrix,
class Vector>
1067 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1071 static const int restartGMRes = getParamFromGroup<int>(this->
paramGroup(),
"LinearSolver.GMResRestart");
1072 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->
residReduction(), restartGMRes,
1075 solver.apply(x, bTmp, result_);
1077 return result_.converged;
1080 const Dune::InverseOperatorResult&
result()
const
1086 {
return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1089 Dune::InverseOperatorResult result_;
1096template<
class M,
class X,
class Y,
int blockLevel = 2>
1099 template<std::
size_t i>
1100 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1102 template<std::
size_t i>
1103 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1105 template<std::
size_t i>
1106 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1108 template<std::
size_t i>
1109 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1111 template<std::
size_t i>
1112 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1114 template<std::
size_t i>
1115 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1117 using AMGTuple =
typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1136 template<
class LOP,
class Criterion,
class SmootherArgs>
1140 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
1161 using namespace Dune::Hybrid;
1162 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1164 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1181 using namespace Dune::Hybrid;
1182 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1184 std::get<decltype(i)::value>(amg_).apply(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_).post(v[i]);
1209 return Dune::SolverCategory::sequential;
1213 template<
class LOP,
class Criterion,
class SmootherArgs, std::size_t... Is>
1215 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1230 template<
class M, std::
size_t i>
1231 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1233 template<
class X, std::
size_t i>
1234 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1236 template<
class M,
class X, std::
size_t i>
1237 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1239 template<
class M,
class X, std::
size_t i>
1240 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1242 template<
class M, std::
size_t i>
1243 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1245 template<
class M,
class X, std::
size_t i>
1246 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1252 template<
class Matrix,
class Vector>
1253 bool solve(
const Matrix& m, Vector& x,
const Vector& b)
1257 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1258 params.setDefaultValuesIsotropic(3);
1259 params.setDebugLevel(this->
verbosity());
1261 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1262 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1264 using namespace Dune::Hybrid;
1265 forEach(std::make_index_sequence<Matrix::N()>{}, [&](
const auto i)
1267 auto& args = std::get<decltype(i)::value>(smootherArgs);
1268 args->iterations = 1;
1269 args->relaxationFactor = 1;
1272 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1277 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1280 solver.apply(x, bTmp, result_);
1282 return result_.converged;
1285 const Dune::InverseOperatorResult&
result()
const
1291 {
return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1294 template<
template<
class M, std::
size_t i>
class Criterion,
class Matrix,
class Params, std::size_t... Is>
1295 auto makeCriterion_(
const Params& p, std::index_sequence<Is...>)
1297 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1300 template<
template<
class M,
class X, std::
size_t i>
class SmootherArgs,
class Matrix,
class Vector, std::size_t... Is>
1301 auto makeSmootherArgs_(std::index_sequence<Is...>)
1303 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1306 template<
template<
class M,
class X, std::
size_t i>
class LinearOperator,
class Matrix,
class Vector, std::size_t... Is>
1307 auto makeLinearOperator_(
const Matrix& m, std::index_sequence<Is...>)
1309 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1312 Dune::InverseOperatorResult result_;
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1229
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1285
std::string name() const
Definition: seqsolverbackend.hh:1290
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1253
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1098
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1197
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1121
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1179
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1207
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1125
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1159
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1127
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1137
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1123
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1023
std::string name() const
Definition: seqsolverbackend.hh:1046
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1041
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1029
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:952
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:972
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:966
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:970
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:986
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:968
void post(X &) final
Definition: seqsolverbackend.hh:997
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1000
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:988
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:980
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1061
std::string name() const
Definition: seqsolverbackend.hh:1085
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1080
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1067
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:742
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:747
std::string name() const
Definition: seqsolverbackend.hh:758
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:307
std::string name() const
Definition: seqsolverbackend.hh:321
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:312
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:492
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:497
std::string name() const
Definition: seqsolverbackend.hh:506
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:603
std::string name() const
Definition: seqsolverbackend.hh:617
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:608
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:639
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:644
std::string name() const
Definition: seqsolverbackend.hh:653
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:676
std::string name() const
Definition: seqsolverbackend.hh:690
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:681
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:193
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:198
std::string name() const
Definition: seqsolverbackend.hh:207
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:381
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:386
std::string name() const
Definition: seqsolverbackend.hh:395
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:714
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:719
std::string name() const
Definition: seqsolverbackend.hh:728
A general solver backend allowing arbitrary preconditioners and solvers.
Definition: seqsolverbackend.hh:54
static bool solveWithILU0PrecGMRes(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:121
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:101
static bool solve(const SolverInterface &s, const Matrix &A, Vector &x, const Vector &b, const std::string &modelParamGroup="")
Definition: seqsolverbackend.hh:58
static bool solveWithParamTree(const Matrix &A, Vector &x, const Vector &b, const Dune::ParameterTree ¶ms)
Definition: seqsolverbackend.hh:144
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:344
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:349
std::string name() const
Definition: seqsolverbackend.hh:358
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:528
std::string name() const
Definition: seqsolverbackend.hh:542
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:533
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 ¶mGroup="")
Construct the solver.
Definition: solver.hh:43
int verbosity() const
the verbosity level
Definition: solver.hh:82
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:47
Adapter to turn a multi-type matrix into a thread-parallel linear operator. Adapts a matrix to the as...
Definition: parallelmatrixadapter.hh:28
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:231
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:236
std::string name() const
Definition: seqsolverbackend.hh:245
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:418
std::string name() const
Definition: seqsolverbackend.hh:432
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:423
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:269
std::string name() const
Definition: seqsolverbackend.hh:283
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:274
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:455
std::string name() const
Definition: seqsolverbackend.hh:469
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:460
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:566
std::string name() const
Definition: seqsolverbackend.hh:580
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:571
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:927
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:932
std::string name() const
Definition: seqsolverbackend.hh:940
constexpr std::size_t preconditionerBlockLevel() noexcept
Returns the block level for the preconditioner for a given matrix.
Definition: seqsolverbackend.hh:170
Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
Type traits to be used with matrix types.
constexpr std::size_t blockSize()
Definition: dunevectors.hh:27
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.