1
// This file is part of Eigen, a lightweight C++ template library
4
// Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5
// Copyright (C) 2012 Kolja Brix <brix@igpm.rwth-aaachen.de>
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.
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.
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.
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/>.
34
* Generalized Minimal Residual Algorithm based on the
35
* Arnoldi algorithm implemented with Householder reflections.
38
* \param mat matrix of linear system of equations
39
* \param Rhs right hand side vector of linear system of equations
40
* \param x on input: initial guess, on output: solution
41
* \param precond preconditioner used
42
* \param iters on input: maximum number of iterations to perform
43
* on output: number of iterations performed
44
* \param restart number of iterations for a restart
45
* \param tol_error on input: residual tolerance
46
* on output: residuum achieved
48
* \sa IterativeMethods::bicgstab()
51
* For references, please see:
53
* Saad, Y. and Schultz, M. H.
54
* GMRES: A Generalized Minimal Residual Algorithm for Solving Nonsymmetric Linear Systems.
55
* SIAM J.Sci.Stat.Comp. 7, 1986, pp. 856 - 869.
58
* Iterative Methods for Sparse Linear Systems.
59
* Society for Industrial and Applied Mathematics, Philadelphia, 2003.
62
* Implementations of the GMRES method.
63
* Comput.Phys.Comm. 53, 1989, pp. 311 - 320.
66
* Implementation of the GMRES Method using Householder Transformations.
67
* SIAM J.Sci.Stat.Comp. 9, 1988, pp. 152 - 163.
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) {
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;
83
RealScalar tol = tol_error;
84
const int maxIters = iters;
87
const int m = mat.rows();
89
VectorType p0 = rhs - mat*x;
90
VectorType r0 = precond.solve(p0);
91
// RealScalar r0_sqnorm = r0.squaredNorm();
93
VectorType w = VectorType::Zero(restart + 1);
95
FMatrixType H = FMatrixType::Zero(m, restart + 1);
96
VectorType tau = VectorType::Zero(restart + 1);
97
std::vector < JacobiRotation < Scalar > > G(restart);
99
// generate first Householder vector
102
r0.makeHouseholder(e, tau.coeffRef(0), beta);
104
H.bottomLeftCorner(m - 1, 1) = e;
106
for (int k = 1; k <= restart; ++k) {
110
VectorType v = VectorType::Unit(m, k - 1), workspace(m);
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());
117
// apply matrix M to v: v = mat * v;
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());
126
if (v.tail(m - k).norm() != 0.0) {
130
// generate new Householder vector
131
VectorType e(m - k - 1);
133
v.tail(m - k).makeHouseholder(e, tau.coeffRef(k), beta);
134
H.col(k).tail(m - k - 1) = e;
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());
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());
149
if (k<m && v(k) != (Scalar) 0) {
150
// determine next Givens rotation
151
G[k - 1].makeGivens(v(k - 1), v(k));
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());
159
// insert coefficients into upper matrix triangle
160
H.col(k - 1).head(k) = v.head(k);
162
bool stop=(k==m || abs(w(k)) < tol || iters == maxIters);
164
if (stop || k == restart) {
166
// solve upper triangular system
167
VectorType y = w.head(k);
168
H.topLeftCorner(k, k).template triangularView < Eigen::Upper > ().solveInPlace(y);
170
// use Horner-like scheme to calculate solution vector
171
VectorType x_new = y(k - 1) * VectorType::Unit(m, k - 1);
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());
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());
189
// reset data for a restart r0 = rhs - mat * x;
191
VectorType p1=precond.solve(p0);
193
// r0_sqnorm = r0.squaredNorm();
194
w = VectorType::Zero(restart + 1);
195
H = FMatrixType::Zero(m, restart + 1);
196
tau = VectorType::Zero(restart + 1);
198
// generate first Householder vector
200
r0.makeHouseholder(e, tau.coeffRef(0), beta);
202
H.bottomLeftCorner(m - 1, 1) = e;
218
template< typename _MatrixType,
219
typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
224
template< typename _MatrixType, typename _Preconditioner>
225
struct traits<GMRES<_MatrixType,_Preconditioner> >
227
typedef _MatrixType MatrixType;
228
typedef _Preconditioner Preconditioner;
233
/** \ingroup IterativeLinearSolvers_Module
234
* \brief A GMRES solver for sparse square problems
236
* This class allows to solve for A.x = b sparse linear problems using a generalized minimal
237
* residual method. The vectors x and b can be either dense or sparse.
239
* \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
240
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
242
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
243
* and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
244
* and NumTraits<Scalar>::epsilon() for the tolerance.
246
* This class can be used as the direct solver classes. Here is a typical usage example:
249
* VectorXd x(n), b(n);
250
* SparseMatrix<double> A(n,n);
252
* GMRES<SparseMatrix<double> > solver(A);
253
* x = solver.solve(b);
254
* std::cout << "#iterations: " << solver.iterations() << std::endl;
255
* std::cout << "estimated error: " << solver.error() << std::endl;
256
* // update b, and solve again
257
* x = solver.solve(b);
260
* By default the iterations start with x=0 as an initial guess of the solution.
261
* One can control the start using the solveWithGuess() method. Here is a step by
262
* step execution example starting with a random guess and printing the evolution
263
* of the estimated error:
265
* x = VectorXd::Random(n);
266
* solver.setMaxIterations(1);
269
* x = solver.solveWithGuess(b,x);
270
* std::cout << i << " : " << solver.error() << std::endl;
272
* } while (solver.info()!=Success && i<100);
274
* Note that such a step by step excution is slightly slower.
276
* \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
278
template< typename _MatrixType, typename _Preconditioner>
279
class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
281
typedef IterativeSolverBase<GMRES> Base;
282
using Base::mp_matrix;
284
using Base::m_iterations;
286
using Base::m_isInitialized;
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;
300
/** Default constructor. */
301
GMRES() : Base(), m_restart(30) {}
303
/** Initialize the solver with matrix \a A for further \c Ax=b solving.
305
* This constructor is a shortcut for the default constructor followed
306
* by a call to compute().
308
* \warning this class stores a reference to the matrix A as well as some
309
* precomputed values that depend on it. Therefore, if \a A is changed
310
* this class becomes invalid. Call compute() to update it with the new
311
* matrix A, or modify a copy of A.
313
GMRES(const MatrixType& A) : Base(A), m_restart(30) {}
317
/** Get the number of iterations after that a restart is performed.
319
int get_restart() { return m_restart; }
321
/** Set the number of iterations after that a restart is performed.
322
* \param restart number of iterations for a restarti, default is 30.
324
void set_restart(const int restart) { m_restart=restart; }
326
/** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A
327
* \a x0 as an initial solution.
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
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);
343
template<typename Rhs,typename Dest>
344
void _solveWithGuess(const Rhs& b, Dest& x) const
347
for(int j=0; j<b.cols(); ++j)
349
m_iterations = Base::maxIterations();
350
m_error = Base::m_tolerance;
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))
356
m_info = failed ? NumericalIssue
357
: m_error <= Base::m_tolerance ? Success
359
m_isInitialized = true;
363
template<typename Rhs,typename Dest>
364
void _solve(const Rhs& b, Dest& x) const
367
_solveWithGuess(b,x);
377
template<typename _MatrixType, typename _Preconditioner, typename Rhs>
378
struct solve_retval<GMRES<_MatrixType, _Preconditioner>, Rhs>
379
: solve_retval_base<GMRES<_MatrixType, _Preconditioner>, Rhs>
381
typedef GMRES<_MatrixType, _Preconditioner> Dec;
382
EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
384
template<typename Dest> void evalTo(Dest& dst) const
386
dec()._solve(rhs(),dst);
390
} // end namespace internal
392
} // end namespace Eigen
394
#endif // EIGEN_GMRES_H