10 #ifndef EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
11 #define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
13 #include "./InternalHeaderCheck.h"
28 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
30 void least_square_conjugate_gradient(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
31 const Preconditioner& precond,
Index& iters,
32 typename Dest::RealScalar& tol_error)
36 typedef typename Dest::RealScalar RealScalar;
37 typedef typename Dest::Scalar Scalar;
38 typedef Matrix<Scalar,Dynamic,1> VectorType;
40 RealScalar tol = tol_error;
41 Index maxIters = iters;
43 Index m = mat.rows(), n = mat.cols();
45 VectorType residual = rhs - mat * x;
46 VectorType normal_residual = mat.adjoint() * residual;
48 RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
56 RealScalar threshold = tol*tol*rhsNorm2;
57 RealScalar residualNorm2 = normal_residual.squaredNorm();
58 if (residualNorm2 < threshold)
61 tol_error =
sqrt(residualNorm2 / rhsNorm2);
66 p = precond.solve(normal_residual);
68 VectorType z(n), tmp(m);
69 RealScalar absNew = numext::real(normal_residual.dot(p));
73 tmp.noalias() = mat * p;
75 Scalar alpha = absNew / tmp.squaredNorm();
77 residual -= alpha * tmp;
78 normal_residual.noalias() = mat.adjoint() * residual;
80 residualNorm2 = normal_residual.squaredNorm();
81 if(residualNorm2 < threshold)
84 z = precond.solve(normal_residual);
86 RealScalar absOld = absNew;
87 absNew = numext::real(normal_residual.dot(z));
88 RealScalar beta = absNew / absOld;
92 tol_error =
sqrt(residualNorm2 / rhsNorm2);
98 template<
typename MatrixType_,
99 typename Preconditioner_ = LeastSquareDiagonalPreconditioner<typename MatrixType_::Scalar> >
100 class LeastSquaresConjugateGradient;
104 template<
typename MatrixType_,
typename Preconditioner_>
105 struct traits<LeastSquaresConjugateGradient<MatrixType_,Preconditioner_> >
107 typedef MatrixType_ MatrixType;
108 typedef Preconditioner_ Preconditioner;
150 template<
typename MatrixType_,
typename Preconditioner_>
156 using Base::m_iterations;
158 using Base::m_isInitialized;
160 typedef MatrixType_ MatrixType;
161 typedef typename MatrixType::Scalar Scalar;
162 typedef typename MatrixType::RealScalar RealScalar;
163 typedef Preconditioner_ Preconditioner;
180 template<
typename MatrixDerived>
186 template<
typename Rhs,
typename Dest>
187 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const
190 m_error = Base::m_tolerance;
192 internal::least_square_conjugate_gradient(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:146
Index maxIterations() const
Definition: IterativeSolverBase.h:286
A conjugate gradient solver for sparse (or dense) least-square problems.
Definition: LeastSquareConjugateGradient.h:152
LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:168
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition: LeastSquareConjugateGradient.h:181
@ Success
Definition: Constants.h:444
@ NoConvergence
Definition: Constants.h:448
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Definition: EigenBase.h:32