10 #ifndef EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H    11 #define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H    26 template<
typename MatrixType, 
typename Rhs, 
typename Dest, 
typename Preconditioner>
    29                                      const Preconditioner& precond, 
Index& iters,
    30                                      typename Dest::RealScalar& tol_error)
    34   typedef typename Dest::RealScalar RealScalar;
    35   typedef typename Dest::Scalar Scalar;
    38   RealScalar tol = tol_error;
    39   Index maxIters = iters;
    43   VectorType residual        = rhs - mat * x;
    44   VectorType normal_residual = mat.adjoint() * residual;
    46   RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
    54   RealScalar threshold = tol*tol*rhsNorm2;
    55   RealScalar residualNorm2 = normal_residual.squaredNorm();
    56   if (residualNorm2 < threshold)
    59     tol_error = 
sqrt(residualNorm2 / rhsNorm2);
    64   p = precond.solve(normal_residual);                         
    66   VectorType z(n), tmp(m);
    67   RealScalar absNew = 
numext::real(normal_residual.dot(p));  
    71     tmp.noalias() = mat * p;
    73     Scalar alpha = absNew / tmp.squaredNorm();      
    75     residual -= alpha * tmp;                        
    76     normal_residual = mat.adjoint() * residual;     
    78     residualNorm2 = normal_residual.squaredNorm();
    79     if(residualNorm2 < threshold)
    82     z = precond.solve(normal_residual);             
    84     RealScalar absOld = absNew;
    86     RealScalar beta = absNew / absOld;              
    90   tol_error = 
sqrt(residualNorm2 / rhsNorm2);
    96 template< 
typename _MatrixType,
   102 template< 
typename _MatrixType, 
typename _Preconditioner>
   148 template< 
typename _MatrixType, 
typename _Preconditioner>
   154   using Base::m_iterations;
   156   using Base::m_isInitialized;
   159   typedef typename MatrixType::Scalar 
Scalar;
   178   template<
typename MatrixDerived>
   184   template<
typename Rhs,
typename Dest>
   187     m_iterations = Base::maxIterations();
   188     m_error = Base::m_tolerance;
   190     for(
Index j=0; j<b.cols(); ++j)
   192       m_iterations = Base::maxIterations();
   193       m_error = Base::m_tolerance;
   199     m_isInitialized = 
true;
   204   using Base::_solve_impl;
   205   template<
typename Rhs,
typename Dest>
   209     _solve_with_guess_impl(b.derived(),x);
   216 #endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
_Preconditioner Preconditioner
void _solve_with_guess_impl(const Rhs &b, Dest &x) const
EIGEN_DEVICE_FUNC RealReturnType real() const
IterativeSolverBase< LeastSquaresConjugateGradient > Base
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1, !IsRowMajor > ColXpr
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
A conjugate gradient solver for sparse (or dense) least-square problems. 
~LeastSquaresConjugateGradient()
_Preconditioner Preconditioner
Jacobi preconditioner for LeastSquaresConjugateGradient. 
#define EIGEN_DONT_INLINE
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API. 
MatrixType::RealScalar RealScalar
MatrixType::Scalar Scalar
EIGEN_DONT_INLINE void least_square_conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
void _solve_impl(const MatrixBase< Rhs > &b, Dest &x) const
LeastSquaresConjugateGradient()
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
Base class for linear iterative solvers. 
Base class for all dense matrices, vectors, and expressions.