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,
97 typename _Preconditioner = LeastSquareDiagonalPreconditioner<typename _MatrixType::Scalar> >
98 class LeastSquaresConjugateGradient;
102 template<
typename _MatrixType,
typename _Preconditioner>
148 template<
typename _MatrixType,
typename _Preconditioner>
159 typedef typename MatrixType::Scalar
Scalar;
178 template<
typename MatrixDerived>
184 template<
typename Rhs,
typename Dest>
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
Definition: ArrayCwiseUnaryOps.h:52
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition: ArrayCwiseUnaryOps.h:187
EIGEN_DEVICE_FUNC RealReturnType real() const
Definition: CommonCwiseUnaryOps.h:100
#define EIGEN_DONT_INLINE
Definition: Macros.h:940
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:144
Index maxIterations() const
Definition: IterativeSolverBase.h:281
ComputationInfo m_info
Definition: IterativeSolverBase.h:438
RealScalar m_error
Definition: IterativeSolverBase.h:436
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:431
Index m_iterations
Definition: IterativeSolverBase.h:437
bool m_isInitialized
Definition: SparseSolverBase.h:119
RealScalar m_tolerance
Definition: IterativeSolverBase.h:434
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & derived()
Definition: SparseSolverBase.h:79
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:419
A conjugate gradient solver for sparse (or dense) least-square problems.
Definition: LeastSquareConjugateGradient.h:150
ComputationInfo m_info
Definition: IterativeSolverBase.h:438
~LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:181
_MatrixType MatrixType
Definition: LeastSquareConjugateGradient.h:158
MatrixType::RealScalar RealScalar
Definition: LeastSquareConjugateGradient.h:160
RealScalar m_error
Definition: IterativeSolverBase.h:436
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition: LeastSquareConjugateGradient.h:179
MatrixType::Scalar Scalar
Definition: LeastSquareConjugateGradient.h:159
Index m_iterations
Definition: IterativeSolverBase.h:437
LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:166
_Preconditioner Preconditioner
Definition: LeastSquareConjugateGradient.h:161
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: LeastSquareConjugateGradient.h:185
IterativeSolverBase< LeastSquaresConjugateGradient > Base
Definition: LeastSquareConjugateGradient.h:151
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:419
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:143
@ Success
Definition: Constants.h:442
@ NoConvergence
Definition: Constants.h:446
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)
Definition: LeastSquareConjugateGradient.h:28
Namespace containing all symbols from the Eigen library.
Definition: LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Definition: document.h:416
Definition: EigenBase.h:30
_Preconditioner Preconditioner
Definition: LeastSquareConjugateGradient.h:106
_MatrixType MatrixType
Definition: LeastSquareConjugateGradient.h:105
Definition: ForwardDeclarations.h:17