GMRES.h
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, 2014 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  Index &iters, const Index &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 Index maxIters = iters;
69  iters = 0;
70 
71  const Index m = mat.rows();
72 
73  // residual and preconditioned residual
74  VectorType p0 = rhs - mat*x;
75  VectorType r0 = precond.solve(p0);
76 
77  const RealScalar r0Norm = r0.norm();
78 
79  // is initial guess already good enough?
80  if(r0Norm == 0)
81  {
82  tol_error = 0;
83  return true;
84  }
85 
86  // storage for Hessenberg matrix and Householder data
87  FMatrixType H = FMatrixType::Zero(m, restart + 1);
88  VectorType w = VectorType::Zero(restart + 1);
89  VectorType tau = VectorType::Zero(restart + 1);
90 
91  // storage for Jacobi rotations
92  std::vector < JacobiRotation < Scalar > > G(restart);
93 
94  // storage for temporaries
95  VectorType t(m), v(m), workspace(m), x_new(m);
96 
97  // generate first Householder vector
98  Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
99  RealScalar beta;
100  r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
101  w(0) = Scalar(beta);
102 
103  for (Index k = 1; k <= restart; ++k)
104  {
105  ++iters;
106 
107  v = VectorType::Unit(m, k - 1);
108 
109  // apply Householder reflections H_{1} ... H_{k-1} to v
110  // TODO: use a HouseholderSequence
111  for (Index i = k - 1; i >= 0; --i) {
112  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
113  }
114 
115  // apply matrix M to v: v = mat * v;
116  t.noalias() = mat * v;
117  v = precond.solve(t);
118 
119  // apply Householder reflections H_{k-1} ... H_{1} to v
120  // TODO: use a HouseholderSequence
121  for (Index i = 0; i < k; ++i) {
122  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
123  }
124 
125  if (v.tail(m - k).norm() != 0.0)
126  {
127  if (k <= restart)
128  {
129  // generate new Householder vector
130  Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
131  v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
132 
133  // apply Householder reflection H_{k} to v
134  v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
135  }
136  }
137 
138  if (k > 1)
139  {
140  for (Index i = 0; i < k - 1; ++i)
141  {
142  // apply old Givens rotations to v
143  v.applyOnTheLeft(i, i + 1, G[i].adjoint());
144  }
145  }
146 
147  if (k<m && v(k) != (Scalar) 0)
148  {
149  // determine next Givens rotation
150  G[k - 1].makeGivens(v(k - 1), v(k));
151 
152  // apply Givens rotation to v and w
153  v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
154  w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
155  }
156 
157  // insert coefficients into upper matrix triangle
158  H.col(k-1).head(k) = v.head(k);
159 
160  bool stop = (k==m || abs(w(k)) < tol * r0Norm || iters == maxIters);
161 
162  if (stop || k == restart)
163  {
164  // solve upper triangular system
165  Ref<VectorType> y = w.head(k);
166  H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
167 
168  // use Horner-like scheme to calculate solution vector
169  x_new.setZero();
170  for (Index i = k - 1; i >= 0; --i)
171  {
172  x_new(i) += y(i);
173  // apply Householder reflection H_{i} to x_new
174  x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
175  }
176 
177  x += x_new;
178 
179  if(stop)
180  {
181  return true;
182  }
183  else
184  {
185  k=0;
186 
187  // reset data for restart
188  p0.noalias() = rhs - mat*x;
189  r0 = precond.solve(p0);
190 
191  // clear Hessenberg matrix and Householder data
192  H.setZero();
193  w.setZero();
194  tau.setZero();
195 
196  // generate first Householder vector
197  r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
198  w(0) = Scalar(beta);
199  }
200  }
201  }
202 
203  return false;
204 
205 }
206 
207 }
208 
209 template< typename _MatrixType,
210  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
211 class GMRES;
212 
213 namespace internal {
214 
215 template< typename _MatrixType, typename _Preconditioner>
216 struct traits<GMRES<_MatrixType,_Preconditioner> >
217 {
218  typedef _MatrixType MatrixType;
219  typedef _Preconditioner Preconditioner;
220 };
221 
222 }
223 
258 template< typename _MatrixType, typename _Preconditioner>
259 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
260 {
261  typedef IterativeSolverBase<GMRES> Base;
262  using Base::matrix;
263  using Base::m_error;
264  using Base::m_iterations;
265  using Base::m_info;
266  using Base::m_isInitialized;
267 
268 private:
269  Index m_restart;
270 
271 public:
272  using Base::_solve_impl;
273  typedef _MatrixType MatrixType;
274  typedef typename MatrixType::Scalar Scalar;
275  typedef typename MatrixType::RealScalar RealScalar;
276  typedef _Preconditioner Preconditioner;
277 
278 public:
279 
281  GMRES() : Base(), m_restart(30) {}
282 
293  template<typename MatrixDerived>
294  explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
295 
296  ~GMRES() {}
297 
300  Index get_restart() { return m_restart; }
301 
305  void set_restart(const Index restart) { m_restart=restart; }
306 
308  template<typename Rhs,typename Dest>
309  void _solve_with_guess_impl(const Rhs& b, Dest& x) const
310  {
311  bool failed = false;
312  for(Index j=0; j<b.cols(); ++j)
313  {
314  m_iterations = Base::maxIterations();
315  m_error = Base::m_tolerance;
316 
317  typename Dest::ColXpr xj(x,j);
318  if(!internal::gmres(matrix(), b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
319  failed = true;
320  }
321  m_info = failed ? NumericalIssue
322  : m_error <= Base::m_tolerance ? Success
323  : NoConvergence;
324  m_isInitialized = true;
325  }
326 
328  template<typename Rhs,typename Dest>
329  void _solve_impl(const Rhs& b, MatrixBase<Dest> &x) const
330  {
331  x = b;
332  if(x.squaredNorm() == 0) return; // Check Zero right hand side
333  _solve_with_guess_impl(b,x.derived());
334  }
335 
336 protected:
337 
338 };
339 
340 } // end namespace Eigen
341 
342 #endif // EIGEN_GMRES_H
GMRES()
Definition: GMRES.h:281
Index get_restart()
Definition: GMRES.h:300
Namespace containing all symbols from the Eigen library.
Definition: CXX11Meta.h:13
GMRES(const EigenBase< MatrixDerived > &A)
Definition: GMRES.h:294
void set_restart(const Index restart)
Definition: GMRES.h:305
A GMRES solver for sparse square problems.
Definition: GMRES.h:211