10 #ifndef EIGEN_CONJUGATE_GRADIENT_H 11 #define EIGEN_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;
45 RealScalar rhsNorm2 = rhs.squaredNorm();
53 RealScalar threshold = tol*tol*rhsNorm2;
54 RealScalar residualNorm2 = residual.squaredNorm();
55 if (residualNorm2 < threshold)
58 tol_error =
sqrt(residualNorm2 / rhsNorm2);
63 p = precond.solve(residual);
65 VectorType
z(n), tmp(n);
70 tmp.noalias() = mat * p;
72 Scalar alpha = absNew / p.dot(tmp);
74 residual -= alpha * tmp;
76 residualNorm2 = residual.squaredNorm();
77 if(residualNorm2 < threshold)
80 z = precond.solve(residual);
82 RealScalar absOld = absNew;
84 RealScalar beta = absNew / absOld;
88 tol_error =
sqrt(residualNorm2 / rhsNorm2);
94 template<
typename _MatrixType,
int _UpLo=
Lower,
100 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
156 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
162 using Base::m_iterations;
164 using Base::m_isInitialized;
167 typedef typename MatrixType::Scalar
Scalar;
190 template<
typename MatrixDerived>
196 template<
typename Rhs,
typename Dest>
202 TransposeInput = (!MatrixWrapper::MatrixFree)
204 && (!MatrixType::IsRowMajor)
212 >::type SelfAdjointWrapper;
213 m_iterations = Base::maxIterations();
214 m_error = Base::m_tolerance;
216 for(
Index j=0; j<b.cols(); ++j)
218 m_iterations = Base::maxIterations();
219 m_error = Base::m_tolerance;
222 RowMajorWrapper row_mat(matrix());
226 m_isInitialized =
true;
231 using Base::_solve_impl;
232 template<
typename Rhs,
typename Dest>
236 _solve_with_guess_impl(b.derived(),x);
245 #endif // EIGEN_CONJUGATE_GRADIENT_H A preconditioner based on the digonal entries.
EIGEN_DEVICE_FUNC RealReturnType real() const
void _solve_with_guess_impl(const Rhs &b, Dest &x) const
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ColXpr
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
#define EIGEN_IMPLIES(a, b)
#define EIGEN_DONT_INLINE
_Preconditioner Preconditioner
MatrixWrapper::ActualMatrixType ActualMatrixType
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
TFSIMD_FORCE_INLINE const tfScalar & z() const
EIGEN_DONT_INLINE void conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
MatrixType::RealScalar RealScalar
void _solve_impl(const MatrixBase< Rhs > &b, Dest &x) const
MatrixType::Scalar Scalar
IterativeSolverBase< ConjugateGradient > Base
ConjugateGradient(const EigenBase< MatrixDerived > &A)
Base class for linear iterative solvers.
EIGEN_DEVICE_FUNC const Scalar & b
Base class for all dense matrices, vectors, and expressions.
_Preconditioner Preconditioner