DUNE-FEM (unstable)

parameter.hh
1#ifndef DUNE_FEM_SOLVERPARAMETER_HH
2#define DUNE_FEM_SOLVERPARAMETER_HH
3
4#include <dune/fem/io/parameter.hh>
5#include <dune/fem/common/staticlistofint.hh>
6
7namespace Dune
8{
9
10 namespace Fem
11 {
12 namespace LinearSolver
13 {
14 struct ToleranceCriteria {
15 static const int absolute = 0;
16 static const int relative = 1;
17 static const int residualReduction = 2;
18 };
19 }
20
21 struct SolverParameter
22#ifndef DOXYGEN
23 : public LocalParameter< SolverParameter, SolverParameter >
24#endif
25 {
26 protected:
27 // key prefix, default is fem.solver (can be overloaded by user)
28 const std::string keyPrefix_;
29
30 ParameterReader parameter_;
31
32 public:
33 // identifier for Fem, ISTL and Petsc solvers
34 LIST_OF_INT_FORWARDED(Solvers,
35 cg = 1, // CG
36 bicgstab = 2, // BiCGStab
37 gmres = 3, // GMRES
38 minres = 4, // MinRes
39 gradient = 5, // GradientSolver
40 loop = 6, // LoopSolver
41 superlu = 7, // SuperLUSolver
42 bicg = 8, // BiCG
43 preonly = 9, // only preconder
44 fgmres = 10);// Flexible GMRes
45
46 static const std::string solverMethodTable(int id)
47 {
48 return Solvers::to_string( id );
49 }
50
51 // identifier for Fem, ISTL and Petsc preconditioner
52 // Note: underscores in variable names will be replaced with dash in the string of the name
53 // e.g. amg_jacobi --> amg-jacobi as string
54 LIST_OF_INT_FORWARDED(Preconditioners,
55 none = 1 , // no preconditioner
56 ssor = 2 , // SSOR preconditioner
57 sor = 3 , // SOR preconditioner
58 ilu = 4 , // ILU preconditioner
59 gauss_seidel = 5 , // Gauss-Seidel preconditioner
60 jacobi = 6 , // Jacobi preconditioner
61 amg_ilu = 7 , // AMG with ILU-0 smoother (deprecated)
62 amg_jacobi = 8 , // AMG with Jacobi smoother
63 ildl = 9 , // ILDL from istl
64 oas = 10, // Overlapping Additive Schwarz
65 icc = 11);// Incomplete Cholesky factorization
66 static const std::string preconditionMethodTable(int id)
67 {
68 return Preconditioners::to_string( id );
69 }
70
71 SolverParameter ( const ParameterReader &parameter = Parameter::container() )
72 : keyPrefix_( "fem.solver." ), parameter_( parameter )
73 {}
74
75 explicit SolverParameter ( const std::string keyPrefix, const ParameterReader &parameter = Parameter::container() )
76 : keyPrefix_( keyPrefix ), parameter_( parameter )
77 {}
78
79 const std::string& keyPrefix() const { return keyPrefix_; }
80
81 const ParameterReader& parameter() const { return parameter_; }
82
83 virtual void reset()
84 {
85 verbose_ = -1;
86 absoluteTol_ = -1;
87 reductionTol_ = -1;
88 maxIterations_ = -1;
89 }
90
91 virtual bool verbose() const
92 {
93 if( verbose_ < 0 )
94 {
95 verbose_ = parameter_.getValue< bool >( keyPrefix_ + "verbose", false ) ? 1 : 0;
96 }
97 return bool(verbose_);
98 }
99
100 virtual void setVerbose( const bool verb )
101 {
102 verbose_ = verb ? 1 : 0;
103 }
104
105 // default value is 'absolute' except when using Eisenstat-Walker
106 // in that case the 'default' should be 'residualreduction'
107 int defaultErrorMeasure = 0;
108 virtual void setDefaultErrorMeasure(int def) { defaultErrorMeasure = def; }
109 virtual int errorMeasure() const
110 {
111 const std::string errorTypeTable[] =
112 { "absolute", "relative", "residualreduction" };
113 const int errorType = parameter_.getEnum( keyPrefix_ + "errormeasure", errorTypeTable,
114 defaultErrorMeasure );
115 return errorType ;
116 }
117
118 virtual double tolerance ( ) const
119 {
120 if( tolerance_ < 0 )
121 {
122 double defaultTol = 1e-8;
123 if(parameter_.exists(keyPrefix_ + "tolerance"))
124 tolerance_ = parameter_.getValue< double >( keyPrefix_ + "tolerance", defaultTol );
125 else
126 {
127 if( parameter_.exists(keyPrefix_ + "absolutetol") ||
128 parameter_.exists(keyPrefix_ + "reductiontol") )
129 // new parameter not used but old parameters exist
130 {
131 int measure = errorMeasure();
132 if (measure == 0)
133 tolerance_ = absoluteTol__();
134 else
135 tolerance_ = reductionTol__();
136 }
137 else tolerance_ = defaultTol; // no parameter set
138 }
139 }
140 return tolerance_;
141 }
142
143 virtual void setTolerance ( const double eps )
144 {
145 assert( eps >= 0.0 );
146 tolerance_ = eps;
147 }
148
149 virtual int maxIterations () const
150 {
151 if( maxIterations_ < 0 )
152 {
153 if(parameter_.exists(keyPrefix_ + "maxlineariterations"))
154 {
155 std::cout << "WARNING: Parameter " + keyPrefix_ +"maxlineariterations is deprecated. Please use " + keyPrefix_ + "maxiterations instead." << std::endl;
156 maxIterations_ = parameter_.getValue< double >(keyPrefix_ + "maxlineariterations");
157 }
158 else
159 maxIterations_ = parameter_.getValue< int >( keyPrefix_ + "maxiterations", std::numeric_limits< int >::max() );
160 }
161 return maxIterations_;
162 }
163
164 virtual void setMaxIterations ( const int maxIter )
165 {
166 assert( maxIter >= 0 );
167 maxIterations_ = maxIter;
168 }
169
170 virtual int solverMethod(
171 const std::vector<int> standardMethods,
172 const std::vector<std::string> &additionalMethods = {},
173 int defaultMethod = 0 // this is the first method passed in
174 ) const
175 {
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];
181 std::size_t method;
182 if( parameter_.exists( keyPrefix_ + "method" ) ||
183 !parameter_.exists( "method" ) )
184 method = parameter_.getEnum( keyPrefix_ + "method", methodTable, defaultMethod );
185 else
186 {
187 method = parameter_.getEnum( "krylovmethod", methodTable, defaultMethod );
188 std::cout << "WARNING: using old parameter name 'krylovmethod' "
189 << "please switch to '" << keyPrefix_ << "method'\n";
190 }
191 if (method < standardMethods.size())
192 return standardMethods[method];
193 else
194 return -(method-standardMethods.size()); // return in [ 0,-additionalMethods.size() )
195 }
196
197 virtual int gmresRestart() const
198 {
199 int defaultRestart = 20;
200 return parameter_.getValue< int >( keyPrefix_ + "gmres.restart", defaultRestart );
201 }
202
203 virtual int preconditionMethod(
204 const std::vector<int> standardMethods,
205 const std::vector<std::string> &additionalMethods = {},
206 int defaultMethod = 0 // this is the first method passed in
207 ) const
208 {
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];
217 else
218 return -(method-standardMethods.size()); // return in [ 0,-additionalMethods.size() )
219 }
220
221 virtual double relaxation () const
222 {
223 return parameter_.getValue< double >( keyPrefix_ + "preconditioning.relaxation", 1.1 );
224 }
225
226 virtual int preconditionerIteration () const
227 {
228 const int defaultIter = 1;
229 if(parameter_.exists(keyPrefix_ + "preconditioning.iteration"))
230 {
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 );
233 }
234 return parameter_.getValue< int >( keyPrefix_ + "preconditioning.iterations", defaultIter );
235 }
236
237 virtual int preconditionerLevel () const
238 {
239 return parameter_.getValue< int >( keyPrefix_ + "preconditioning.level", 0 );
240 }
241
242 virtual bool threading () const
243 {
244 return parameter_.getValue< bool >( keyPrefix_ + "threading", true );
245 }
246
247 virtual bool knollTrick() const
248 {
249 return parameter_.getValue< bool >( keyPrefix_ + "knolltrick", false );
250 }
251
252 private:
253 virtual double absoluteTol__ ( ) const
254 {
255 if( absoluteTol_ < 0 )
256 {
257 if(parameter_.exists(keyPrefix_ + "linabstol"))
258 {
259 std::cout << "WARNING: Parameter " + keyPrefix_ + "linabstol is deprecated. Please use " + keyPrefix_ + "absolutetol instead." << std::endl;
260 absoluteTol_ = parameter_.getValue< double >(keyPrefix_ + "linabstol");
261 }
262 else
263 {
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 );
266 }
267 }
268 return absoluteTol_;
269 }
270
271 virtual double reductionTol__ ( ) const
272 {
273 if( reductionTol_ < 0 )
274 {
275 if(parameter_.exists(keyPrefix_ + "linreduction"))
276 {
277 std::cout << "WARNING: Parameter " + keyPrefix_ +"linreduction is deprecated. Please use " + keyPrefix_ + "reductiontol instead." << std::endl;
278 reductionTol_ = parameter_.getValue< double >(keyPrefix_ + "linreduction");
279 }
280 else
281 {
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 );
284 }
285 }
286 return reductionTol_;
287 }
288
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.;
294 };
295 }
296}
297
298#endif // #ifndef DUNE_FEM_SOLVERPARAMETER_HH
@ 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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (May 31, 22:31, 2026)