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;
155 template<
class Preconditioner,
class Solver,
class Matrix,
class Vector>
157 const Dune::ParameterTree& params)
160 using MatrixAdapter = Dune::MatrixAdapter<Matrix, Vector, Vector>;
161 const auto linearOperator = std::make_shared<MatrixAdapter>(A);
163 auto precond = std::make_shared<Preconditioner>(linearOperator, params.sub(
"preconditioner"));
164 Solver solver(linearOperator, precond, params);
168 Dune::InverseOperatorResult result;
169 solver.apply(x, bTmp, result);
171 return result.converged;
209 template<
class Matrix,
class Vector>
210 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
212 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
213 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
214 using Solver = Dune::BiCGSTABSolver<Vector>;
216 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
221 return "ILUn preconditioned BiCGSTAB solver";
247 template<
class Matrix,
class Vector>
248 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
250 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
251 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
252 using Solver = Dune::BiCGSTABSolver<Vector>;
254 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
259 return "SOR preconditioned BiCGSTAB solver";
285 template<
class Matrix,
class Vector>
286 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
288 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
289 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
290 using Solver = Dune::BiCGSTABSolver<Vector>;
292 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
297 return "SSOR preconditioned BiCGSTAB solver";
323 template<
class Matrix,
class Vector>
324 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
326 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
327 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
328 using Solver = Dune::BiCGSTABSolver<Vector>;
330 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
335 return "SSOR preconditioned BiCGSTAB solver";
360 template<
class Matrix,
class Vector>
361 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
363 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
364 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
365 using Solver = Dune::BiCGSTABSolver<Vector>;
367 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
372 return "Jac preconditioned BiCGSTAB solver";
397 template<
class Matrix,
class Vector>
398 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
400 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
401 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
402 using Solver = Dune::CGSolver<Vector>;
404 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
409 return "ILUn preconditioned CG solver";
434 template<
class Matrix,
class Vector>
435 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
437 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
438 using Preconditioner = Dune::SeqSOR<Matrix, Vector, Vector, precondBlockLevel>;
439 using Solver = Dune::CGSolver<Vector>;
441 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
446 return "SOR preconditioned CG solver";
471 template<
class Matrix,
class Vector>
472 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
474 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
475 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
476 using Solver = Dune::CGSolver<Vector>;
478 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
483 return "SSOR preconditioned CG solver";
508 template<
class Matrix,
class Vector>
509 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
511 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
512 using Preconditioner = Dune::SeqGS<Matrix, Vector, Vector, precondBlockLevel>;
513 using Solver = Dune::CGSolver<Vector>;
515 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
520 return "GS preconditioned CG solver";
544 template<
class Matrix,
class Vector>
545 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
547 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
548 using Preconditioner = Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel>;
549 using Solver = Dune::CGSolver<Vector>;
551 return IterativePreconditionedSolverImpl::template solve<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
556 return "GS preconditioned CG solver";
582 template<
class Matrix,
class Vector>
583 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
585 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
586 using Preconditioner = Dune::SeqSSOR<Matrix, Vector, Vector, precondBlockLevel>;
587 using Solver = Dune::RestartedGMResSolver<Vector>;
589 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
594 return "SSOR preconditioned GMRes solver";
619 template<
class Matrix,
class Vector>
620 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
622 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
623 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
624 using Solver = Dune::BiCGSTABSolver<Vector>;
626 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
631 return "ILU0 preconditioned BiCGSTAB solver";
655 template<
class Matrix,
class Vector>
656 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
658 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
659 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
660 using Solver = Dune::CGSolver<Vector>;
662 return IterativePreconditionedSolverImpl::template solveWithILU0Prec<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
667 return "ILU0 preconditioned BiCGSTAB solver";
692 template<
class Matrix,
class Vector>
693 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
695 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
696 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
697 using Solver = Dune::RestartedGMResSolver<Vector>;
699 return IterativePreconditionedSolverImpl::template solveWithILU0PrecGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
704 return "ILU0 preconditioned BiCGSTAB solver";
730 template<
class Matrix,
class Vector>
731 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
733 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
734 using Preconditioner = Dune::SeqILU<Matrix, Vector, Vector, precondBlockLevel>;
735 using Solver = Dune::RestartedGMResSolver<Vector>;
737 return IterativePreconditionedSolverImpl::template solveWithGMRes<Preconditioner, Solver>(*
this, A, x, b, this->
paramGroup());
742 return "ILUn preconditioned GMRes solver";
758 template<
class Matrix,
class Vector>
759 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
762 constexpr auto precondBlockLevel = preconditionerBlockLevel<Matrix>();
763 Dune::SeqJac<Matrix, Vector, Vector, precondBlockLevel> jac(A, 1, 1.0);
772 return "Explicit diagonal matrix solver";
785class SuperLUBackend :
public LinearSolver
790 template<
class Matrix,
class Vector>
791 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
793 static_assert(isBCRSMatrix<Matrix>::value,
"SuperLU only works with BCRS matrices!");
794 using BlockType =
typename Matrix::block_type;
795 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
796 constexpr auto blockSize = BlockType::rows;
798 Dune::SuperLU<Matrix> solver(A, this->verbosity() > 0);
801 solver.apply(x, bTmp, result_);
804 for (
int i = 0; i < size; i++)
810 if (isnan(x[i][j]) || isinf(x[i][j]))
812 result_.converged =
false;
818 return result_.converged;
821 std::string name()
const
823 return "SuperLU solver";
826 const Dune::InverseOperatorResult& result()
const
832 Dune::InverseOperatorResult result_;
859class UMFPackBackend :
public LinearSolver
863 UMFPackBackend(
const std::string& paramGroup =
"")
864 : LinearSolver(paramGroup)
866 ordering_ = getParamFromGroup<int>(this->paramGroup(),
"LinearSolver.UMFPackOrdering", 1);
870 void setOrdering(
int i)
875 {
return ordering_; }
877 template<
class Matrix,
class Vector>
878 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
880 static_assert(isBCRSMatrix<Matrix>::value,
"UMFPack only works with BCRS matrices!");
881 using BlockType =
typename Matrix::block_type;
882 static_assert(BlockType::rows == BlockType::cols,
"Matrix block must be quadratic!");
883 constexpr auto blockSize = BlockType::rows;
885 Dune::UMFPack<Matrix> solver;
886 solver.setVerbosity(this->verbosity() > 0);
887 solver.setOption(UMFPACK_ORDERING, ordering_);
891 solver.apply(x, bTmp, result_);
894 for (
int i = 0; i < size; i++)
900 if (isnan(x[i][j]) || isinf(x[i][j]))
902 result_.converged =
false;
908 return result_.converged;
911 std::string name()
const
913 return "UMFPack solver";
916 const Dune::InverseOperatorResult& result()
const
922 Dune::InverseOperatorResult result_;
937template <
class LinearSolverTraits>
943 template<
class Matrix,
class Vector>
944 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
947 using Solver = Dune::BiCGSTABSolver<Vector>;
949 return IterativePreconditionedSolverImpl::template solveWithParamTree<Preconditioner, Solver>(A, x, b, solverParams);
954 return "Uzawa preconditioned BiCGSTAB solver";
962template<
class M,
class X,
class Y,
int blockLevel = 2>
965 template<std::
size_t i>
966 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
968 template<std::
size_t i>
969 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
971 template<std::
size_t i>
972 using BlockILU = Dune::SeqILU<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>, blockLevel-1>;
974 using ILUTuple =
typename makeFromIndexedType<std::tuple, BlockILU, std::make_index_sequence<M::N()> >::type;
995 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
998 void pre (X& v, Y& d)
final {}
1002 using namespace Dune::Hybrid;
1003 forEach(integralRange(Dune::Hybrid::size(ilu_)), [&](
const auto i)
1005 std::get<decltype(i)::value>(ilu_).apply(v[i], d[i]);
1014 return Dune::SolverCategory::sequential;
1018 template<std::size_t... Is>
1020 : ilu_(std::make_tuple(BlockILU<Is>(m[
Dune::index_constant<Is>{}][Dune::index_constant<Is>{}], w)...))
1040 template<
class Matrix,
class Vector>
1041 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1044 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1045 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1048 solver.apply(x, bTmp, result_);
1050 return result_.converged;
1053 const Dune::InverseOperatorResult&
result()
const
1059 {
return "block-diagonal ILU0-preconditioned BiCGSTAB solver"; }
1062 Dune::InverseOperatorResult result_;
1078 template<
int precondBlockLevel = 2,
class Matrix,
class Vector>
1079 bool solve(
const Matrix& M, Vector& x,
const Vector& b)
1082 Dune::MatrixAdapter<Matrix, Vector, Vector> op(M);
1083 static const int restartGMRes = getParamFromGroup<int>(this->
paramGroup(),
"LinearSolver.GMResRestart");
1084 Dune::RestartedGMResSolver<Vector> solver(op, preconditioner, this->
residReduction(), restartGMRes,
1087 solver.apply(x, bTmp, result_);
1089 return result_.converged;
1092 const Dune::InverseOperatorResult&
result()
const
1098 {
return "block-diagonal ILU0-preconditioned restarted GMRes solver"; }
1101 Dune::InverseOperatorResult result_;
1108template<
class M,
class X,
class Y,
int blockLevel = 2>
1111 template<std::
size_t i>
1112 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1114 template<std::
size_t i>
1115 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1117 template<std::
size_t i>
1118 using Smoother = Dune::SeqSSOR<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1120 template<std::
size_t i>
1121 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<i>, VecBlockType<i>, VecBlockType<i>>;
1123 template<std::
size_t i>
1124 using ScalarProduct = Dune::SeqScalarProduct<VecBlockType<i>>;
1126 template<std::
size_t i>
1127 using BlockAMG = Dune::Amg::AMG<LinearOperator<i>, VecBlockType<i>, Smoother<i>>;
1129 using AMGTuple =
typename makeFromIndexedType<std::tuple, BlockAMG, std::make_index_sequence<M::N()> >::type;
1148 template<
class LOP,
class Criterion,
class SmootherArgs>
1152 static_assert(blockLevel >= 2,
"Only makes sense for MultiTypeBlockMatrix!");
1173 using namespace Dune::Hybrid;
1174 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1176 std::get<decltype(i)::value>(amg_).pre(v[i], d[i]);
1193 using namespace Dune::Hybrid;
1194 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1196 std::get<decltype(i)::value>(amg_).apply(v[i], d[i]);
1211 using namespace Dune::Hybrid;
1212 forEach(integralRange(Dune::Hybrid::size(amg_)), [&](
const auto i)
1214 std::get<decltype(i)::value>(amg_).post(v[i]);
1221 return Dune::SolverCategory::sequential;
1225 template<
class LOP,
class Criterion,
class SmootherArgs, std::size_t... Is>
1227 : amg_(std::make_tuple(BlockAMG<Is>(*std::get<Is>(lop), *std::get<Is>(c), *std::get<Is>(sa))...))
1242 template<
class M, std::
size_t i>
1243 using DiagBlockType = std::decay_t<decltype(std::declval<M>()[Dune::index_constant<i>{}][Dune::index_constant<i>{}])>;
1245 template<
class X, std::
size_t i>
1246 using VecBlockType = std::decay_t<decltype(std::declval<X>()[Dune::index_constant<i>{}])>;
1248 template<
class M,
class X, std::
size_t i>
1249 using Smoother = Dune::SeqSSOR<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1251 template<
class M,
class X, std::
size_t i>
1252 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother<M, X, i>>::Arguments;
1254 template<
class M, std::
size_t i>
1255 using Criterion = Dune::Amg::CoarsenCriterion<Dune::Amg::SymmetricCriterion<DiagBlockType<M, i>, Dune::Amg::FirstDiagonal>>;
1257 template<
class M,
class X, std::
size_t i>
1258 using LinearOperator = Dune::MatrixAdapter<DiagBlockType<M, i>, VecBlockType<X, i>, VecBlockType<X, i>>;
1264 template<
class Matrix,
class Vector>
1265 bool solve(
const Matrix& m, Vector& x,
const Vector& b)
1269 Dune::Amg::Parameters params(15, 2000, 1.2, 1.6, Dune::Amg::atOnceAccu);
1270 params.setDefaultValuesIsotropic(3);
1271 params.setDebugLevel(this->
verbosity());
1273 auto criterion = makeCriterion_<Criterion, Matrix>(params, std::make_index_sequence<Matrix::N()>{});
1274 auto smootherArgs = makeSmootherArgs_<SmootherArgs, Matrix, Vector>(std::make_index_sequence<Matrix::N()>{});
1276 using namespace Dune::Hybrid;
1277 forEach(std::make_index_sequence<Matrix::N()>{}, [&](
const auto i)
1279 auto& args = std::get<decltype(i)::value>(smootherArgs);
1280 args->iterations = 1;
1281 args->relaxationFactor = 1;
1284 auto linearOperator = makeLinearOperator_<LinearOperator, Matrix, Vector>(m, std::make_index_sequence<Matrix::N()>{});
1288 Dune::MatrixAdapter<Matrix, Vector, Vector> op(m);
1289 Dune::BiCGSTABSolver<Vector> solver(op, preconditioner, this->
residReduction(),
1292 solver.apply(x, bTmp, result_);
1294 return result_.converged;
1297 const Dune::InverseOperatorResult&
result()
const
1303 {
return "block-diagonal AMG-preconditioned BiCGSTAB solver"; }
1306 template<
template<
class M, std::
size_t i>
class Criterion,
class Matrix,
class Params, std::size_t... Is>
1307 auto makeCriterion_(
const Params& p, std::index_sequence<Is...>)
1309 return std::make_tuple(std::make_shared<Criterion<Matrix, Is>>(p)...);
1312 template<
template<
class M,
class X, std::
size_t i>
class SmootherArgs,
class Matrix,
class Vector, std::size_t... Is>
1313 auto makeSmootherArgs_(std::index_sequence<Is...>)
1315 return std::make_tuple(std::make_shared<SmootherArgs<Matrix, Vector, Is>>()...);
1318 template<
template<
class M,
class X, std::
size_t i>
class LinearOperator,
class Matrix,
class Vector, std::size_t... Is>
1319 auto makeLinearOperator_(
const Matrix& m, std::index_sequence<Is...>)
1321 return std::make_tuple(std::make_shared<LinearOperator<Matrix, Vector, Is>>(m[Dune::index_constant<Is>{}][Dune::index_constant<Is>{}])...);
1324 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:182
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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: deprecated.hh:149
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:49
A preconditioner based on the Uzawa algorithm for saddle-point problems of the form .
Definition: preconditioners.hh:79
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
static bool solveWithParamTree(const Matrix &A, Vector &x, const Vector &b, const Dune::ParameterTree ¶ms)
Definition: seqsolverbackend.hh:156
Sequential ILU(n)-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:205
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:210
std::string name() const
Definition: seqsolverbackend.hh:219
Sequential SOR-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:243
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:248
std::string name() const
Definition: seqsolverbackend.hh:257
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:281
std::string name() const
Definition: seqsolverbackend.hh:295
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:286
Sequential GS-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:319
std::string name() const
Definition: seqsolverbackend.hh:333
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:324
Sequential Jacobi-preconditioned BiCSTAB solver.
Definition: seqsolverbackend.hh:356
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:361
std::string name() const
Definition: seqsolverbackend.hh:370
Sequential ILU(n)-preconditioned CG solver.
Definition: seqsolverbackend.hh:393
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:398
std::string name() const
Definition: seqsolverbackend.hh:407
Sequential SOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:430
std::string name() const
Definition: seqsolverbackend.hh:444
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:435
Sequential SSOR-preconditioned CG solver.
Definition: seqsolverbackend.hh:467
std::string name() const
Definition: seqsolverbackend.hh:481
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:472
Sequential GS-preconditioned CG solver.
Definition: seqsolverbackend.hh:504
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:509
std::string name() const
Definition: seqsolverbackend.hh:518
Sequential Jacobi-preconditioned CG solver.
Definition: seqsolverbackend.hh:540
std::string name() const
Definition: seqsolverbackend.hh:554
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:545
Sequential SSOR-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:578
std::string name() const
Definition: seqsolverbackend.hh:592
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:583
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:615
std::string name() const
Definition: seqsolverbackend.hh:629
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:620
Sequential ILU(0)-preconditioned CG solver.
Definition: seqsolverbackend.hh:651
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:656
std::string name() const
Definition: seqsolverbackend.hh:665
Sequential ILU0-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:688
std::string name() const
Definition: seqsolverbackend.hh:702
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:693
Sequential ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:726
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:731
std::string name() const
Definition: seqsolverbackend.hh:740
Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes)
Definition: seqsolverbackend.hh:754
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:759
std::string name() const
Definition: seqsolverbackend.hh:770
A Uzawa preconditioned BiCGSTAB solver for saddle-point problems.
Definition: seqsolverbackend.hh:939
bool solve(const Matrix &A, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:944
std::string name() const
Definition: seqsolverbackend.hh:952
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:964
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:984
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:978
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:982
void pre(X &v, Y &d) final
Definition: seqsolverbackend.hh:998
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:980
void post(X &) final
Definition: seqsolverbackend.hh:1009
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1012
void apply(X &v, const Y &d) final
Definition: seqsolverbackend.hh:1000
BlockDiagILU0Preconditioner(const M &m, double w=1.0)
Constructor.
Definition: seqsolverbackend.hh:992
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1035
std::string name() const
Definition: seqsolverbackend.hh:1058
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1053
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1041
A simple ilu0 block diagonal preconditioned RestartedGMResSolver.
Definition: seqsolverbackend.hh:1073
std::string name() const
Definition: seqsolverbackend.hh:1097
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1092
bool solve(const Matrix &M, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1079
A simple ilu0 block diagonal preconditioner.
Definition: seqsolverbackend.hh:1110
void post(X &v) final
Clean up.
Definition: seqsolverbackend.hh:1209
typename std::decay_t< M > matrix_type
The matrix type the preconditioner is for.
Definition: seqsolverbackend.hh:1133
void apply(X &v, const Y &d) final
Apply one step of the preconditioner to the system A(v)=d.
Definition: seqsolverbackend.hh:1191
Dune::SolverCategory::Category category() const final
Category of the preconditioner (see SolverCategory::Category)
Definition: seqsolverbackend.hh:1219
Y range_type
The range type of the preconditioner.
Definition: seqsolverbackend.hh:1137
void pre(X &v, Y &d) final
Prepare the preconditioner.
Definition: seqsolverbackend.hh:1171
typename X::field_type field_type
The field type of the preconditioner.
Definition: seqsolverbackend.hh:1139
BlockDiagAMGPreconditioner(const LOP &lop, const Criterion &c, const SmootherArgs &sa)
Constructor.
Definition: seqsolverbackend.hh:1149
X domain_type
The domain type of the preconditioner.
Definition: seqsolverbackend.hh:1135
A simple ilu0 block diagonal preconditioned BiCGSTABSolver.
Definition: seqsolverbackend.hh:1241
const Dune::InverseOperatorResult & result() const
Definition: seqsolverbackend.hh:1297
std::string name() const
Definition: seqsolverbackend.hh:1302
bool solve(const Matrix &m, Vector &x, const Vector &b)
Definition: seqsolverbackend.hh:1265
Base class for linear solvers.
Definition: solver.hh:37
Scalar 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="")
Construct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:83