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