10 #ifndef EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
11 #define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
28 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
31 const Preconditioner& precond,
Index& iters,
32 typename Dest::RealScalar& tol_error)
36 typedef typename Dest::RealScalar RealScalar;
37 typedef typename Dest::Scalar Scalar;
40 RealScalar tol = tol_error;
41 Index maxIters = iters;
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);
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;
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_> >
108 typedef Preconditioner_ Preconditioner;
150 template<
typename MatrixType_,
typename Preconditioner_>
180 template<
typename MatrixDerived>
186 template<
typename Rhs,
typename Dest>
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
RealReturnType real() const
#define EIGEN_DONT_INLINE
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
internal::traits< Derived >::Scalar Scalar
Base class for linear iterative solvers.
Index maxIterations() const
Preconditioner m_preconditioner
LeastSquaresConjugateGradient< MatrixType_, Preconditioner_ > & derived()
const ActualMatrixType & matrix() const
A conjugate gradient solver for sparse (or dense) least-square problems.
IterativeSolverBase< LeastSquaresConjugateGradient > Base
~LeastSquaresConjugateGradient()
MatrixType::RealScalar RealScalar
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Preconditioner_ Preconditioner
MatrixType::Scalar Scalar
LeastSquaresConjugateGradient()
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
const ActualMatrixType & matrix() const
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)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)