10#ifndef EIGEN_CONJUGATE_GRADIENT_H
11#define EIGEN_CONJUGATE_GRADIENT_H
14#include "./InternalHeaderCheck.h"
29template <
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
30EIGEN_DONT_INLINE
void conjugate_gradient(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
const Preconditioner& precond,
31 Index& iters,
typename Dest::RealScalar& tol_error) {
32 typedef typename Dest::RealScalar RealScalar;
33 typedef typename Dest::Scalar Scalar;
34 typedef Matrix<Scalar, Dynamic, 1> VectorType;
36 RealScalar tol = tol_error;
37 Index maxIters = iters;
41 VectorType residual = rhs - mat * x;
43 RealScalar rhsNorm2 = rhs.squaredNorm();
50 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
51 RealScalar threshold = numext::maxi(RealScalar(tol * tol * rhsNorm2), considerAsZero);
52 RealScalar residualNorm2 = residual.squaredNorm();
53 if (residualNorm2 < threshold) {
55 tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
60 p = precond.solve(residual);
62 VectorType z(n), tmp(n);
63 RealScalar absNew = numext::real(residual.dot(p));
65 while (i < maxIters) {
66 tmp.noalias() = mat * p;
68 Scalar alpha = absNew / p.dot(tmp);
70 residual -= alpha * tmp;
72 residualNorm2 = residual.squaredNorm();
73 if (residualNorm2 < threshold)
break;
75 z = precond.solve(residual);
77 RealScalar absOld = absNew;
78 absNew = numext::real(residual.dot(z));
79 RealScalar beta = absNew / absOld;
83 tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
89template <
typename MatrixType_,
int UpLo_ =
Lower,
90 typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
91class ConjugateGradient;
95template <
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
96struct traits<ConjugateGradient<MatrixType_, UpLo_, Preconditioner_> > {
97 typedef MatrixType_ MatrixType;
98 typedef Preconditioner_ Preconditioner;
151template <
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
156 using Base::m_isInitialized;
157 using Base::m_iterations;
161 typedef MatrixType_ MatrixType;
162 typedef typename MatrixType::Scalar Scalar;
163 typedef typename MatrixType::RealScalar RealScalar;
164 typedef Preconditioner_ Preconditioner;
166 enum { UpLo = UpLo_ };
182 template <
typename MatrixDerived>
188 template <
typename Rhs,
typename Dest>
189 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const {
190 typedef typename Base::MatrixWrapper MatrixWrapper;
191 typedef typename Base::ActualMatrixType ActualMatrixType;
193 TransposeInput = (!MatrixWrapper::MatrixFree) && (UpLo == (
Lower |
Upper)) && (!MatrixType::IsRowMajor) &&
194 (!NumTraits<Scalar>::IsComplex)
196 typedef std::conditional_t<TransposeInput, Transpose<const ActualMatrixType>, ActualMatrixType
const&>
198 EIGEN_STATIC_ASSERT(internal::check_implication(MatrixWrapper::MatrixFree, UpLo == (
Lower |
Upper)),
199 MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
200 typedef std::conditional_t<UpLo == (
Lower |
Upper), RowMajorWrapper,
201 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type>
205 m_error = Base::m_tolerance;
207 RowMajorWrapper row_mat(matrix());
208 internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b, x, Base::m_preconditioner, m_iterations, m_error);
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition ConjugateGradient.h:152
ConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition ConjugateGradient.h:183
ConjugateGradient()
Definition ConjugateGradient.h:170
Base class for linear iterative solvers.
Definition IterativeSolverBase.h:124
Index maxIterations() const
Definition IterativeSolverBase.h:251
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Namespace containing all symbols from the Eigen library.
Definition Core:137
Definition EigenBase.h:33