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 Kolja Brix <brix@igpm.rwth-aaachen.de>
6 //
7 // Eigen is free software; you can redistribute it and/or
8 // modify it under the terms of the GNU Lesser General Public
9 // License as published by the Free Software Foundation; either
10 // version 3 of the License, or (at your option) any later version.
11 //
12 // Alternatively, you can redistribute it and/or
13 // modify it under the terms of the GNU General Public License as
14 // published by the Free Software Foundation; either version 2 of
15 // the License, or (at your option) any later version.
16 //
17 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20 // GNU General Public License for more details.
21 //
22 // You should have received a copy of the GNU Lesser General Public
23 // License and a copy of the GNU General Public License along with
24 // Eigen. If not, see <http://www.gnu.org/licenses/>.
25 
26 #ifndef EIGEN_GMRES_H
27 #define EIGEN_GMRES_H
28 
29 namespace Eigen {
30 
31 namespace internal {
32 
70 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
71 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
72  int &iters, const int &restart, typename Dest::RealScalar & tol_error) {
73 
74  using std::sqrt;
75  using std::abs;
76 
77  typedef typename Dest::RealScalar RealScalar;
78  typedef typename Dest::Scalar Scalar;
79  typedef Matrix < RealScalar, Dynamic, 1 > RealVectorType;
80  typedef Matrix < Scalar, Dynamic, 1 > VectorType;
81  typedef Matrix < Scalar, Dynamic, Dynamic > FMatrixType;
82 
83  RealScalar tol = tol_error;
84  const int maxIters = iters;
85  iters = 0;
86 
87  const int m = mat.rows();
88 
89  VectorType p0 = rhs - mat*x;
90  VectorType r0 = precond.solve(p0);
91 // RealScalar r0_sqnorm = r0.squaredNorm();
92 
93  VectorType w = VectorType::Zero(restart + 1);
94 
95  FMatrixType H = FMatrixType::Zero(m, restart + 1);
96  VectorType tau = VectorType::Zero(restart + 1);
97  std::vector < JacobiRotation < Scalar > > G(restart);
98 
99  // generate first Householder vector
100  VectorType e;
101  RealScalar beta;
102  r0.makeHouseholder(e, tau.coeffRef(0), beta);
103  w(0)=(Scalar) beta;
104  H.bottomLeftCorner(m - 1, 1) = e;
105 
106  for (int k = 1; k <= restart; ++k) {
107 
108  ++iters;
109 
110  VectorType v = VectorType::Unit(m, k - 1), workspace(m);
111 
112  // apply Householder reflections H_{1} ... H_{k-1} to v
113  for (int i = k - 1; i >= 0; --i) {
114  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
115  }
116 
117  // apply matrix M to v: v = mat * v;
118  VectorType t=mat*v;
119  v=precond.solve(t);
120 
121  // apply Householder reflections H_{k-1} ... H_{1} to v
122  for (int i = 0; i < k; ++i) {
123  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
124  }
125 
126  if (v.tail(m - k).norm() != 0.0) {
127 
128  if (k <= restart) {
129 
130  // generate new Householder vector
131  VectorType e(m - k - 1);
132  RealScalar beta;
133  v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
134  H.col(k).tail(m - k - 1) = e;
135 
136  // apply Householder reflection H_{k} to v
137  v.tail(m - k).applyHouseholderOnTheLeft(H.col(k).tail(m - k - 1), tau.coeffRef(k), workspace.data());
138 
139  }
140  }
141 
142  if (k > 1) {
143  for (int i = 0; i < k - 1; ++i) {
144  // apply old Givens rotations to v
145  v.applyOnTheLeft(i, i + 1, G[i].adjoint());
146  }
147  }
148 
149  if (k<m && v(k) != (Scalar) 0) {
150  // determine next Givens rotation
151  G[k - 1].makeGivens(v(k - 1), v(k));
152 
153  // apply Givens rotation to v and w
154  v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
155  w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
156 
157  }
158 
159  // insert coefficients into upper matrix triangle
160  H.col(k - 1).head(k) = v.head(k);
161 
162  bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
163 
164  if (stop || k == restart) {
165 
166  // solve upper triangular system
167  VectorType y = w.head(k);
168  H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
169 
170  // use Horner-like scheme to calculate solution vector
171  VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
172 
173  // apply Householder reflection H_{k} to x_new
174  x_new.tail(m - k + 1).applyHouseholderOnTheLeft(H.col(k - 1).tail(m - k), tau.coeffRef(k - 1), workspace.data());
175 
176  for (int i = k - 2; i >= 0; --i) {
177  x_new += y(i) * VectorType::Unit(m, i);
178  // apply Householder reflection H_{i} to x_new
179  x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
180  }
181 
182  x += x_new;
183 
184  if (stop) {
185  return true;
186  } else {
187  k=0;
188 
189  // reset data for a restart r0 = rhs - mat * x;
190  VectorType p0=mat*x;
191  VectorType p1=precond.solve(p0);
192  r0 = rhs - p1;
193 // r0_sqnorm = r0.squaredNorm();
194  w = VectorType::Zero(restart + 1);
195  H = FMatrixType::Zero(m, restart + 1);
196  tau = VectorType::Zero(restart + 1);
197 
198  // generate first Householder vector
199  RealScalar beta;
200  r0.makeHouseholder(e, tau.coeffRef(0), beta);
201  w(0)=(Scalar) beta;
202  H.bottomLeftCorner(m - 1, 1) = e;
203 
204  }
205 
206  }
207 
208 
209 
210  }
211 
212  return false;
213 
214 }
215 
216 }
217 
218 template< typename _MatrixType,
219  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
220 class GMRES;
221 
222 namespace internal {
223 
224 template< typename _MatrixType, typename _Preconditioner>
225 struct traits<GMRES<_MatrixType,_Preconditioner> >
226 {
227  typedef _MatrixType MatrixType;
228  typedef _Preconditioner Preconditioner;
229 };
230 
231 }
232 
278 template< typename _MatrixType, typename _Preconditioner>
279 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
280 {
281  typedef IterativeSolverBase<GMRES> Base;
282  using Base::mp_matrix;
283  using Base::m_error;
284  using Base::m_iterations;
285  using Base::m_info;
286  using Base::m_isInitialized;
287 
288 private:
289  int m_restart;
290 
291 public:
292  typedef _MatrixType MatrixType;
293  typedef typename MatrixType::Scalar Scalar;
294  typedef typename MatrixType::Index Index;
295  typedef typename MatrixType::RealScalar RealScalar;
296  typedef _Preconditioner Preconditioner;
297 
298 public:
299 
301  GMRES() : Base(), m_restart(30) {}
302 
313  GMRES(const MatrixType& A) : Base(A), m_restart(30) {}
314 
315  ~GMRES() {}
316 
319  int get_restart() { return m_restart; }
320 
324  void set_restart(const int restart) { m_restart=restart; }
325 
331  template<typename Rhs,typename Guess>
332  inline const internal::solve_retval_with_guess<GMRES, Rhs, Guess>
333  solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
334  {
335  eigen_assert(m_isInitialized && "GMRES is not initialized.");
336  eigen_assert(Base::rows()==b.rows()
337  && "GMRES::solve(): invalid number of rows of the right hand side matrix b");
338  return internal::solve_retval_with_guess
339  <GMRES, Rhs, Guess>(*this, b.derived(), x0);
340  }
341 
343  template<typename Rhs,typename Dest>
344  void _solveWithGuess(const Rhs& b, Dest& x) const
345  {
346  bool failed = false;
347  for(int j=0; j<b.cols(); ++j)
348  {
349  m_iterations = Base::maxIterations();
350  m_error = Base::m_tolerance;
351 
352  typename Dest::ColXpr xj(x,j);
353  if(!internal::gmres(*mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error))
354  failed = true;
355  }
356  m_info = failed ? NumericalIssue
357  : m_error <= Base::m_tolerance ? Success
358  : NoConvergence;
359  m_isInitialized = true;
360  }
361 
363  template<typename Rhs,typename Dest>
364  void _solve(const Rhs& b, Dest& x) const
365  {
366  x.setZero();
367  _solveWithGuess(b,x);
368  }
369 
370 protected:
371 
372 };
373 
374 
375 namespace internal {
376 
377  template<typename _MatrixType, typename _Preconditioner, typename Rhs>
378 struct solve_retval<GMRES<_MatrixType, _Preconditioner>, Rhs>
379  : solve_retval_base<GMRES<_MatrixType, _Preconditioner>, Rhs>
380 {
381  typedef GMRES<_MatrixType, _Preconditioner> Dec;
382  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
383 
384  template<typename Dest> void evalTo(Dest& dst) const
385  {
386  dec()._solve(rhs(),dst);
387  }
388 };
389 
390 } // end namespace internal
391 
392 } // end namespace Eigen
393 
394 #endif // EIGEN_GMRES_H