10 #ifndef EIGEN_CONJUGATE_GRADIENT_H
11 #define EIGEN_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)
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;
45 RealScalar rhsNorm2 = rhs.squaredNorm();
54 RealScalar threshold =
numext::maxi(RealScalar(tol*tol*rhsNorm2),considerAsZero);
55 RealScalar residualNorm2 = residual.squaredNorm();
56 if (residualNorm2 < threshold)
64 p = precond.solve(residual);
66 VectorType z(
n), tmp(
n);
71 tmp.noalias() =
mat *
p;
73 Scalar alpha = absNew /
p.dot(tmp);
75 residual -= alpha * tmp;
77 residualNorm2 = residual.squaredNorm();
78 if(residualNorm2 < threshold)
81 z = precond.solve(residual);
83 RealScalar absOld = absNew;
85 RealScalar beta = absNew / absOld;
95 template<
typename MatrixType_,
int UpLo_=
Lower,
96 typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
97 class ConjugateGradient;
101 template<
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
102 struct traits<ConjugateGradient<MatrixType_,UpLo_,Preconditioner_> >
105 typedef Preconditioner_ Preconditioner;
157 template<
typename MatrixType_,
int UpLo_,
typename Preconditioner_>
191 template<
typename MatrixDerived>
197 template<
typename Rhs,
typename Dest>
203 TransposeInput = (!MatrixWrapper::MatrixFree)
205 && (!MatrixType::IsRowMajor)
208 typedef std::conditional_t<TransposeInput,Transpose<const ActualMatrixType>,
ActualMatrixType const&> RowMajorWrapper;
213 > SelfAdjointWrapper;
218 RowMajorWrapper row_mat(
matrix());
RealReturnType real() const
#define EIGEN_DONT_INLINE
#define EIGEN_STATIC_ASSERT(X, MSG)
Matrix< float, 1, Dynamic > MatrixType
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
IterativeSolverBase< ConjugateGradient > Base
MatrixType::Scalar Scalar
MatrixType::RealScalar RealScalar
ConjugateGradient(const EigenBase< MatrixDerived > &A)
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
const ActualMatrixType & matrix() const
Preconditioner_ Preconditioner
NumTraits< Scalar >::Real RealScalar
internal::traits< Derived >::Scalar Scalar
Base class for linear iterative solvers.
internal::generic_matrix_wrapper< MatrixType > MatrixWrapper
Index maxIterations() const
MatrixWrapper::ActualMatrixType ActualMatrixType
Preconditioner m_preconditioner
ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ > & derived()
const ActualMatrixType & matrix() const
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
constexpr bool check_implication(bool a, bool b)
EIGEN_DONT_INLINE void conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
EIGEN_ALWAYS_INLINE float sqrt(const float &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.