ConjugateGradient.h 8.7 KB

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