10#ifndef EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
11#define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
14#include "./InternalHeaderCheck.h"
29template <
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
30EIGEN_DONT_INLINE
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) {
35 typedef typename Dest::RealScalar RealScalar;
36 typedef typename Dest::Scalar Scalar;
37 typedef Matrix<Scalar, Dynamic, 1> VectorType;
39 RealScalar tol = tol_error;
40 Index maxIters = iters;
42 Index m = mat.rows(), n = mat.cols();
44 VectorType residual = rhs - mat * x;
45 VectorType normal_residual = mat.adjoint() * residual;
47 RealScalar rhsNorm2 = (mat.adjoint() * rhs).squaredNorm();
54 RealScalar threshold = tol * tol * rhsNorm2;
55 RealScalar residualNorm2 = normal_residual.squaredNorm();
56 if (residualNorm2 < threshold) {
58 tol_error =
sqrt(residualNorm2 / rhsNorm2);
63 p = precond.solve(normal_residual);
65 VectorType z(n), tmp(m);
66 RealScalar absNew = numext::real(normal_residual.dot(p));
68 while (i < maxIters) {
69 tmp.noalias() = mat * p;
71 Scalar alpha = absNew / tmp.squaredNorm();
73 residual -= alpha * tmp;
74 normal_residual.noalias() = mat.adjoint() * residual;
76 residualNorm2 = normal_residual.squaredNorm();
77 if (residualNorm2 < threshold)
break;
79 z = precond.solve(normal_residual);
81 RealScalar absOld = absNew;
82 absNew = numext::real(normal_residual.dot(z));
83 RealScalar beta = absNew / absOld;
87 tol_error =
sqrt(residualNorm2 / rhsNorm2);
93template <
typename MatrixType_,
94 typename Preconditioner_ = LeastSquareDiagonalPreconditioner<typename MatrixType_::Scalar> >
95class LeastSquaresConjugateGradient;
99template <
typename MatrixType_,
typename Preconditioner_>
100struct traits<LeastSquaresConjugateGradient<MatrixType_, Preconditioner_> > {
101 typedef MatrixType_ MatrixType;
102 typedef Preconditioner_ Preconditioner;
145template <
typename MatrixType_,
typename Preconditioner_>
147 :
public IterativeSolverBase<LeastSquaresConjugateGradient<MatrixType_, Preconditioner_> > {
151 using Base::m_isInitialized;
152 using Base::m_iterations;
156 typedef MatrixType_ MatrixType;
157 typedef typename MatrixType::Scalar Scalar;
158 typedef typename MatrixType::RealScalar RealScalar;
159 typedef Preconditioner_ Preconditioner;
175 template <
typename MatrixDerived>
181 template <
typename Rhs,
typename Dest>
182 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const {
184 m_error = Base::m_tolerance;
186 internal::least_square_conjugate_gradient(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
Base class for linear iterative solvers.
Definition IterativeSolverBase.h:124
Index maxIterations() const
Definition IterativeSolverBase.h:251
A conjugate gradient solver for sparse (or dense) least-square problems.
Definition LeastSquareConjugateGradient.h:147
LeastSquaresConjugateGradient()
Definition LeastSquareConjugateGradient.h:163
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition LeastSquareConjugateGradient.h:176
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Definition EigenBase.h:33