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,
46 static const std::string solverMethodTable(
int id)
48 return Solvers::to_string(
id );
54 LIST_OF_INT_FORWARDED(Preconditioners,
66 static const std::string preconditionMethodTable(
int id)
68 return Preconditioners::to_string(
id );
71 SolverParameter (
const ParameterReader ¶meter = Parameter::container() )
72 : keyPrefix_(
"fem.solver." ), parameter_( parameter )
75 explicit SolverParameter (
const std::string keyPrefix,
const ParameterReader ¶meter = Parameter::container() )
76 : keyPrefix_( keyPrefix ), parameter_( parameter )
79 const std::string& keyPrefix()
const {
return keyPrefix_; }
81 const ParameterReader& parameter()
const {
return parameter_; }
91 virtual bool verbose()
const
95 verbose_ = parameter_.getValue<
bool >( keyPrefix_ +
"verbose", false ) ? 1 : 0;
97 return bool(verbose_);
100 virtual void setVerbose(
const bool verb )
102 verbose_ = verb ? 1 : 0;
107 int defaultErrorMeasure = 0;
108 virtual void setDefaultErrorMeasure(
int def) { defaultErrorMeasure = def; }
109 virtual int errorMeasure()
const
111 const std::string errorTypeTable[] =
112 {
"absolute",
"relative",
"residualreduction" };
113 const int errorType = parameter_.getEnum( keyPrefix_ +
"errormeasure", errorTypeTable,
114 defaultErrorMeasure );
118 virtual double tolerance ( )
const
122 double defaultTol = 1e-8;
123 if(parameter_.exists(keyPrefix_ +
"tolerance"))
124 tolerance_ = parameter_.getValue<
double >( keyPrefix_ +
"tolerance", defaultTol );
127 if( parameter_.exists(keyPrefix_ +
"absolutetol") ||
128 parameter_.exists(keyPrefix_ +
"reductiontol") )
131 int measure = errorMeasure();
133 tolerance_ = absoluteTol__();
135 tolerance_ = reductionTol__();
137 else tolerance_ = defaultTol;
143 virtual void setTolerance (
const double eps )
145 assert( eps >= 0.0 );
149 virtual int maxIterations ()
const
151 if( maxIterations_ < 0 )
153 if(parameter_.exists(keyPrefix_ +
"maxlineariterations"))
155 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"maxlineariterations is deprecated. Please use " + keyPrefix_ +
"maxiterations instead." << std::endl;
156 maxIterations_ = parameter_.getValue<
double >(keyPrefix_ +
"maxlineariterations");
161 return maxIterations_;
164 virtual void setMaxIterations (
const int maxIter )
166 assert( maxIter >= 0 );
167 maxIterations_ = maxIter;
170 virtual int solverMethod(
171 const std::vector<int> standardMethods,
172 const std::vector<std::string> &additionalMethods = {},
173 int defaultMethod = 0
176 std::vector<std::string> methodTable(standardMethods.size()+additionalMethods.size());
177 for (std::size_t i=0;i<standardMethods.size();++i)
178 methodTable[i] = solverMethodTable(standardMethods[i]);
179 for (std::size_t i=0;i<additionalMethods.size();++i)
180 methodTable[standardMethods.size()+i] = additionalMethods[i];
182 if( parameter_.exists( keyPrefix_ +
"method" ) ||
183 !parameter_.exists(
"method" ) )
184 method = parameter_.getEnum( keyPrefix_ +
"method", methodTable, defaultMethod );
187 method = parameter_.getEnum(
"krylovmethod", methodTable, defaultMethod );
188 std::cout <<
"WARNING: using old parameter name 'krylovmethod' "
189 <<
"please switch to '" << keyPrefix_ <<
"method'\n";
191 if (method < standardMethods.size())
192 return standardMethods[method];
194 return -(method-standardMethods.size());
197 virtual int gmresRestart()
const
199 int defaultRestart = 20;
200 return parameter_.getValue<
int >( keyPrefix_ +
"gmres.restart", defaultRestart );
203 virtual int preconditionMethod(
204 const std::vector<int> standardMethods,
205 const std::vector<std::string> &additionalMethods = {},
206 int defaultMethod = 0
209 std::vector<std::string> methodTable(standardMethods.size()+additionalMethods.size());
210 for (std::size_t i=0;i<standardMethods.size();++i)
211 methodTable[i] = preconditionMethodTable(standardMethods[i]);
212 for (std::size_t i=0;i<additionalMethods.size();++i)
213 methodTable[standardMethods.size()+i] = additionalMethods[i];
214 std::size_t method = parameter_.getEnum( keyPrefix_ +
"preconditioning.method", methodTable, defaultMethod );
215 if (method < standardMethods.size())
216 return standardMethods[method];
218 return -(method-standardMethods.size());
221 virtual double relaxation ()
const
223 return parameter_.getValue<
double >( keyPrefix_ +
"preconditioning.relaxation", 1.1 );
226 virtual int preconditionerIteration ()
const
228 const int defaultIter = 1;
229 if(parameter_.exists(keyPrefix_ +
"preconditioning.iteration"))
231 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"preconditioning.iteration is deprecated. Please use " + keyPrefix_ +
"preconditioning.iterations instead." << std::endl;
232 return parameter_.getValue<
int >( keyPrefix_ +
"preconditioning.iteration", defaultIter );
234 return parameter_.getValue<
int >( keyPrefix_ +
"preconditioning.iterations", defaultIter );
237 virtual int preconditionerLevel ()
const
239 return parameter_.getValue<
int >( keyPrefix_ +
"preconditioning.level", 0 );
242 virtual bool threading ()
const
244 return parameter_.getValue<
bool >( keyPrefix_ +
"threading", true );
247 virtual bool knollTrick()
const
249 return parameter_.getValue<
bool >( keyPrefix_ +
"knolltrick", false );
253 virtual double absoluteTol__ ( )
const
255 if( absoluteTol_ < 0 )
257 if(parameter_.exists(keyPrefix_ +
"linabstol"))
259 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"linabstol is deprecated. Please use " + keyPrefix_ +
"absolutetol instead." << std::endl;
260 absoluteTol_ = parameter_.getValue<
double >(keyPrefix_ +
"linabstol");
264 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"absolutetol is deprecated. Please use " + keyPrefix_ +
"tolerance instead." << std::endl;
265 absoluteTol_ = parameter_.getValue<
double >(keyPrefix_ +
"absolutetol", 1e-8 );
271 virtual double reductionTol__ ( )
const
273 if( reductionTol_ < 0 )
275 if(parameter_.exists(keyPrefix_ +
"linreduction"))
277 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"linreduction is deprecated. Please use " + keyPrefix_ +
"reductiontol instead." << std::endl;
278 reductionTol_ = parameter_.getValue<
double >(keyPrefix_ +
"linreduction");
282 std::cout <<
"WARNING: Parameter " + keyPrefix_ +
"reductiontol is deprecated. Please use " + keyPrefix_ +
"tolerance instead." << std::endl;
283 reductionTol_ = parameter_.getValue<
double >( keyPrefix_ +
"reductiontol", 1e-8 );
286 return reductionTol_;
289 mutable int verbose_ = -1;
290 mutable int maxIterations_ = -1;
291 mutable double absoluteTol_ = -1.;
292 mutable double reductionTol_ = -1.;
293 mutable double tolerance_ = -1.;
@ absolute
|a-b| <= epsilon
Definition: float_cmp.hh:113
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:489
Dune namespace
Definition: alignedallocator.hh:13