ConjugateGradient.h (8887B)
1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2011-2014 Gael Guennebaud <gael.guennebaud@inria.fr> 5 // 6 // This Source Code Form is subject to the terms of the Mozilla 7 // Public License v. 2.0. If a copy of the MPL was not distributed 8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 9 10 #ifndef EIGEN_CONJUGATE_GRADIENT_H 11 #define EIGEN_CONJUGATE_GRADIENT_H 12 13 namespace Eigen { 14 15 namespace internal { 16 17 /** \internal Low-level conjugate gradient algorithm 18 * \param mat The matrix A 19 * \param rhs The right hand side vector b 20 * \param x On input and initial solution, on output the computed solution. 21 * \param precond A preconditioner being able to efficiently solve for an 22 * approximation of Ax=b (regardless of b) 23 * \param iters On input the max number of iteration, on output the number of performed iterations. 24 * \param tol_error On input the tolerance error, on output an estimation of the relative error. 25 */ 26 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner> 27 EIGEN_DONT_INLINE 28 void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x, 29 const Preconditioner& precond, Index& iters, 30 typename Dest::RealScalar& tol_error) 31 { 32 using std::sqrt; 33 using std::abs; 34 typedef typename Dest::RealScalar RealScalar; 35 typedef typename Dest::Scalar Scalar; 36 typedef Matrix<Scalar,Dynamic,1> VectorType; 37 38 RealScalar tol = tol_error; 39 Index maxIters = iters; 40 41 Index n = mat.cols(); 42 43 VectorType residual = rhs - mat * x; //initial residual 44 45 RealScalar rhsNorm2 = rhs.squaredNorm(); 46 if(rhsNorm2 == 0) 47 { 48 x.setZero(); 49 iters = 0; 50 tol_error = 0; 51 return; 52 } 53 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)(); 54 RealScalar threshold = numext::maxi(RealScalar(tol*tol*rhsNorm2),considerAsZero); 55 RealScalar residualNorm2 = residual.squaredNorm(); 56 if (residualNorm2 < threshold) 57 { 58 iters = 0; 59 tol_error = sqrt(residualNorm2 / rhsNorm2); 60 return; 61 } 62 63 VectorType p(n); 64 p = precond.solve(residual); // initial search direction 65 66 VectorType z(n), tmp(n); 67 RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM 68 Index i = 0; 69 while(i < maxIters) 70 { 71 tmp.noalias() = mat * p; // the bottleneck of the algorithm 72 73 Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir 74 x += alpha * p; // update solution 75 residual -= alpha * tmp; // update residual 76 77 residualNorm2 = residual.squaredNorm(); 78 if(residualNorm2 < threshold) 79 break; 80 81 z = precond.solve(residual); // approximately solve for "A z = residual" 82 83 RealScalar absOld = absNew; 84 absNew = numext::real(residual.dot(z)); // update the absolute value of r 85 RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction 86 p = z + beta * p; // update search direction 87 i++; 88 } 89 tol_error = sqrt(residualNorm2 / rhsNorm2); 90 iters = i; 91 } 92 93 } 94 95 template< typename _MatrixType, int _UpLo=Lower, 96 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> > 97 class ConjugateGradient; 98 99 namespace internal { 100 101 template< typename _MatrixType, int _UpLo, typename _Preconditioner> 102 struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> > 103 { 104 typedef _MatrixType MatrixType; 105 typedef _Preconditioner Preconditioner; 106 }; 107 108 } 109 110 /** \ingroup IterativeLinearSolvers_Module 111 * \brief A conjugate gradient solver for sparse (or dense) self-adjoint problems 112 * 113 * This class allows to solve for A.x = b linear problems using an iterative conjugate gradient algorithm. 114 * The matrix A must be selfadjoint. The matrix A and the vectors x and b can be either dense or sparse. 115 * 116 * \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix. 117 * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower, 118 * \c Upper, or \c Lower|Upper in which the full matrix entries will be considered. 119 * Default is \c Lower, best performance is \c Lower|Upper. 120 * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner 121 * 122 * \implsparsesolverconcept 123 * 124 * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() 125 * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations 126 * and NumTraits<Scalar>::epsilon() for the tolerance. 127 * 128 * The tolerance corresponds to the relative residual error: |Ax-b|/|b| 129 * 130 * \b Performance: Even though the default value of \c _UpLo is \c Lower, significantly higher performance is 131 * achieved when using a complete matrix and \b Lower|Upper as the \a _UpLo template parameter. Moreover, in this 132 * case multi-threading can be exploited if the user code is compiled with OpenMP enabled. 133 * See \ref TopicMultiThreading for details. 134 * 135 * This class can be used as the direct solver classes. Here is a typical usage example: 136 \code 137 int n = 10000; 138 VectorXd x(n), b(n); 139 SparseMatrix<double> A(n,n); 140 // fill A and b 141 ConjugateGradient<SparseMatrix<double>, Lower|Upper> cg; 142 cg.compute(A); 143 x = cg.solve(b); 144 std::cout << "#iterations: " << cg.iterations() << std::endl; 145 std::cout << "estimated error: " << cg.error() << std::endl; 146 // update b, and solve again 147 x = cg.solve(b); 148 \endcode 149 * 150 * By default the iterations start with x=0 as an initial guess of the solution. 151 * One can control the start using the solveWithGuess() method. 152 * 153 * ConjugateGradient can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink. 154 * 155 * \sa class LeastSquaresConjugateGradient, class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner 156 */ 157 template< typename _MatrixType, int _UpLo, typename _Preconditioner> 158 class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> > 159 { 160 typedef IterativeSolverBase<ConjugateGradient> Base; 161 using Base::matrix; 162 using Base::m_error; 163 using Base::m_iterations; 164 using Base::m_info; 165 using Base::m_isInitialized; 166 public: 167 typedef _MatrixType MatrixType; 168 typedef typename MatrixType::Scalar Scalar; 169 typedef typename MatrixType::RealScalar RealScalar; 170 typedef _Preconditioner Preconditioner; 171 172 enum { 173 UpLo = _UpLo 174 }; 175 176 public: 177 178 /** Default constructor. */ 179 ConjugateGradient() : Base() {} 180 181 /** Initialize the solver with matrix \a A for further \c Ax=b solving. 182 * 183 * This constructor is a shortcut for the default constructor followed 184 * by a call to compute(). 185 * 186 * \warning this class stores a reference to the matrix A as well as some 187 * precomputed values that depend on it. Therefore, if \a A is changed 188 * this class becomes invalid. Call compute() to update it with the new 189 * matrix A, or modify a copy of A. 190 */ 191 template<typename MatrixDerived> 192 explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {} 193 194 ~ConjugateGradient() {} 195 196 /** \internal */ 197 template<typename Rhs,typename Dest> 198 void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const 199 { 200 typedef typename Base::MatrixWrapper MatrixWrapper; 201 typedef typename Base::ActualMatrixType ActualMatrixType; 202 enum { 203 TransposeInput = (!MatrixWrapper::MatrixFree) 204 && (UpLo==(Lower|Upper)) 205 && (!MatrixType::IsRowMajor) 206 && (!NumTraits<Scalar>::IsComplex) 207 }; 208 typedef typename internal::conditional<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType const&>::type RowMajorWrapper; 209 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(MatrixWrapper::MatrixFree,UpLo==(Lower|Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY); 210 typedef typename internal::conditional<UpLo==(Lower|Upper), 211 RowMajorWrapper, 212 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type 213 >::type SelfAdjointWrapper; 214 215 m_iterations = Base::maxIterations(); 216 m_error = Base::m_tolerance; 217 218 RowMajorWrapper row_mat(matrix()); 219 internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b, x, Base::m_preconditioner, m_iterations, m_error); 220 m_info = m_error <= Base::m_tolerance ? Success : NoConvergence; 221 } 222 223 protected: 224 225 }; 226 227 } // end namespace Eigen 228 229 #endif // EIGEN_CONJUGATE_GRADIENT_H