FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
gmres.hpp
1// FEAT3: Finite Element Analysis Toolbox, Version 3
2// Copyright (C) 2010 by Stefan Turek & the FEAT group
3// FEAT3 is released under the GNU General Public License version 3,
4// see the file 'copyright.txt' in the top level directory for details.
5
6#pragma once
7
8// includes, FEAT
9#include <kernel/solver/iterative.hpp>
10
11// includes, FEAT
12#include <vector>
13
14namespace FEAT
15{
16 namespace Solver
17 {
73 template<
74 typename Matrix_,
75 typename Filter_>
76 class GMRES :
77 public PreconditionedIterativeSolver<typename Matrix_::VectorTypeR>
78 {
79 public:
80 typedef Matrix_ MatrixType;
81 typedef Filter_ FilterType;
82 typedef typename MatrixType::VectorTypeR VectorType;
83 typedef typename MatrixType::DataType DataType;
85
87
88 protected:
90 const MatrixType& _system_matrix;
92 const FilterType& _system_filter;
98 std::vector<VectorType> _vec_v;
99 VectorType _vec_w;
101 std::vector<DataType> _c, _s, _q;
103 std::vector<std::vector<DataType>> _h;
104
105 public:
124 explicit GMRES(const MatrixType& matrix, const FilterType& filter, Index krylov_dim,
125 DataType inner_res_scale = DataType(0), std::shared_ptr<PrecondType> precond = nullptr) :
126 BaseClass("GMRES(" + stringify(krylov_dim) + ")", precond),
127 _system_matrix(matrix),
128 _system_filter(filter),
129 _krylov_dim(krylov_dim)
130 {
131 // we always need to compute defect norms
132 this->_force_def_norm_calc = true;
133
134 // set communicator by system matrix
135 this->_set_comm_by_matrix(matrix);
136 set_inner_res_scale(inner_res_scale);
137 }
138
139 explicit GMRES(const String& section_name, const PropertyMap* section,
140 const MatrixType& matrix, const FilterType& filter, std::shared_ptr<PrecondType> precond = nullptr) :
141 BaseClass("GMRES", section_name, section, precond),
142 _system_matrix(matrix),
143 _system_filter(filter),
144 _inner_res_scale(DataType(0))
145 {
146 // we always need to compute defect norms
147 this->_force_def_norm_calc = true;
148
149 // set communicator by system matrix
150 this->_set_comm_by_matrix(matrix);
151
152 // Set _inner_res_scale by parameter or use default value
153 auto inner_res_scale_p = section->query("inner_res_scale");
154 if(inner_res_scale_p.second && (!inner_res_scale_p.first.parse(this->_inner_res_scale) || (this->_inner_res_scale < DataType(0))))
155 throw ParseError(section_name + ".inner_res_scale", inner_res_scale_p.first, "a non-negative float");
156
157 // Check if we have set _krylov_vim
158 auto krylov_dim_p = section->query("krylov_dim");
159 if(!krylov_dim_p.second)
160 throw ParseError("GMRES config section is missing the mandatory krylov_dim!");
161
162 if(!krylov_dim_p.first.parse(this->_krylov_dim) || (this->_krylov_dim <= Index(0)))
163 throw ParseError(section_name + ".krylov_dim", krylov_dim_p.first, "a positive integer");
164
165 this->set_plot_name("GMRES("+stringify(_krylov_dim)+")");
166 }
167
171 virtual ~GMRES()
172 {
173 }
174
176 virtual String name() const override
177 {
178 return "GMRES";
179 }
180
187 virtual void force_defect_norm_calc(bool DOXY(force)) override
188 {
189 // force is already applied in constructor
190 }
191
192 virtual void init_symbolic() override
193 {
195
196 _c.reserve(_krylov_dim);
197 _s.reserve(_krylov_dim);
198 _q.reserve(_krylov_dim + 1);
199 _h.resize(_krylov_dim);
200
201 for(Index i(0); i < _krylov_dim; ++i)
202 {
203 _h.at(i).resize(i+1);
204 }
205
206 _vec_v.push_back(this->_system_matrix.create_vector_r());
207 _vec_w = this->_system_matrix.create_vector_r();
208 for(Index i(0); i < _krylov_dim; ++i)
209 {
210 _vec_v.push_back(this->_vec_v.front().clone(LAFEM::CloneMode::Layout));
211 }
212 }
213
214 virtual void done_symbolic() override
215 {
216 _vec_v.clear();
217 _vec_w.clear();
219 }
220
228 virtual void set_krylov_dim(Index krylov_dim)
229 {
230 XASSERT(krylov_dim > Index(0));
231 _krylov_dim = krylov_dim;
232 }
233
240 virtual void set_inner_res_scale(DataType inner_res_scale)
241 {
242 XASSERT(inner_res_scale >= DataType(0));
243 _inner_res_scale = inner_res_scale;
244 }
245
247 virtual Status apply(VectorType& vec_sol, const VectorType& vec_rhs) override
248 {
249
250 // BEGIN: Test auskommentiert
251 // save input rhs vector as initial defect
252 this->_vec_v.at(0).copy(vec_rhs);
253
254 // clear solution vector
255 vec_sol.format();
256
257 // apply
258 this->_status = _apply_intern(vec_sol, vec_rhs);
259 this->plot_summary();
260 return this->_status;
261 }
262
264 virtual Status correct(VectorType& vec_sol, const VectorType& vec_rhs) override
265 {
266 // compute initial defect
267 this->_system_matrix.apply(this->_vec_v.at(0), vec_sol, vec_rhs, -DataType(1));
268 this->_system_filter.filter_def(this->_vec_v.at(0));
269
270 // apply
271 this->_status = _apply_intern(vec_sol, vec_rhs);
272 this->plot_summary();
273 return this->_status;
274 }
275
276 protected:
277 virtual Status _apply_intern(VectorType& vec_sol, const VectorType& vec_rhs)
278 {
279 IterationStats pre_iter(*this);
280 Statistics::add_solver_expression(std::make_shared<ExpressionStartSolve>(this->name()));
281 const MatrixType& matrix(this->_system_matrix);
282 const FilterType& filter(this->_system_filter);
283
284 // compute initial defect
285 Status status = this->_set_initial_defect(this->_vec_v.at(0), vec_sol);
286
287 pre_iter.destroy();
288
289 // outer GMRES loop
290 while(status == Status::progress)
291 {
292 IterationStats stat(*this);
293
294 _q.clear();
295 _s.clear();
296 _c.clear();
297 _q.push_back(this->_def_cur);
298
299 // normalize v[0]
300 this->_vec_v.at(0).scale(this->_vec_v.at(0), DataType(1) / _q.back());
301
302 // inner GMRES loop
303 Index i(0);
304 while(i < this->_krylov_dim)
305 {
306 // apply preconditioner
307 if(!this->_apply_precond(this->_vec_w, this->_vec_v.at(i), filter))
308 {
309 stat.destroy();
310 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
311 return Status::aborted;
312 }
313 //filter.filter_cor(this->_vec_z.at(i));
314
315 // v[i+1] := A*w
316 matrix.apply(this->_vec_v.at(i+1), this->_vec_w);
317 filter.filter_def(this->_vec_v.at(i+1));
318
319 // Gram-Schmidt process
320 for(Index k(0); k <= i; ++k)
321 {
322 this->_h.at(i).at(k) = this->_vec_v.at(i+1).dot(this->_vec_v.at(k));
323 this->_vec_v.at(i+1).axpy(this->_vec_v.at(k), -this->_h.at(i).at(k));
324 }
325
326 // normalize v[i+1]
327 DataType alpha = this->_vec_v.at(i+1).norm2();
328 this->_vec_v.at(i+1).scale(this->_vec_v.at(i+1), DataType(1) / alpha);
329
330 // apply Givens rotations
331 for(Index k(0); k < i; ++k)
332 {
333 DataType t(this->_h.at(i).at(k));
334 this->_h.at(i).at(k ) = this->_c.at(k) * t + this->_s.at(k) * this->_h.at(i).at(k+1);
335 this->_h.at(i).at(k+1) = this->_s.at(k) * t - this->_c.at(k) * this->_h.at(i).at(k+1);
336 }
337
338 // compute beta
339 DataType beta = Math::sqrt(Math::sqr(this->_h.at(i).at(i)) + Math::sqr(alpha));
340
341 // compute next plane rotation
342 _s.push_back(alpha / beta);
343 _c.push_back(this->_h.at(i).at(i) / beta);
344
345 this->_h.at(i).at(i) = beta;
346 this->_q.push_back(this->_s.back() * this->_q.at(i));
347 this->_q.at(i) *= this->_c.back();
348
349 // push our new defect
350 if(++i < this->_krylov_dim)
351 {
352 // get the absolute defect
353 DataType def_cur = Math::abs(this->_q.back());
354
355 // did we diverge?
356 if((def_cur > this->_div_abs) || (def_cur > (this->_div_rel * this->_def_init)))
357 break;
358
359 // minimum number of iterations performed?
360 if(!(this->_num_iter < this->_min_iter))
361 {
362 // did we converge?
363 if((def_cur <= _inner_res_scale * this->_tol_abs) &&
364 ((def_cur <= _inner_res_scale * (this->_tol_rel * this->_def_init)) || (def_cur <= _inner_res_scale * this->_tol_abs_low) ))
365 break;
366
367 // maximum number of iterations performed?
368 if(this->_num_iter >= this->_max_iter)
369 break;
370 }
371
372 // set our pseudo defect
373 // increase iteration count
374 ++this->_num_iter;
375
376 // compute new defect
377 this->_def_cur = def_cur;
378
379 // plot?
380 if(this->_plot_iter())
381 {
382 String msg = this->_plot_name
383 + "* " + stringify(this->_num_iter).pad_front(this->_iter_digits)
384 + " : " + stringify_fp_sci(this->_def_cur);
385 this->_print_line(msg);
386 }
387 }
388 }
389
390 if(this->_plot_iter())
391 {
392 /* DEBUG */
393 String msg = this->_plot_name
394 + "* " + stringify(this->_num_iter).pad_front(this->_iter_digits)
395 + " : " + stringify_fp_sci(Math::abs(this->_q.back()))
396 + " <= " + stringify_fp_sci(this->_inner_res_scale * this->_tol_abs)
397 + " && " + stringify_fp_sci(this->_inner_res_scale * this->_tol_rel);
398 this->_print_line(msg);
399 }
400
401 Index n = Math::min(i, this->_krylov_dim);
402
403 // solve H*q = q
404 for(Index k(n); k > 0;)
405 {
406 --k;
407 this->_q.at(k) /= this->_h.at(k).at(k);
408 for(Index j(k); j > 0;)
409 {
410 --j;
411 this->_q.at(j) -= this->_h.at(k).at(j) * this->_q.at(k);
412 }
413 }
414
415 // update solution
416 this->_vec_w.format();
417 for(Index k(0); k < n; ++k)
418 this->_vec_w.axpy(this->_vec_v.at(k), this->_q.at(k));
419
420 // re-use v[0] for preconditioned correction
421 if(!this->_apply_precond(this->_vec_v.front(), this->_vec_w, filter))
422 {
423 stat.destroy();
424 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), Status::aborted, this->get_num_iter()));
425 return Status::aborted;
426 }
427 vec_sol.axpy(this->_vec_v.front());
428
429 // compute "real" residual
430 matrix.apply(this->_vec_v.at(0), vec_sol, vec_rhs, -DataType(1));
431 filter.filter_def(this->_vec_v.at(0));
432
433 // set the current defect
434 status = this->_set_new_defect(this->_vec_v.at(0), vec_sol);
435
436 }
437
438 // finished
439 Statistics::add_solver_expression(std::make_shared<ExpressionEndSolve>(this->name(), status, this->get_num_iter()));
440 return status;
441 }
442 }; // class GMRES<...>
443
466 template<typename Matrix_, typename Filter_>
467 inline std::shared_ptr<GMRES<Matrix_, Filter_>> new_gmres(
468 const Matrix_& matrix, const Filter_& filter, Index krylov_dim,
469 typename Matrix_::DataType inner_res_scale = typename Matrix_::DataType(0),
470 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
471 {
472 return std::make_shared<GMRES<Matrix_, Filter_>>(matrix, filter, krylov_dim, inner_res_scale, precond);
473 }
474
496 template<typename Matrix_, typename Filter_>
497 inline std::shared_ptr<GMRES<Matrix_, Filter_>> new_gmres(
498 const String& section_name, const PropertyMap* section,
499 const Matrix_& matrix, const Filter_& filter,
500 std::shared_ptr<SolverBase<typename Matrix_::VectorTypeL>> precond = nullptr)
501 {
502 return std::make_shared<GMRES<Matrix_, Filter_>>(section_name, section, matrix, filter, precond);
503 }
504 } // namespace Solver
505} // namespace FEAT
#define XASSERT(expr)
Assertion macro definition.
Definition: assertion.hpp:262
Class for parser related errors.
Definition: exception.hpp:132
A class organizing a tree of key-value pairs.
std::pair< String, bool > query(String key_path) const
Queries a value by its key path.
GMRES(k) solver implementation.
Definition: gmres.hpp:78
const FilterType & _system_filter
the filter for the solver
Definition: gmres.hpp:92
virtual String name() const override
Returns a descriptive string.
Definition: gmres.hpp:176
virtual void init_symbolic() override
Symbolic initialization method.
Definition: gmres.hpp:192
virtual void set_inner_res_scale(DataType inner_res_scale)
Sets the inner residual scale.
Definition: gmres.hpp:240
GMRES(const MatrixType &matrix, const FilterType &filter, Index krylov_dim, DataType inner_res_scale=DataType(0), std::shared_ptr< PrecondType > precond=nullptr)
Constructor.
Definition: gmres.hpp:124
virtual Status correct(VectorType &vec_sol, const VectorType &vec_rhs) override
Definition: gmres.hpp:264
DataType _inner_res_scale
inner pseudo-residual scaling factor, see the documentation of this class for details
Definition: gmres.hpp:96
virtual Status apply(VectorType &vec_sol, const VectorType &vec_rhs) override
Solver application method.
Definition: gmres.hpp:247
std::vector< VectorType > _vec_v
krylov basis vector
Definition: gmres.hpp:98
std::vector< std::vector< DataType > > _h
Hessenberg matrix.
Definition: gmres.hpp:103
virtual void set_krylov_dim(Index krylov_dim)
Sets the inner Krylov space dimension.
Definition: gmres.hpp:228
const MatrixType & _system_matrix
the matrix for the solver
Definition: gmres.hpp:90
virtual ~GMRES()
Empty virtual destructor.
Definition: gmres.hpp:171
virtual void done_symbolic() override
Symbolic finalization method.
Definition: gmres.hpp:214
std::vector< DataType > _c
Givens rotation coefficients.
Definition: gmres.hpp:101
virtual void force_defect_norm_calc(bool force) override
Forces the calculation of defect norms in every iteration (overridden)
Definition: gmres.hpp:187
Index _krylov_dim
krylov dimension
Definition: gmres.hpp:94
Helper class for iteration statistics collection.
Definition: base.hpp:392
bool _force_def_norm_calc
whether skipping of defect computation is allowed or not
Definition: iterative.hpp:249
void set_plot_name(const String &plot_name)
Sets the plot name of the solver.
Definition: iterative.hpp:517
DataType _def_init
initial defect
Definition: iterative.hpp:237
Index _min_iter
minimum number of iterations
Definition: iterative.hpp:227
Index _iter_digits
iteration count digits for plotting
Definition: iterative.hpp:243
Index get_num_iter() const
Returns number of performed iterations.
Definition: iterative.hpp:462
Status _status
current status of the solver
Definition: iterative.hpp:213
void _set_comm_by_matrix(const Matrix_ &matrix)
Sets the communicator for the solver from a matrix.
Definition: iterative.hpp:680
void _print_line(const String &line) const
Prints a line.
Definition: iterative.hpp:752
String _plot_name
name of the solver in plots
Definition: iterative.hpp:211
DataType _div_rel
relative divergence parameter
Definition: iterative.hpp:221
DataType _div_abs
absolute divergence parameter
Definition: iterative.hpp:223
virtual Status _set_new_defect(const VectorType &vec_def, const VectorType &vec_sol)
Internal function: sets the new (next) defect vector.
Definition: iterative.hpp:924
Index _num_iter
number of performed iterations
Definition: iterative.hpp:231
bool _plot_iter(Status st=Status::progress) const
Plot the current iteration?
Definition: iterative.hpp:720
DataType _tol_rel
relative tolerance parameter
Definition: iterative.hpp:215
DataType _def_cur
current defect
Definition: iterative.hpp:239
virtual void plot_summary() const
Plot a summary of the last solver run.
Definition: iterative.hpp:627
Index _max_iter
maximum number of iterations
Definition: iterative.hpp:229
DataType _tol_abs_low
absolute low tolerance parameter
Definition: iterative.hpp:219
virtual Status _set_initial_defect(const VectorType &vec_def, const VectorType &vec_sol)
Internal function: sets the initial defect vector.
Definition: iterative.hpp:802
DataType _tol_abs
absolute tolerance parameter
Definition: iterative.hpp:217
Abstract base-class for preconditioned iterative solvers.
Definition: iterative.hpp:1027
virtual void init_symbolic() override
Symbolic initialization method.
Definition: iterative.hpp:1086
bool _apply_precond(VectorType &vec_cor, const VectorType &vec_def, const Filter_ &filter)
Applies the preconditioner onto a defect vector.
Definition: iterative.hpp:1140
virtual void done_symbolic() override
Symbolic finalization method.
Definition: iterative.hpp:1110
Polymorphic solver interface.
Definition: base.hpp:183
String class implementation.
Definition: string.hpp:47
String pad_front(size_type len, char c=' ') const
Pads the front of the string up to a desired length.
Definition: string.hpp:393
T_ sqrt(T_ x)
Returns the square-root of a value.
Definition: math.hpp:300
T_ abs(T_ x)
Returns the absolute value.
Definition: math.hpp:275
T_ sqr(T_ x)
Returns the square of a value.
Definition: math.hpp:95
T_ min(T_ a, T_ b)
Returns the minimum of two values.
Definition: math.hpp:123
Status
Solver status return codes enumeration.
Definition: base.hpp:47
@ progress
continue iteration (internal use only)
@ aborted
premature abort (solver aborted due to internal errors or preconditioner failure)
std::shared_ptr< GMRES< Matrix_, Filter_ > > new_gmres(const Matrix_ &matrix, const Filter_ &filter, Index krylov_dim, typename Matrix_::DataType inner_res_scale=typename Matrix_::DataType(0), std::shared_ptr< SolverBase< typename Matrix_::VectorTypeL > > precond=nullptr)
Creates a new GMRES solver object.
Definition: gmres.hpp:467
FEAT namespace.
Definition: adjactor.hpp:12
String stringify(const T_ &item)
Converts an item into a String.
Definition: string.hpp:993
String stringify_fp_sci(DataType_ value, int precision=0, int width=0, bool sign=false)
Prints a floating point value to a string in scientific notation.
Definition: string.hpp:1137
std::uint64_t Index
Index data type.