10#ifndef EIGEN_CONJUGATE_GRADIENT_H
11#define EIGEN_CONJUGATE_GRADIENT_H
26template<
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();
54 RealScalar threshold =
numext::maxi(RealScalar(tol*tol*rhsNorm2),considerAsZero);
55 RealScalar residualNorm2 = residual.squaredNorm();
56 if (residualNorm2 < threshold)
59 tol_error =
sqrt(residualNorm2 / rhsNorm2);
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;
89 tol_error =
sqrt(residualNorm2 / rhsNorm2);
95template<
typename _MatrixType,
int _UpLo=
Lower,
96 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
97class ConjugateGradient;
101template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
157template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
168 typedef typename MatrixType::Scalar
Scalar;
191 template<
typename MatrixDerived>
197 template<
typename Rhs,
typename Dest>
203 TransposeInput = (!MatrixWrapper::MatrixFree)
205 && (!MatrixType::IsRowMajor)
213 >
::type SelfAdjointWrapper;
218 RowMajorWrapper row_mat(matrix());
EIGEN_DEVICE_FUNC RealReturnType real() const
Definition: CommonCwiseUnaryOps.h:100
#define EIGEN_DONT_INLINE
Definition: Macros.h:950
#define EIGEN_IMPLIES(a, b)
Definition: Macros.h:1325
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:127
constexpr common_return_t< T1, T2 > beta(const T1 a, const T2 b) noexcept
Compile-time beta function.
Definition: beta.hpp:36
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition: ConjugateGradient.h:159
_Preconditioner Preconditioner
Definition: ConjugateGradient.h:170
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: ConjugateGradient.h:198
MatrixType::Scalar Scalar
Definition: ConjugateGradient.h:168
ConjugateGradient()
Default constructor.
Definition: ConjugateGradient.h:179
~ConjugateGradient()
Definition: ConjugateGradient.h:194
MatrixType::RealScalar RealScalar
Definition: ConjugateGradient.h:169
ConjugateGradient(const EigenBase< MatrixDerived > &A)
Initialize the solver with matrix A for further Ax=b solving.
Definition: ConjugateGradient.h:192
_MatrixType MatrixType
Definition: ConjugateGradient.h:167
@ UpLo
Definition: ConjugateGradient.h:173
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:144
internal::generic_matrix_wrapper< MatrixType > MatrixWrapper
Definition: IterativeSolverBase.h:416
Index maxIterations() const
Definition: IterativeSolverBase.h:281
ComputationInfo m_info
Definition: IterativeSolverBase.h:438
MatrixWrapper::ActualMatrixType ActualMatrixType
Definition: IterativeSolverBase.h:417
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
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived()
Definition: SparseSolverBase.h:79
RealScalar m_tolerance
Definition: IterativeSolverBase.h:434
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:419
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:145
Definition: IterativeSolverBase.h:47
type
Definition: core.h:575
UnitType abs(const UnitType x) noexcept
Compute absolute value.
Definition: math.h:721
auto sqrt(const UnitType &value) noexcept -> unit_t< square_root< typename units::traits::unit_t_traits< UnitType >::unit_type >, typename units::traits::unit_t_traits< UnitType >::underlying_type, linear_scale >
computes the square root of value
Definition: math.h:483
@ Lower
View matrix as a lower triangular matrix.
Definition: Constants.h:209
@ Upper
View matrix as an upper triangular matrix.
Definition: Constants.h:211
@ Success
Computation was successful.
Definition: Constants.h:442
@ NoConvergence
Iterative procedure did not converge.
Definition: Constants.h:446
constexpr common_t< T1, T2 > min(const T1 x, const T2 y) noexcept
Compile-time pairwise minimum function.
Definition: min.hpp:35
Type
Definition: Constants.h:471
EIGEN_DONT_INLINE void conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
Definition: ConjugateGradient.h:28
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:1091
Namespace containing all symbols from the Eigen library.
Definition: MatrixExponential.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Definition: Eigen_Colamd.h:50
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:30
Holds information about the various numeric (i.e.
Definition: NumTraits.h:233
_MatrixType MatrixType
Definition: ConjugateGradient.h:104
_Preconditioner Preconditioner
Definition: ConjugateGradient.h:105
Definition: ForwardDeclarations.h:17