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 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 
58 template <typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
59 bool gmres(const MatrixType& mat, const Rhs& rhs, Dest& x, const Preconditioner& precond, Index& iters,
60  const Index& restart, typename Dest::RealScalar& tol_error) {
61  using std::abs;
62  using std::sqrt;
63 
64  typedef typename Dest::RealScalar RealScalar;
65  typedef typename Dest::Scalar Scalar;
68 
69  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
70 
71  if (rhs.norm() <= considerAsZero) {
72  x.setZero();
73  tol_error = 0;
74  return true;
75  }
76 
77  RealScalar tol = tol_error;
78  const Index maxIters = iters;
79  iters = 0;
80 
81  const Index m = mat.rows();
82 
83  // residual and preconditioned residual
84  VectorType p0 = rhs - mat * x;
85  VectorType r0 = precond.solve(p0);
86 
87  const RealScalar r0Norm = r0.norm();
88 
89  // is initial guess already good enough?
90  if (r0Norm == 0) {
91  tol_error = 0;
92  return true;
93  }
94 
95  // storage for Hessenberg matrix and Householder data
96  FMatrixType H = FMatrixType::Zero(m, restart + 1);
99 
100  // storage for Jacobi rotations
101  std::vector<JacobiRotation<Scalar> > G(restart);
102 
103  // storage for temporaries
104  VectorType t(m), v(m), workspace(m), x_new(m);
105 
106  // generate first Householder vector
107  Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
109  r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
110  w(0) = Scalar(beta);
111 
112  for (Index k = 1; k <= restart; ++k) {
113  ++iters;
114 
115  v = VectorType::Unit(m, k - 1);
116 
117  // apply Householder reflections H_{1} ... H_{k-1} to v
118  // TODO: use a HouseholderSequence
119  for (Index i = k - 1; i >= 0; --i) {
120  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
121  }
122 
123  // apply matrix M to v: v = mat * v;
124  t.noalias() = mat * v;
125  v = precond.solve(t);
126 
127  // apply Householder reflections H_{k-1} ... H_{1} to v
128  // TODO: use a HouseholderSequence
129  for (Index i = 0; i < k; ++i) {
130  v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
131  }
132 
133  if (v.tail(m - k).norm() != 0.0) {
134  if (k <= restart) {
135  // generate new Householder vector
136  Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
137  v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
138 
139  // apply Householder reflection H_{k} to v
140  v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
141  }
142  }
143 
144  if (k > 1) {
145  for (Index i = 0; i < k - 1; ++i) {
146  // apply old Givens rotations to v
147  v.applyOnTheLeft(i, i + 1, G[i].adjoint());
148  }
149  }
150 
151  if (k < m && v(k) != (Scalar)0) {
152  // determine next Givens rotation
153  G[k - 1].makeGivens(v(k - 1), v(k));
154 
155  // apply Givens rotation to v and w
156  v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
157  w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
158  }
159 
160  // insert coefficients into upper matrix triangle
161  H.col(k - 1).head(k) = v.head(k);
162 
163  tol_error = abs(w(k)) / r0Norm;
164  bool stop = (k == m || tol_error < tol || iters == maxIters);
165 
166  if (stop || k == restart) {
167  // solve upper triangular system
168  Ref<VectorType> y = w.head(k);
169  H.topLeftCorner(k, k).template triangularView<Upper>().solveInPlace(y);
170 
171  // use Horner-like scheme to calculate solution vector
172  x_new.setZero();
173  for (Index i = k - 1; i >= 0; --i) {
174  x_new(i) += y(i);
175  // apply Householder reflection H_{i} to x_new
176  x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
177  }
178 
179  x += x_new;
180 
181  if (stop) {
182  return true;
183  } else {
184  k = 0;
185 
186  // reset data for restart
187  p0.noalias() = rhs - mat * x;
188  r0 = precond.solve(p0);
189 
190  // clear Hessenberg matrix and Householder data
191  H.setZero();
192  w.setZero();
193  tau.setZero();
194 
195  // generate first Householder vector
196  r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
197  w(0) = Scalar(beta);
198  }
199  }
200  }
201 
202  return false;
203 }
204 
205 } // namespace internal
206 
207 template <typename MatrixType_, typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
208 class GMRES;
209 
210 namespace internal {
211 
212 template <typename MatrixType_, typename Preconditioner_>
213 struct traits<GMRES<MatrixType_, Preconditioner_> > {
214  typedef MatrixType_ MatrixType;
215  typedef Preconditioner_ Preconditioner;
216 };
217 
218 } // namespace internal
219 
254 template <typename MatrixType_, typename Preconditioner_>
255 class GMRES : public IterativeSolverBase<GMRES<MatrixType_, Preconditioner_> > {
257  using Base::m_error;
258  using Base::m_info;
259  using Base::m_isInitialized;
260  using Base::m_iterations;
261  using Base::matrix;
262 
263  private:
265 
266  public:
267  using Base::_solve_impl;
268  typedef MatrixType_ MatrixType;
269  typedef typename MatrixType::Scalar Scalar;
271  typedef Preconditioner_ Preconditioner;
272 
273  public:
275  GMRES() : Base(), m_restart(30) {}
276 
287  template <typename MatrixDerived>
288  explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
289 
290  ~GMRES() {}
291 
295 
300 
302  template <typename Rhs, typename Dest>
303  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const {
308  }
309 
310  protected:
311 };
312 
313 } // end namespace Eigen
314 
315 #endif // EIGEN_GMRES_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
MatrixXf H
Definition: HessenbergDecomposition_matrixH.cpp:4
JacobiRotation< float > G
Definition: Jacobi_makeGivens.cpp:2
Vector3f p0
Definition: MatrixBase_all.cpp:2
RowVector3d w
Definition: Matrix_resize_int.cpp:3
void adjoint(const MatrixType &m)
Definition: adjoint.cpp:85
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
A GMRES solver for sparse square problems.
Definition: GMRES.h:255
MatrixType::RealScalar RealScalar
Definition: GMRES.h:270
~GMRES()
Definition: GMRES.h:290
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
MatrixType_ MatrixType
Definition: GMRES.h:268
GMRES(const EigenBase< MatrixDerived > &A)
Definition: GMRES.h:288
Preconditioner_ Preconditioner
Definition: GMRES.h:271
IterativeSolverBase< GMRES > Base
Definition: GMRES.h:256
Index m_restart
Definition: GMRES.h:264
RealScalar m_error
Definition: IterativeSolverBase.h:387
Index m_iterations
Definition: IterativeSolverBase.h:388
Index get_restart()
Definition: GMRES.h:294
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: GMRES.h:303
MatrixType::Scalar Scalar
Definition: GMRES.h:269
GMRES()
Definition: GMRES.h:275
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
void set_restart(const Index restart)
Definition: GMRES.h:299
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:124
Index maxIterations() const
Definition: IterativeSolverBase.h:251
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
RealScalar m_error
Definition: IterativeSolverBase.h:387
void _solve_impl(const Rhs &b, Dest &x) const
Definition: IterativeSolverBase.h:357
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:382
Index m_iterations
Definition: IterativeSolverBase.h:388
bool m_isInitialized
Definition: SparseSolverBase.h:110
RealScalar m_tolerance
Definition: IterativeSolverBase.h:385
GMRES< MatrixType_, Preconditioner_ > & derived()
Definition: SparseSolverBase.h:76
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:264
Index rows() const
Definition: SparseMatrix.h:159
Definition: restart2.cpp:8
#define min(a, b)
Definition: datatypes.h:22
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
Eigen::DenseIndex ret
Definition: level1_cplx_impl.h:43
int * m
Definition: level2_cplx_impl.h:294
Scalar beta
Definition: level2_cplx_impl.h:36
char char char int int * k
Definition: level2_impl.h:374
@ Rhs
Definition: TensorContractionMapper.h:20
const Scalar & y
Definition: RandomImpl.h:36
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:59
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
unsigned Preconditioner
----------------------—Domain Properties------------------------—
Definition: space_time_oscillating_cylinder.cc:725
Definition: Eigen_Colamd.h:49
double Zero
Definition: pseudosolid_node_update_elements.cc:35
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: EigenBase.h:33
Definition: ForwardDeclarations.h:21
Definition: fft_test_shared.h:66