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 namespace Eigen { 15 16 namespace internal { 17 18 /** 19 * 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 */ 55 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner> 56 bool gmres(const MatrixType & mat, const Rhs & rhs, Dest & x, const Preconditioner & precond, 57 Index &iters, const Index &restart, typename Dest::RealScalar & tol_error) { 58 59 using std::sqrt; 60 using std::abs; 61 62 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; 66 67 RealScalar tol = tol_error; 68 const Index maxIters = iters; 69 iters = 0; 70 71 const Index m = mat.rows(); 72 73 // residual and preconditioned residual 74 VectorType p0 = rhs - mat*x; 75 VectorType r0 = precond.solve(p0); 76 77 const RealScalar r0Norm = r0.norm(); 78 79 // is initial guess already good enough? 80 if(r0Norm == 0) 81 { 82 tol_error = 0; 83 return true; 84 } 85 86 // storage for Hessenberg matrix and Householder data 87 FMatrixType H = FMatrixType::Zero(m, restart + 1); 88 VectorType w = VectorType::Zero(restart + 1); 89 VectorType tau = VectorType::Zero(restart + 1); 90 91 // storage for Jacobi rotations 92 std::vector < JacobiRotation < Scalar > > G(restart); 93 94 // storage for temporaries 95 VectorType t(m), v(m), workspace(m), x_new(m); 96 97 // generate first Householder vector 98 Ref<VectorType> H0_tail = H.col(0).tail(m - 1); 99 RealScalar beta; 100 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta); 101 w(0) = Scalar(beta); 102 103 for (Index k = 1; k <= restart; ++k) 104 { 105 ++iters; 106 107 v = VectorType::Unit(m, k - 1); 108 109 // apply Householder reflections H_{1} ... H_{k-1} to v 110 // TODO: use a HouseholderSequence 111 for (Index i = k - 1; i >= 0; --i) { 112 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data()); 113 } 114 115 // apply matrix M to v: v = mat * v; 116 t.noalias() = mat * v; 117 v = precond.solve(t); 118 119 // apply Householder reflections H_{k-1} ... H_{1} to v 120 // TODO: use a HouseholderSequence 121 for (Index i = 0; i < k; ++i) { 122 v.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data()); 123 } 124 125 if (v.tail(m - k).norm() != 0.0) 126 { 127 if (k <= restart) 128 { 129 // generate new Householder vector 130 Ref<VectorType> Hk_tail = H.col(k).tail(m - k - 1); 131 v.tail(m - k).makeHouseholder(Hk_tail, tau.coeffRef(k), beta); 132 133 // apply Householder reflection H_{k} to v 134 v.tail(m - k).applyHouseholderOnTheLeft(Hk_tail, tau.coeffRef(k), workspace.data()); 135 } 136 } 137 138 if (k > 1) 139 { 140 for (Index i = 0; i < k - 1; ++i) 141 { 142 // apply old Givens rotations to v 143 v.applyOnTheLeft(i, i + 1, G[i].adjoint()); 144 } 145 } 146 147 if (k<m && v(k) != (Scalar) 0) 148 { 149 // determine next Givens rotation 150 G[k - 1].makeGivens(v(k - 1), v(k)); 151 152 // apply Givens rotation to v and w 153 v.applyOnTheLeft(k - 1, k, G[k - 1].adjoint()); 154 w.applyOnTheLeft(k - 1, k, G[k - 1].adjoint()); 155 } 156 157 // insert coefficients into upper matrix triangle 158 H.col(k-1).head(k) = v.head(k); 159 160 tol_error = abs(w(k)) / r0Norm; 161 bool stop = (k==m || tol_error < tol || iters == maxIters); 162 163 if (stop || k == restart) 164 { 165 // solve upper triangular system 166 Ref<VectorType> y = w.head(k); 167 H.topLeftCorner(k, k).template triangularView <Upper>().solveInPlace(y); 168 169 // use Horner-like scheme to calculate solution vector 170 x_new.setZero(); 171 for (Index i = k - 1; i >= 0; --i) 172 { 173 x_new(i) += y(i); 174 // apply Householder reflection H_{i} to x_new 175 x_new.tail(m - i).applyHouseholderOnTheLeft(H.col(i).tail(m - i - 1), tau.coeffRef(i), workspace.data()); 176 } 177 178 x += x_new; 179 180 if(stop) 181 { 182 return true; 183 } 184 else 185 { 186 k=0; 187 188 // reset data for restart 189 p0.noalias() = rhs - mat*x; 190 r0 = precond.solve(p0); 191 192 // clear Hessenberg matrix and Householder data 193 H.setZero(); 194 w.setZero(); 195 tau.setZero(); 196 197 // generate first Householder vector 198 r0.makeHouseholder(H0_tail, tau.coeffRef(0), beta); 199 w(0) = Scalar(beta); 200 } 201 } 202 } 203 204 return false; 205 206 } 207 208 } 209 210 template< typename _MatrixType, 211 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> > 212 class GMRES; 213 214 namespace internal { 215 216 template< typename _MatrixType, typename _Preconditioner> 217 struct traits<GMRES<_MatrixType,_Preconditioner> > 218 { 219 typedef _MatrixType MatrixType; 220 typedef _Preconditioner Preconditioner; 221 }; 222 223 } 224 225 /** \ingroup IterativeLinearSolvers_Module 226 * \brief A GMRES solver for sparse square problems 227 * 228 * This class allows to solve for A.x = b sparse linear problems using a generalized minimal 229 * residual method. The vectors x and b can be either dense or sparse. 230 * 231 * \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix. 232 * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner 233 * 234 * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() 235 * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations 236 * and NumTraits<Scalar>::epsilon() for the tolerance. 237 * 238 * This class can be used as the direct solver classes. Here is a typical usage example: 239 * \code 240 * int n = 10000; 241 * VectorXd x(n), b(n); 242 * SparseMatrix<double> A(n,n); 243 * // fill A and b 244 * GMRES<SparseMatrix<double> > solver(A); 245 * x = solver.solve(b); 246 * std::cout << "#iterations: " << solver.iterations() << std::endl; 247 * std::cout << "estimated error: " << solver.error() << std::endl; 248 * // update b, and solve again 249 * x = solver.solve(b); 250 * \endcode 251 * 252 * By default the iterations start with x=0 as an initial guess of the solution. 253 * One can control the start using the solveWithGuess() method. 254 * 255 * GMRES can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink. 256 * 257 * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner 258 */ 259 template< typename _MatrixType, typename _Preconditioner> 260 class GMRES : public IterativeSolverBase<GMRES<_MatrixType,_Preconditioner> > 261 { 262 typedef IterativeSolverBase<GMRES> Base; 263 using Base::matrix; 264 using Base::m_error; 265 using Base::m_iterations; 266 using Base::m_info; 267 using Base::m_isInitialized; 268 269 private: 270 Index m_restart; 271 272 public: 273 using Base::_solve_impl; 274 typedef _MatrixType MatrixType; 275 typedef typename MatrixType::Scalar Scalar; 276 typedef typename MatrixType::RealScalar RealScalar; 277 typedef _Preconditioner Preconditioner; 278 279 public: 280 281 /** Default constructor. */ 282 GMRES() : Base(), m_restart(30) {} 283 284 /** Initialize the solver with matrix \a A for further \c Ax=b solving. 285 * 286 * This constructor is a shortcut for the default constructor followed 287 * by a call to compute(). 288 * 289 * \warning this class stores a reference to the matrix A as well as some 290 * precomputed values that depend on it. Therefore, if \a A is changed 291 * this class becomes invalid. Call compute() to update it with the new 292 * matrix A, or modify a copy of A. 293 */ 294 template<typename MatrixDerived> 295 explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {} 296 297 ~GMRES() {} 298 299 /** Get the number of iterations after that a restart is performed. 300 */ 301 Index get_restart() { return m_restart; } 302 303 /** Set the number of iterations after that a restart is performed. 304 * \param restart number of iterations for a restarti, default is 30. 305 */ 306 void set_restart(const Index restart) { m_restart=restart; } 307 308 /** \internal */ 309 template<typename Rhs,typename Dest> 310 void _solve_with_guess_impl(const Rhs& b, Dest& x) const 311 { 312 bool failed = false; 313 for(Index j=0; j<b.cols(); ++j) 314 { 315 m_iterations = Base::maxIterations(); 316 m_error = Base::m_tolerance; 317 318 typename Dest::ColXpr xj(x,j); 319 if(!internal::gmres(matrix(), b.col(j), xj, Base::m_preconditioner, m_iterations, m_restart, m_error)) 320 failed = true; 321 } 322 m_info = failed ? NumericalIssue 323 : m_error <= Base::m_tolerance ? Success 324 : NoConvergence; 325 m_isInitialized = true; 326 } 327 328 /** \internal */ 329 template<typename Rhs,typename Dest> 330 void _solve_impl(const Rhs& b, MatrixBase<Dest> &x) const 331 { 332 x = b; 333 if(x.squaredNorm() == 0) return; // Check Zero right hand side 334 _solve_with_guess_impl(b,x.derived()); 335 } 336 337 protected: 338 339 }; 340 341 } // end namespace Eigen 342 343 #endif // EIGEN_GMRES_H 344