1#ifndef DUNE_FEM_SOLVERPARAMETER_HH
2#define DUNE_FEM_SOLVERPARAMETER_HH
4#include <dune/fem/io/parameter.hh>
5#include <dune/fem/common/staticlistofint.hh>
12 namespace LinearSolver
14 struct ToleranceCriteria {
16 static const int relative = 1;
17 static const int residualReduction = 2;
21 struct SolverParameter
23 :
public LocalParameter< SolverParameter, SolverParameter >
28 const std::string keyPrefix_;
30 ParameterReader parameter_;
34 LIST_OF_INT_FORWARDED(Solvers,
45 static const std::string solverMethodTable(
int id)
47 return Solvers::to_string(
id );
53 LIST_OF_INT_FORWARDED(Preconditioners,
65 static const std::string preconditionMethodTable(
int id)
67 return Preconditioners::to_string(
id );
70 SolverParameter (
const ParameterReader ¶meter = Parameter::container() )
71 : keyPrefix_(
"fem.solver." ), parameter_( parameter )
74 explicit SolverParameter (
const std::string keyPrefix,
const ParameterReader ¶meter = Parameter::container() )
75 : keyPrefix_( keyPrefix ), parameter_( parameter )
78 const std::string& keyPrefix()
const {
return keyPrefix_; }
80 const ParameterReader& parameter()
const {
return parameter_; }
90 virtual bool verbose()
const
94 verbose_ = parameter_.getValue<
bool >( keyPrefix_ +
"verbose", false ) ? 1 : 0;
96 return bool(verbose_);
99 virtual void setVerbose(
const bool verb )
101 verbose_ = verb ? 1 : 0;
104 virtual int errorMeasure()
const
106 const std::string errorTypeTable[] =
107 {
"absolute",
"relative",
"residualreduction" };
108 const int errorType = parameter_.getEnum( keyPrefix_ +
"errormeasure", errorTypeTable, 0 );
112 virtual double tolerance ( )
const
116 double defaultTol = 1e-8;
117 if(parameter_.exists(keyPrefix_ +
"tolerance"))
118 tolerance_ = parameter_.getValue<
double >( keyPrefix_ +
"tolerance", defaultTol );
121 if( parameter_.exists(keyPrefix_ +
"absolutetol") ||
122 parameter_.exists(keyPrefix_ +
"reductiontol") )
125 int measure = errorMeasure();
127 tolerance_ = absoluteTol__();
129 tolerance_ = reductionTol__();
131 else tolerance_ = defaultTol;
137 virtual void setTolerance (
const double eps )
139 assert( eps >= 0.0 );
143 virtual int maxIterations ()
const
145 if( maxIterations_ < 0 )
147 if(parameter_.exists(keyPrefix_ +
"maxlineariterations"))
149 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"maxlineariterations is deprecated. Please use " + keyPrefix_ +
"maxiterations instead." << std::endl;
150 maxIterations_ = parameter_.getValue<
double >(keyPrefix_ +
"maxlineariterations");
155 return maxIterations_;
158 virtual void setMaxIterations (
const int maxIter )
160 assert( maxIter >= 0 );
161 maxIterations_ = maxIter;
164 virtual int solverMethod(
165 const std::vector<int> standardMethods,
166 const std::vector<std::string> &additionalMethods = {},
167 int defaultMethod = 0
170 std::vector<std::string> methodTable(standardMethods.size()+additionalMethods.size());
171 for (std::size_t i=0;i<standardMethods.size();++i)
172 methodTable[i] = solverMethodTable(standardMethods[i]);
173 for (std::size_t i=0;i<additionalMethods.size();++i)
174 methodTable[standardMethods.size()+i] = additionalMethods[i];
176 if( parameter_.exists( keyPrefix_ +
"method" ) ||
177 !parameter_.exists(
"method" ) )
178 method = parameter_.getEnum( keyPrefix_ +
"method", methodTable, defaultMethod );
181 method = parameter_.getEnum(
"krylovmethod", methodTable, defaultMethod );
182 std::cout <<
"WARNING: using old parameter name 'krylovmethod' "
183 <<
"please switch to '" << keyPrefix_ <<
"method'\n";
185 if (method < standardMethods.size())
186 return standardMethods[method];
188 return -(method-standardMethods.size());
191 virtual int gmresRestart()
const
193 int defaultRestart = 20;
194 return parameter_.getValue<
int >( keyPrefix_ +
"gmres.restart", defaultRestart );
197 virtual int preconditionMethod(
198 const std::vector<int> standardMethods,
199 const std::vector<std::string> &additionalMethods = {},
200 int defaultMethod = 0
203 std::vector<std::string> methodTable(standardMethods.size()+additionalMethods.size());
204 for (std::size_t i=0;i<standardMethods.size();++i)
205 methodTable[i] = preconditionMethodTable(standardMethods[i]);
206 for (std::size_t i=0;i<additionalMethods.size();++i)
207 methodTable[standardMethods.size()+i] = additionalMethods[i];
208 std::size_t method = parameter_.getEnum( keyPrefix_ +
"preconditioning.method", methodTable, defaultMethod );
209 if (method < standardMethods.size())
210 return standardMethods[method];
212 return -(method-standardMethods.size());
215 virtual double relaxation ()
const
217 return parameter_.getValue<
double >( keyPrefix_ +
"preconditioning.relaxation", 1.1 );
220 virtual int preconditionerIteration ()
const
222 return parameter_.getValue<
int >( keyPrefix_ +
"preconditioning.iterations", 1 );
225 virtual int preconditionerLevel ()
const
227 return parameter_.getValue<
int >( keyPrefix_ +
"preconditioning.level", 0 );
230 virtual bool threading ()
const
232 return parameter_.getValue<
bool >( keyPrefix_ +
"threading", true );
235 virtual bool knollTrick()
const
237 return parameter_.getValue<
bool >( keyPrefix_ +
"knolltrick", false );
241 virtual double absoluteTol__ ( )
const
243 if( absoluteTol_ < 0 )
245 if(parameter_.exists(keyPrefix_ +
"linabstol"))
247 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"linabstol is deprecated. Please use " + keyPrefix_ +
"absolutetol instead." << std::endl;
248 absoluteTol_ = parameter_.getValue<
double >(keyPrefix_ +
"linabstol");
252 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"absolutetol is deprecated. Please use " + keyPrefix_ +
"tolerance instead." << std::endl;
253 absoluteTol_ = parameter_.getValue<
double >(keyPrefix_ +
"absolutetol", 1e-8 );
259 virtual double reductionTol__ ( )
const
261 if( reductionTol_ < 0 )
263 if(parameter_.exists(keyPrefix_ +
"linreduction"))
265 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"linreduction is deprecated. Please use " + keyPrefix_ +
"reductiontol instead." << std::endl;
266 reductionTol_ = parameter_.getValue<
double >(keyPrefix_ +
"linreduction");
270 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"reductiontol is deprecated. Please use " + keyPrefix_ +
"tolerance instead." << std::endl;
271 reductionTol_ = parameter_.getValue<
double >( keyPrefix_ +
"reductiontol", 1e-8 );
274 return reductionTol_;
277 mutable int verbose_ = -1;
278 mutable int maxIterations_ = -1;
279 mutable double absoluteTol_ = -1.;
280 mutable double reductionTol_ = -1.;
281 mutable double tolerance_ = -1.;
@ absolute
|a-b| <= epsilon
Definition: float_cmp.hh:110
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:471
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
Dune namespace.
Definition: alignedallocator.hh:13