GMRES.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2012 Kolja Brix <brix@igpm.rwth-aaachen.de>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_GMRES_H
12 #define EIGEN_GMRES_H
13 
14 namespace Eigen {
15 
16 namespace internal {
17 
55 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
56 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
57  int &iters, const int &restart, typename Dest::RealScalar & tol_error) {
58 
59  using std::sqrt;
60  using std::abs;
61 
62  typedef typename Dest::RealScalar RealScalar;
63  typedef typename Dest::Scalar Scalar;
64  typedef Matrix < Scalar, Dynamic, 1 > VectorType;
65  typedef Matrix < Scalar, Dynamic, Dynamic > FMatrixType;
66 
67  RealScalar tol = tol_error;
68  const int maxIters = iters;
69  iters = 0;
70 
71  const int m = mat.rows();
72 
73  VectorType p0 = rhs - mat*x;
74  VectorType r0 = precond.solve(p0);
75 // RealScalar r0_sqnorm = r0.squaredNorm();
76 
77  VectorType w = VectorType::Zero(restart + 1);
78 
79  FMatrixType H = FMatrixType::Zero(m, restart + 1);
80  VectorType tau = VectorType::Zero(restart + 1);
81  std::vector < JacobiRotation < Scalar > > G(restart);
82 
83  // generate first Householder vector
84  VectorType e;
85  RealScalar beta;
86  r0.makeHouseholder(e, tau.coeffRef(0), beta);
87  w(0)=(Scalar) beta;
88  H.bottomLeftCorner(m - 1, 1) = e;
89 
90  for (int k = 1; k <= restart; ++k) {
91 
92  ++iters;
93 
94  VectorType v = VectorType::Unit(m, k - 1), workspace(m);
95 
96  // apply Householder reflections H_{1} ... H_{k-1} to v
97  for (int i = k - 1; i >= 0; --i) {
98  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
99  }
100 
101  // apply matrix M to v: v = mat * v;
102  VectorType t=mat*v;
103  v=precond.solve(t);
104 
105  // apply Householder reflections H_{k-1} ... H_{1} to v
106  for (int i = 0; i < k; ++i) {
107  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
108  }
109 
110  if (v.tail(m - k).norm() != 0.0) {
111 
112  if (k <= restart) {
113 
114  // generate new Householder vector
115  VectorType e(m - k - 1);
116  RealScalar beta;
117  v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
118  H.col(k).tail(m - k - 1) = e;
119 
120  // apply Householder reflection H_{k} to v
121  v.tail(m - k).applyHouseholderOnTheLeft(H.col(k).tail(m - k - 1), tau.coeffRef(k), workspace.data());
122 
123  }
124  }
125 
126  if (k > 1) {
127  for (int i = 0; i < k - 1; ++i) {
128  // apply old Givens rotations to v
129  v.applyOnTheLeft(i, i + 1, G[i].adjoint());
130  }
131  }
132 
133  if (k<m && v(k) != (Scalar) 0) {
134  // determine next Givens rotation
135  G[k - 1].makeGivens(v(k - 1), v(k));
136 
137  // apply Givens rotation to v and w
138  v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
139  w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
140 
141  }
142 
143  // insert coefficients into upper matrix triangle
144  H.col(k - 1).head(k) = v.head(k);
145 
146  bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
147 
148  if (stop || k == restart) {
149 
150  // solve upper triangular system
151  VectorType y = w.head(k);
152  H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
153 
154  // use Horner-like scheme to calculate solution vector
155  VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
156 
157  // apply Householder reflection H_{k} to x_new
158  x_new.tail(m - k + 1).applyHouseholderOnTheLeft(H.col(k - 1).tail(m - k), tau.coeffRef(k - 1), workspace.data());
159 
160  for (int i = k - 2; i >= 0; --i) {
161  x_new += y(i) * VectorType::Unit(m, i);
162  // apply Householder reflection H_{i} to x_new
163  x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
164  }
165 
166  x += x_new;
167 
168  if (stop) {
169  return true;
170  } else {
171  k=0;
172 
173  // reset data for a restart r0 = rhs - mat * x;
174  VectorType p0=mat*x;
175  VectorType p1=precond.solve(p0);
176  r0 = rhs - p1;
177 // r0_sqnorm = r0.squaredNorm();
178  w = VectorType::Zero(restart + 1);
179  H = FMatrixType::Zero(m, restart + 1);
180  tau = VectorType::Zero(restart + 1);
181 
182  // generate first Householder vector
183  RealScalar beta;
184  r0.makeHouseholder(e, tau.coeffRef(0), beta);
185  w(0)=(Scalar) beta;
186  H.bottomLeftCorner(m - 1, 1) = e;
187 
188  }
189 
190  }
191 
192 
193 
194  }
195 
196  return false;
197 
198 }
199 
200 }
201 
202 template< typename _MatrixType,
203  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
204 class GMRES;
205 
206 namespace internal {
207 
208 template< typename _MatrixType, typename _Preconditioner>
209 struct traits<GMRES<_MatrixType,_Preconditioner> >
210 {
211  typedef _MatrixType MatrixType;
212  typedef _Preconditioner Preconditioner;
213 };
214 
215 }
216 
262 template< typename _MatrixType, typename _Preconditioner>
263 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
264 {
266  using Base::mp_matrix;
267  using Base::m_error;
268  using Base::m_iterations;
269  using Base::m_info;
270  using Base::m_isInitialized;
271 
272 private:
274 
275 public:
276  typedef _MatrixType MatrixType;
277  typedef typename MatrixType::Scalar Scalar;
278  typedef typename MatrixType::Index Index;
279  typedef typename MatrixType::RealScalar RealScalar;
280  typedef _Preconditioner Preconditioner;
281 
282 public:
283 
285  GMRES() : Base(), m_restart(30) {}
286 
297  GMRES(const MatrixType& A) : Base(A), m_restart(30) {}
298 
299  ~GMRES() {}
300 
303  int get_restart() { return m_restart; }
304 
308  void set_restart(const int restart) { m_restart=restart; }
309 
315  template<typename Rhs,typename Guess>
317  solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
318  {
319  eigen_assert(m_isInitialized && "GMRES is not initialized.");
320  eigen_assert(Base::rows()==b.rows()
321  && "GMRES::solve(): invalid number of rows of the right hand side matrix b");
323  <GMRES, Rhs, Guess>(*this, b.derived(), x0);
324  }
325 
327  template<typename Rhs,typename Dest>
328  void _solveWithGuess(const Rhs& b, Dest& x) const
329  {
330  bool failed = false;
331  for(int j=0; j<b.cols(); ++j)
332  {
333  m_iterations = Base::maxIterations();
334  m_error = Base::m_tolerance;
335 
336  typename Dest::ColXpr xj(x,j);
337  if(!internal::gmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
338  failed = true;
339  }
340  m_info = failed ? NumericalIssue
341  : m_error <= Base::m_tolerance ? Success
342  : NoConvergence;
343  m_isInitialized = true;
344  }
345 
347  template<typename Rhs,typename Dest>
348  void _solve(const Rhs& b, Dest& x) const
349  {
350  x = b;
351  if(x.squaredNorm() == 0) return; // Check Zero right hand side
352  _solveWithGuess(b,x);
353  }
354 
355 protected:
356 
357 };
358 
359 
360 namespace internal {
361 
362  template<typename _MatrixType, typename _Preconditioner, typename Rhs>
363 struct solve_retval<GMRES<_MatrixType, _Preconditioner>, Rhs>
364  : solve_retval_base<GMRES<_MatrixType, _Preconditioner>, Rhs>
365 {
367  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
368 
369  template<typename Dest> void evalTo(Dest& dst) const
370  {
371  dec()._solve(rhs(),dst);
372  }
373 };
374 
375 } // end namespace internal
376 
377 } // end namespace Eigen
378 
379 #endif // EIGEN_GMRES_H
A preconditioner based on the digonal entries.
MatrixType::RealScalar RealScalar
Definition: GMRES.h:279
GMRES(const MatrixType &A)
Definition: GMRES.h:297
IntermediateState sqrt(const Expression &arg)
void _solve(const Rhs &b, Dest &x) const
Definition: GMRES.h:348
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ColXpr
Definition: BlockMethods.h:15
bool gmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, int &iters, const int &restart, typename Dest::RealScalar &tol_error)
Definition: GMRES.h:56
_MatrixType MatrixType
Definition: GMRES.h:276
EIGEN_STRONG_INLINE Index rows() const
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
int get_restart()
Definition: GMRES.h:303
_Preconditioner Preconditioner
Definition: GMRES.h:280
void _solveWithGuess(const Rhs &b, Dest &x) const
Definition: GMRES.h:328
#define v
void rhs(const real_t *x, real_t *f)
IterativeSolverBase< GMRES > Base
Definition: GMRES.h:265
A GMRES solver for sparse square problems.
Definition: GMRES.h:204
int m_restart
Definition: GMRES.h:273
MatrixType::Index Index
Definition: GMRES.h:278
MatrixType::Scalar Scalar
Definition: GMRES.h:277
const internal::solve_retval_with_guess< GMRES, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: GMRES.h:317
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
Definition: Solve.h:61
#define eigen_assert(x)
Base class for linear iterative solvers.
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
void set_restart(const int restart)
Definition: GMRES.h:308


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:39