24#ifndef DUMUX_LINEAR_SOLVER_HH
25#define DUMUX_LINEAR_SOLVER_HH
27#include <dune/common/exceptions.hh>
56 verbosity_ = getParamFromGroup<int>(
paramGroup,
"LinearSolver.Verbosity", 0);
57 maxIter_ = getParamFromGroup<int>(
paramGroup,
"LinearSolver.MaxIterations", 250);
58 residReduction_ = getParamFromGroup<double>(
paramGroup,
"LinearSolver.ResidualReduction", 1e-13);
65 std::cerr <<
"Deprecation warning: parameter LinearSolver.PreconditionerRelaxation is deprecated and will be removed after 3.2. "
66 <<
"Use LinearSolver.Preconditioner.Relaxation instead (Preconditioner subgroup)." << std::endl;
67 relaxation_ = getParamFromGroup<double>(
paramGroup,
"LinearSolver.PreconditionerRelaxation");
70 relaxation_ = getParamFromGroup<double>(
paramGroup,
"LinearSolver.Preconditioner.Relaxation", 1);
74 std::cerr <<
"Deprecation warning: parameter LinearSolver.PreconditionerIterations is deprecated and will be removed after 3.2. "
75 <<
"Use LinearSolver.Preconditioner.Iterations instead (Preconditioner subgroup)." << std::endl;
76 precondIter_ = getParamFromGroup<int>(
paramGroup,
"LinearSolver.PreconditionerIterations");
79 precondIter_ = getParamFromGroup<int>(
paramGroup,
"LinearSolver.Preconditioner.Iterations", 1);
83 std::cerr <<
"Deprecation warning: parameter LinearSolver.PreconditionerVerbosity is deprecated and will be removed after 3.2. "
84 <<
"Use LinearSolver.Preconditioner.Verbosity instead (Preconditioner subgroup)." << std::endl;
85 precondVerbosity_ = getParamFromGroup<int>(
paramGroup,
"LinearSolver.PreconditionerVerbosity");
88 precondVerbosity_ = getParamFromGroup<int>(
paramGroup,
"LinearSolver.Preconditioner.Verbosity", 0);
100 template<
class Matrix,
class Vector>
101 bool solve(
const Matrix& A, Vector& x,
const Vector& b)
103 DUNE_THROW(Dune::NotImplemented,
"Linear solver doesn't implement a solve method!");
108 {
return "unknown solver"; }
112 {
return paramGroup_; }
116 {
return verbosity_; }
132 {
return residReduction_; }
136 { residReduction_ = r; }
140 {
return relaxation_; }
148 {
return precondIter_; }
152 { precondIter_ = i; }
156 {
return precondVerbosity_; }
160 { precondVerbosity_ = verbosityLevel; }
165 double residReduction_;
168 int precondVerbosity_;
169 const std::string paramGroup_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:391
Base class for linear solvers.
Definition: solver.hh:37
bool solve(const Matrix &A, Vector &x, const Vector &b)
Solve the linear system Ax = b.
Definition: solver.hh:101
int precondVerbosity() const
the preconditioner verbosity
Definition: solver.hh:155
void setResidualReduction(double r)
set the linear solver residual reduction
Definition: solver.hh:135
void setMaxIter(int i)
set the maximum number of linear solver iterations
Definition: solver.hh:127
int precondIter() const
the number of preconditioner iterations
Definition: solver.hh:147
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
void setPrecondIter(int i)
set the number of preconditioner iterations
Definition: solver.hh:151
void setVerbosity(int v)
set the verbosity level
Definition: solver.hh:119
void setPrecondVerbosity(int verbosityLevel)
set the preconditioner verbosity
Definition: solver.hh:159
const std::string & paramGroup() const
the parameter group for getting parameter from the parameter tree
Definition: solver.hh:111
void setRelaxation(double r)
set the linear solver relaxation factor
Definition: solver.hh:143
LinearSolver(const std::string ¶mGroup="")
Contruct the solver.
Definition: solver.hh:53
int verbosity() const
the verbosity level
Definition: solver.hh:115
std::string name() const
the name of the linear solver
Definition: solver.hh:107
double Scalar
Definition: solver.hh:41
double relaxation() const
the linear solver relaxation factor
Definition: solver.hh:139