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, 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 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
57 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
58 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
59  Index &iters, const Index &restart, typename Dest::RealScalar & tol_error) {
60 
61  using std::sqrt;
62  using std::abs;
63 
64  typedef typename Dest::RealScalar RealScalar;
65  typedef typename Dest::Scalar Scalar;
66  typedef Matrix < Scalar, Dynamic, 1 > VectorType;
68 
69  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
70 
71  if(rhs.norm() <= considerAsZero)
72  {
73  x.setZero();
74  tol_error = 0;
75  return true;
76  }
77 
78  RealScalar tol = tol_error;
79  const Index maxIters = iters;
80  iters = 0;
81 
82  const Index m = mat.rows();
83 
84  // residual and preconditioned residual
85  VectorType p0 = rhs - mat*x;
86  VectorType r0 = precond.solve(p0);
87 
88  const RealScalar r0Norm = r0.norm();
89 
90  // is initial guess already good enough?
91  if(r0Norm == 0)
92  {
93  tol_error = 0;
94  return true;
95  }
96 
97  // storage for Hessenberg matrix and Householder data
98  FMatrixType H = FMatrixType::Zero(m, restart + 1);
99  VectorType w = VectorType::Zero(restart + 1);
100  VectorType tau = VectorType::Zero(restart + 1);
101 
102  // storage for Jacobi rotations
103  std::vector < JacobiRotation < Scalar > > G(restart);
104 
105  // storage for temporaries
106  VectorType t(m), v(m), workspace(m), x_new(m);
107 
108  // generate first Householder vector
109  Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
110  RealScalar beta;
111  r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
112  w(0) = Scalar(beta);
113 
114  for (Index k = 1; k <= restart; ++k)
115  {
116  ++iters;
117 
118  v = VectorType::Unit(m, k - 1);
119 
120  // apply Householder reflections H_{1} ... H_{k-1} to v
121  // TODO: use a HouseholderSequence
122  for (Index i = k - 1; i >= 0; --i) {
123  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
124  }
125 
126  // apply matrix M to v: v = mat * v;
127  t.noalias() = mat * v;
128  v = precond.solve(t);
129 
130  // apply Householder reflections H_{k-1} ... H_{1} to v
131  // TODO: use a HouseholderSequence
132  for (Index i = 0; i < k; ++i) {
133  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
134  }
135 
136  if (v.tail(m - k).norm() != 0.0)
137  {
138  if (k <= restart)
139  {
140  // generate new Householder vector
141  Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
142  v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
143 
144  // apply Householder reflection H_{k} to v
145  v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
146  }
147  }
148 
149  if (k > 1)
150  {
151  for (Index i = 0; i < k - 1; ++i)
152  {
153  // apply old Givens rotations to v
154  v.applyOnTheLeft(i, i + 1, G[i].adjoint());
155  }
156  }
157 
158  if (k<m && v(k) != (Scalar) 0)
159  {
160  // determine next Givens rotation
161  G[k - 1].makeGivens(v(k - 1), v(k));
162 
163  // apply Givens rotation to v and w
164  v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
165  w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
166  }
167 
168  // insert coefficients into upper matrix triangle
169  H.col(k-1).head(k) = v.head(k);
170 
171  tol_error = abs(w(k)) / r0Norm;
172  bool stop = (k==m || tol_error < tol || iters == maxIters);
173 
174  if (stop || k == restart)
175  {
176  // solve upper triangular system
177  Ref<VectorType> y = w.head(k);
178  H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
179 
180  // use Horner-like scheme to calculate solution vector
181  x_new.setZero();
182  for (Index i = k - 1; i >= 0; --i)
183  {
184  x_new(i) += y(i);
185  // apply Householder reflection H_{i} to x_new
186  x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
187  }
188 
189  x += x_new;
190 
191  if(stop)
192  {
193  return true;
194  }
195  else
196  {
197  k=0;
198 
199  // reset data for restart
200  p0.noalias() = rhs - mat*x;
201  r0 = precond.solve(p0);
202 
203  // clear Hessenberg matrix and Householder data
204  H.setZero();
205  w.setZero();
206  tau.setZero();
207 
208  // generate first Householder vector
209  r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
210  w(0) = Scalar(beta);
211  }
212  }
213  }
214 
215  return false;
216 
217 }
218 
219 }
220 
221 template< typename MatrixType_,
222  typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
223 class GMRES;
224 
225 namespace internal {
226 
227 template< typename MatrixType_, typename Preconditioner_>
228 struct traits<GMRES<MatrixType_,Preconditioner_> >
229 {
230  typedef MatrixType_ MatrixType;
231  typedef Preconditioner_ Preconditioner;
232 };
233 
234 }
235 
270 template< typename MatrixType_, typename Preconditioner_>
271 class GMRES : public IterativeSolverBase<GMRES<MatrixType_,Preconditioner_> >
272 {
274  using Base::matrix;
275  using Base::m_error;
276  using Base::m_iterations;
277  using Base::m_info;
278  using Base::m_isInitialized;
279 
280 private:
282 
283 public:
284  using Base::_solve_impl;
285  typedef MatrixType_ MatrixType;
286  typedef typename MatrixType::Scalar Scalar;
287  typedef typename MatrixType::RealScalar RealScalar;
288  typedef Preconditioner_ Preconditioner;
289 
290 public:
291 
293  GMRES() : Base(), m_restart(30) {}
294 
305  template<typename MatrixDerived>
306  explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
307 
308  ~GMRES() {}
309 
313 
317  void set_restart(const Index restart) { m_restart=restart; }
318 
320  template<typename Rhs,typename Dest>
321  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
322  {
326  m_info = (!ret) ? NumericalIssue
328  : NoConvergence;
329  }
330 
331 protected:
332 
333 };
334 
335 } // end namespace Eigen
336 
337 #endif // EIGEN_GMRES_H
Matrix3f m
Array< int, Dynamic, 1 > v
SparseMatrix< double > A(n, n)
int i
MatrixXf H
JacobiRotation< float > G
Vector3f p0
RowVector3d w
MatrixXf mat
Matrix< float, 1, Dynamic > MatrixType
A GMRES solver for sparse square problems.
Definition: GMRES.h:272
RealScalar m_error
MatrixType::RealScalar RealScalar
Definition: GMRES.h:287
MatrixType_ MatrixType
Definition: GMRES.h:285
GMRES(const EigenBase< MatrixDerived > &A)
Definition: GMRES.h:306
const ActualMatrixType & matrix() const
Preconditioner_ Preconditioner
Definition: GMRES.h:288
IterativeSolverBase< GMRES > Base
Definition: GMRES.h:273
ComputationInfo m_info
Index m_restart
Definition: GMRES.h:281
Index m_iterations
Index get_restart()
Definition: GMRES.h:312
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: GMRES.h:321
MatrixType::Scalar Scalar
Definition: GMRES.h:286
void set_restart(const Index restart)
Definition: GMRES.h:317
void _solve_impl(const Rhs &b, Dest &x) const
Preconditioner m_preconditioner
const ActualMatrixType & matrix() const
NumericalIssue
const Scalar & y
bool gmres(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, const Index &restart, typename Dest::RealScalar &tol_error)
Definition: GMRES.h:58
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Eigen::internal::remove_all_t< DerType >, typename Eigen::internal::traits< Eigen::internal::remove_all_t< DerType >>::Scalar, product) > sqrt(const Eigen::AutoDiffScalar< DerType > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
adouble abs(const adouble &x)
Definition: AdolcForward:74