blob: ff912094f1ae440ea0323dd18de4eaeb3142c712 [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// 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
Austin Schuh189376f2018-12-20 22:11:15 +110014namespace Eigen {
Brian Silverman72890c22015-09-19 14:37:37 -040015
16namespace internal {
17
18/**
Austin Schuh189376f2018-12-20 22:11:15 +110019* Generalized Minimal Residual Algorithm based on the
20* Arnoldi algorithm implemented with Householder reflections.
21*
22* Parameters:
23* \param mat matrix of linear system of equations
24* \param rhs right hand side vector of linear system of equations
25* \param x on input: initial guess, on output: solution
26* \param precond preconditioner used
27* \param iters on input: maximum number of iterations to perform
28* on output: number of iterations performed
29* \param restart number of iterations for a restart
30* \param tol_error on input: relative residual tolerance
31* on output: residuum achieved
32*
33* \sa IterativeMethods::bicgstab()
34*
35*
36* For references, please see:
37*
38* Saad, Y. and Schultz, M. H.
39* GMRES: A Generalized Minimal Residual Algorithm for Solving Nonsymmetric Linear Systems.
40* SIAM J.Sci.Stat.Comp. 7, 1986, pp. 856 - 869.
41*
42* Saad, Y.
43* Iterative Methods for Sparse Linear Systems.
44* Society for Industrial and Applied Mathematics, Philadelphia, 2003.
45*
46* Walker, H. F.
47* Implementations of the GMRES method.
48* Comput.Phys.Comm. 53, 1989, pp. 311 - 320.
49*
50* Walker, H. F.
51* Implementation of the GMRES Method using Householder Transformations.
52* SIAM J.Sci.Stat.Comp. 9, 1988, pp. 152 - 163.
53*
54*/
Brian Silverman72890c22015-09-19 14:37:37 -040055template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
56bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond,
Austin Schuh189376f2018-12-20 22:11:15 +110057 Index &iters, const Index &restart, typename Dest::RealScalar & tol_error) {
Brian Silverman72890c22015-09-19 14:37:37 -040058
Austin Schuh189376f2018-12-20 22:11:15 +110059 using std::sqrt;
60 using std::abs;
Brian Silverman72890c22015-09-19 14:37:37 -040061
Austin Schuh189376f2018-12-20 22:11:15 +110062 typedef typename Dest::RealScalar RealScalar;
63 typedef typename Dest::Scalar Scalar;
64 typedef Matrix < Scalar, Dynamic, 1 > VectorType;
65 typedef Matrix < Scalar, Dynamic, Dynamic, ColMajor> FMatrixType;
Brian Silverman72890c22015-09-19 14:37:37 -040066
Austin Schuhc55b0172022-02-20 17:52:35 -080067 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
68
69 if(rhs.norm() <= considerAsZero)
70 {
71 x.setZero();
72 tol_error = 0;
73 return true;
74 }
75
Austin Schuh189376f2018-12-20 22:11:15 +110076 RealScalar tol = tol_error;
77 const Index maxIters = iters;
78 iters = 0;
Brian Silverman72890c22015-09-19 14:37:37 -040079
Austin Schuh189376f2018-12-20 22:11:15 +110080 const Index m = mat.rows();
Brian Silverman72890c22015-09-19 14:37:37 -040081
Austin Schuh189376f2018-12-20 22:11:15 +110082 // residual and preconditioned residual
83 VectorType p0 = rhs - mat*x;
84 VectorType r0 = precond.solve(p0);
Brian Silverman72890c22015-09-19 14:37:37 -040085
Austin Schuh189376f2018-12-20 22:11:15 +110086 const RealScalar r0Norm = r0.norm();
Brian Silverman72890c22015-09-19 14:37:37 -040087
Austin Schuh189376f2018-12-20 22:11:15 +110088 // is initial guess already good enough?
89 if(r0Norm == 0)
90 {
91 tol_error = 0;
92 return true;
93 }
Brian Silverman72890c22015-09-19 14:37:37 -040094
Austin Schuh189376f2018-12-20 22:11:15 +110095 // storage for Hessenberg matrix and Householder data
96 FMatrixType H = FMatrixType::Zero(m, restart + 1);
97 VectorType w = VectorType::Zero(restart + 1);
98 VectorType tau = VectorType::Zero(restart + 1);
Brian Silverman72890c22015-09-19 14:37:37 -040099
Austin Schuh189376f2018-12-20 22:11:15 +1100100 // 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);
Brian Silverman72890c22015-09-19 14:37:37 -0400105
Austin Schuh189376f2018-12-20 22:11:15 +1100106 // generate first Householder vector
107 Ref<VectorType> H0_tail = H.col(0).tail(m - 1);
108 RealScalar beta;
109 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
110 w(0) = Scalar(beta);
111
112 for (Index k = 1; k <= restart; ++k)
113 {
114 ++iters;
Brian Silverman72890c22015-09-19 14:37:37 -0400115
Austin Schuh189376f2018-12-20 22:11:15 +1100116 v = VectorType::Unit(m, k - 1);
Brian Silverman72890c22015-09-19 14:37:37 -0400117
Austin Schuh189376f2018-12-20 22:11:15 +1100118 // apply Householder reflections H_{1} ... H_{k-1} to v
119 // TODO: use a HouseholderSequence
120 for (Index i = k - 1; i >= 0; --i) {
121 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
122 }
Brian Silverman72890c22015-09-19 14:37:37 -0400123
Austin Schuh189376f2018-12-20 22:11:15 +1100124 // apply matrix M to v: v = mat * v;
125 t.noalias() = mat * v;
126 v = precond.solve(t);
Brian Silverman72890c22015-09-19 14:37:37 -0400127
Austin Schuh189376f2018-12-20 22:11:15 +1100128 // apply Householder reflections H_{k-1} ... H_{1} to v
129 // TODO: use a HouseholderSequence
130 for (Index i = 0; i < k; ++i) {
131 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
132 }
Brian Silverman72890c22015-09-19 14:37:37 -0400133
Austin Schuh189376f2018-12-20 22:11:15 +1100134 if (v.tail(m - k).norm() != 0.0)
135 {
136 if (k <= restart)
137 {
138 // generate new Householder vector
139 Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1);
140 v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta);
Brian Silverman72890c22015-09-19 14:37:37 -0400141
Austin Schuh189376f2018-12-20 22:11:15 +1100142 // apply Householder reflection H_{k} to v
143 v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data());
144 }
145 }
Brian Silverman72890c22015-09-19 14:37:37 -0400146
Austin Schuh189376f2018-12-20 22:11:15 +1100147 if (k > 1)
148 {
149 for (Index i = 0; i < k - 1; ++i)
150 {
151 // apply old Givens rotations to v
152 v.applyOnTheLeft(i, i + 1, G[i].adjoint());
153 }
154 }
Brian Silverman72890c22015-09-19 14:37:37 -0400155
Austin Schuh189376f2018-12-20 22:11:15 +1100156 if (k<m && v(k) != (Scalar) 0)
157 {
158 // determine next Givens rotation
159 G[k - 1].makeGivens(v(k - 1), v(k));
Brian Silverman72890c22015-09-19 14:37:37 -0400160
Austin Schuh189376f2018-12-20 22:11:15 +1100161 // apply Givens rotation to v and w
162 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
163 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint());
164 }
Brian Silverman72890c22015-09-19 14:37:37 -0400165
Austin Schuh189376f2018-12-20 22:11:15 +1100166 // insert coefficients into upper matrix triangle
167 H.col(k-1).head(k) = v.head(k);
Brian Silverman72890c22015-09-19 14:37:37 -0400168
Austin Schuh189376f2018-12-20 22:11:15 +1100169 tol_error = abs(w(k)) / r0Norm;
170 bool stop = (k==m || tol_error < tol || iters == maxIters);
Brian Silverman72890c22015-09-19 14:37:37 -0400171
Austin Schuh189376f2018-12-20 22:11:15 +1100172 if (stop || k == restart)
173 {
174 // solve upper triangular system
175 Ref<VectorType> y = w.head(k);
176 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y);
Brian Silverman72890c22015-09-19 14:37:37 -0400177
Austin Schuh189376f2018-12-20 22:11:15 +1100178 // use Horner-like scheme to calculate solution vector
179 x_new.setZero();
180 for (Index i = k - 1; i >= 0; --i)
181 {
182 x_new(i) += y(i);
183 // apply Householder reflection H_{i} to x_new
184 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data());
185 }
Brian Silverman72890c22015-09-19 14:37:37 -0400186
Austin Schuh189376f2018-12-20 22:11:15 +1100187 x += x_new;
Brian Silverman72890c22015-09-19 14:37:37 -0400188
Austin Schuh189376f2018-12-20 22:11:15 +1100189 if(stop)
190 {
191 return true;
192 }
193 else
194 {
195 k=0;
Brian Silverman72890c22015-09-19 14:37:37 -0400196
Austin Schuh189376f2018-12-20 22:11:15 +1100197 // reset data for restart
198 p0.noalias() = rhs - mat*x;
199 r0 = precond.solve(p0);
Brian Silverman72890c22015-09-19 14:37:37 -0400200
Austin Schuh189376f2018-12-20 22:11:15 +1100201 // clear Hessenberg matrix and Householder data
202 H.setZero();
203 w.setZero();
204 tau.setZero();
Brian Silverman72890c22015-09-19 14:37:37 -0400205
Austin Schuh189376f2018-12-20 22:11:15 +1100206 // generate first Householder vector
207 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta);
208 w(0) = Scalar(beta);
209 }
210 }
211 }
Brian Silverman72890c22015-09-19 14:37:37 -0400212
Austin Schuh189376f2018-12-20 22:11:15 +1100213 return false;
Brian Silverman72890c22015-09-19 14:37:37 -0400214
215}
216
217}
218
219template< typename _MatrixType,
220 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
221class GMRES;
222
223namespace internal {
224
225template< typename _MatrixType, typename _Preconditioner>
226struct traits<GMRES<_MatrixType,_Preconditioner> >
227{
228 typedef _MatrixType MatrixType;
229 typedef _Preconditioner Preconditioner;
230};
231
232}
233
234/** \ingroup IterativeLinearSolvers_Module
235 * \brief A GMRES solver for sparse square problems
236 *
237 * This class allows to solve for A.x = b sparse linear problems using a generalized minimal
238 * residual method. The vectors x and b can be either dense or sparse.
239 *
240 * \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
241 * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
242 *
243 * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
244 * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
245 * and NumTraits<Scalar>::epsilon() for the tolerance.
Austin Schuh189376f2018-12-20 22:11:15 +1100246 *
Brian Silverman72890c22015-09-19 14:37:37 -0400247 * This class can be used as the direct solver classes. Here is a typical usage example:
248 * \code
249 * int n = 10000;
250 * VectorXd x(n), b(n);
251 * SparseMatrix<double> A(n,n);
252 * // fill A and b
253 * GMRES<SparseMatrix<double> > solver(A);
254 * x = solver.solve(b);
255 * std::cout << "#iterations: " << solver.iterations() << std::endl;
256 * std::cout << "estimated error: " << solver.error() << std::endl;
257 * // update b, and solve again
258 * x = solver.solve(b);
259 * \endcode
Austin Schuh189376f2018-12-20 22:11:15 +1100260 *
Brian Silverman72890c22015-09-19 14:37:37 -0400261 * By default the iterations start with x=0 as an initial guess of the solution.
262 * One can control the start using the solveWithGuess() method.
263 *
Austin Schuh189376f2018-12-20 22:11:15 +1100264 * GMRES can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink.
265 *
Brian Silverman72890c22015-09-19 14:37:37 -0400266 * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
267 */
268template< typename _MatrixType, typename _Preconditioner>
269class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> >
270{
271 typedef IterativeSolverBase<GMRES> Base;
Austin Schuh189376f2018-12-20 22:11:15 +1100272 using Base::matrix;
Brian Silverman72890c22015-09-19 14:37:37 -0400273 using Base::m_error;
274 using Base::m_iterations;
275 using Base::m_info;
276 using Base::m_isInitialized;
Austin Schuh189376f2018-12-20 22:11:15 +1100277
Brian Silverman72890c22015-09-19 14:37:37 -0400278private:
Austin Schuh189376f2018-12-20 22:11:15 +1100279 Index m_restart;
280
Brian Silverman72890c22015-09-19 14:37:37 -0400281public:
Austin Schuh189376f2018-12-20 22:11:15 +1100282 using Base::_solve_impl;
Brian Silverman72890c22015-09-19 14:37:37 -0400283 typedef _MatrixType MatrixType;
284 typedef typename MatrixType::Scalar Scalar;
Brian Silverman72890c22015-09-19 14:37:37 -0400285 typedef typename MatrixType::RealScalar RealScalar;
286 typedef _Preconditioner Preconditioner;
287
288public:
289
290 /** Default constructor. */
291 GMRES() : Base(), m_restart(30) {}
292
293 /** Initialize the solver with matrix \a A for further \c Ax=b solving.
Austin Schuh189376f2018-12-20 22:11:15 +1100294 *
Brian Silverman72890c22015-09-19 14:37:37 -0400295 * This constructor is a shortcut for the default constructor followed
296 * by a call to compute().
Austin Schuh189376f2018-12-20 22:11:15 +1100297 *
Brian Silverman72890c22015-09-19 14:37:37 -0400298 * \warning this class stores a reference to the matrix A as well as some
299 * precomputed values that depend on it. Therefore, if \a A is changed
300 * this class becomes invalid. Call compute() to update it with the new
301 * matrix A, or modify a copy of A.
302 */
Austin Schuh189376f2018-12-20 22:11:15 +1100303 template<typename MatrixDerived>
304 explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
Brian Silverman72890c22015-09-19 14:37:37 -0400305
306 ~GMRES() {}
Austin Schuh189376f2018-12-20 22:11:15 +1100307
Brian Silverman72890c22015-09-19 14:37:37 -0400308 /** Get the number of iterations after that a restart is performed.
309 */
Austin Schuh189376f2018-12-20 22:11:15 +1100310 Index get_restart() { return m_restart; }
311
Brian Silverman72890c22015-09-19 14:37:37 -0400312 /** Set the number of iterations after that a restart is performed.
313 * \param restart number of iterations for a restarti, default is 30.
314 */
Austin Schuh189376f2018-12-20 22:11:15 +1100315 void set_restart(const Index restart) { m_restart=restart; }
316
Brian Silverman72890c22015-09-19 14:37:37 -0400317 /** \internal */
318 template<typename Rhs,typename Dest>
Austin Schuhc55b0172022-02-20 17:52:35 -0800319 void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
Austin Schuh189376f2018-12-20 22:11:15 +1100320 {
Austin Schuhc55b0172022-02-20 17:52:35 -0800321 m_iterations = Base::maxIterations();
322 m_error = Base::m_tolerance;
323 bool ret = internal::gmres(matrix(), b, x, Base::m_preconditioner, m_iterations, m_restart, m_error);
324 m_info = (!ret) ? NumericalIssue
Austin Schuh189376f2018-12-20 22:11:15 +1100325 : m_error <= Base::m_tolerance ? Success
326 : NoConvergence;
Brian Silverman72890c22015-09-19 14:37:37 -0400327 }
328
329protected:
330
331};
332
Brian Silverman72890c22015-09-19 14:37:37 -0400333} // end namespace Eigen
334
335#endif // EIGEN_GMRES_H