Home | History | Annotate | Download | only in IterativeLinearSolvers
      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 (at) 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   RealScalar threshold = tol*tol*rhsNorm2;
     54   RealScalar residualNorm2 = residual.squaredNorm();
     55   if (residualNorm2 < threshold)
     56   {
     57     iters = 0;
     58     tol_error = sqrt(residualNorm2 / rhsNorm2);
     59     return;
     60   }
     61 
     62   VectorType p(n);
     63   p = precond.solve(residual);      // initial search direction
     64 
     65   VectorType z(n), tmp(n);
     66   RealScalar absNew = numext::real(residual.dot(p));  // the square of the absolute value of r scaled by invM
     67   Index i = 0;
     68   while(i < maxIters)
     69   {
     70     tmp.noalias() = mat * p;                    // the bottleneck of the algorithm
     71 
     72     Scalar alpha = absNew / p.dot(tmp);         // the amount we travel on dir
     73     x += alpha * p;                             // update solution
     74     residual -= alpha * tmp;                    // update residual
     75 
     76     residualNorm2 = residual.squaredNorm();
     77     if(residualNorm2 < threshold)
     78       break;
     79 
     80     z = precond.solve(residual);                // approximately solve for "A z = residual"
     81 
     82     RealScalar absOld = absNew;
     83     absNew = numext::real(residual.dot(z));     // update the absolute value of r
     84     RealScalar beta = absNew / absOld;          // calculate the Gram-Schmidt value used to create the new search direction
     85     p = z + beta * p;                           // update search direction
     86     i++;
     87   }
     88   tol_error = sqrt(residualNorm2 / rhsNorm2);
     89   iters = i;
     90 }
     91 
     92 }
     93 
     94 template< typename _MatrixType, int _UpLo=Lower,
     95           typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
     96 class ConjugateGradient;
     97 
     98 namespace internal {
     99 
    100 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
    101 struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
    102 {
    103   typedef _MatrixType MatrixType;
    104   typedef _Preconditioner Preconditioner;
    105 };
    106 
    107 }
    108 
    109 /** \ingroup IterativeLinearSolvers_Module
    110   * \brief A conjugate gradient solver for sparse (or dense) self-adjoint problems
    111   *
    112   * This class allows to solve for A.x = b linear problems using an iterative conjugate gradient algorithm.
    113   * The matrix A must be selfadjoint. The matrix A and the vectors x and b can be either dense or sparse.
    114   *
    115   * \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
    116   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower,
    117   *               \c Upper, or \c Lower|Upper in which the full matrix entries will be considered.
    118   *               Default is \c Lower, best performance is \c Lower|Upper.
    119   * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
    120   *
    121   * \implsparsesolverconcept
    122   *
    123   * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
    124   * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
    125   * and NumTraits<Scalar>::epsilon() for the tolerance.
    126   *
    127   * The tolerance corresponds to the relative residual error: |Ax-b|/|b|
    128   *
    129   * \b Performance: Even though the default value of \c _UpLo is \c Lower, significantly higher performance is
    130   * achieved when using a complete matrix and \b Lower|Upper as the \a _UpLo template parameter. Moreover, in this
    131   * case multi-threading can be exploited if the user code is compiled with OpenMP enabled.
    132   * See \ref TopicMultiThreading for details.
    133   *
    134   * This class can be used as the direct solver classes. Here is a typical usage example:
    135     \code
    136     int n = 10000;
    137     VectorXd x(n), b(n);
    138     SparseMatrix<double> A(n,n);
    139     // fill A and b
    140     ConjugateGradient<SparseMatrix<double>, Lower|Upper> cg;
    141     cg.compute(A);
    142     x = cg.solve(b);
    143     std::cout << "#iterations:     " << cg.iterations() << std::endl;
    144     std::cout << "estimated error: " << cg.error()      << std::endl;
    145     // update b, and solve again
    146     x = cg.solve(b);
    147     \endcode
    148   *
    149   * By default the iterations start with x=0 as an initial guess of the solution.
    150   * One can control the start using the solveWithGuess() method.
    151   *
    152   * ConjugateGradient can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink.
    153   *
    154   * \sa class LeastSquaresConjugateGradient, class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
    155   */
    156 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
    157 class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
    158 {
    159   typedef IterativeSolverBase<ConjugateGradient> Base;
    160   using Base::matrix;
    161   using Base::m_error;
    162   using Base::m_iterations;
    163   using Base::m_info;
    164   using Base::m_isInitialized;
    165 public:
    166   typedef _MatrixType MatrixType;
    167   typedef typename MatrixType::Scalar Scalar;
    168   typedef typename MatrixType::RealScalar RealScalar;
    169   typedef _Preconditioner Preconditioner;
    170 
    171   enum {
    172     UpLo = _UpLo
    173   };
    174 
    175 public:
    176 
    177   /** Default constructor. */
    178   ConjugateGradient() : Base() {}
    179 
    180   /** Initialize the solver with matrix \a A for further \c Ax=b solving.
    181     *
    182     * This constructor is a shortcut for the default constructor followed
    183     * by a call to compute().
    184     *
    185     * \warning this class stores a reference to the matrix A as well as some
    186     * precomputed values that depend on it. Therefore, if \a A is changed
    187     * this class becomes invalid. Call compute() to update it with the new
    188     * matrix A, or modify a copy of A.
    189     */
    190   template<typename MatrixDerived>
    191   explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
    192 
    193   ~ConjugateGradient() {}
    194 
    195   /** \internal */
    196   template<typename Rhs,typename Dest>
    197   void _solve_with_guess_impl(const Rhs& b, Dest& x) const
    198   {
    199     typedef typename Base::MatrixWrapper MatrixWrapper;
    200     typedef typename Base::ActualMatrixType ActualMatrixType;
    201     enum {
    202       TransposeInput  =   (!MatrixWrapper::MatrixFree)
    203                       &&  (UpLo==(Lower|Upper))
    204                       &&  (!MatrixType::IsRowMajor)
    205                       &&  (!NumTraits<Scalar>::IsComplex)
    206     };
    207     typedef typename internal::conditional<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType const&>::type RowMajorWrapper;
    208     EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(MatrixWrapper::MatrixFree,UpLo==(Lower|Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
    209     typedef typename internal::conditional<UpLo==(Lower|Upper),
    210                                            RowMajorWrapper,
    211                                            typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
    212                                           >::type SelfAdjointWrapper;
    213     m_iterations = Base::maxIterations();
    214     m_error = Base::m_tolerance;
    215 
    216     for(Index j=0; j<b.cols(); ++j)
    217     {
    218       m_iterations = Base::maxIterations();
    219       m_error = Base::m_tolerance;
    220 
    221       typename Dest::ColXpr xj(x,j);
    222       RowMajorWrapper row_mat(matrix());
    223       internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
    224     }
    225 
    226     m_isInitialized = true;
    227     m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
    228   }
    229 
    230   /** \internal */
    231   using Base::_solve_impl;
    232   template<typename Rhs,typename Dest>
    233   void _solve_impl(const MatrixBase<Rhs>& b, Dest& x) const
    234   {
    235     x.setZero();
    236     _solve_with_guess_impl(b.derived(),x);
    237   }
    238 
    239 protected:
    240 
    241 };
    242 
    243 } // end namespace Eigen
    244 
    245 #endif // EIGEN_CONJUGATE_GRADIENT_H
    246