10#ifndef EIGEN_CONJUGATE_GRADIENT_H
11#define EIGEN_CONJUGATE_GRADIENT_H
26template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
28void conjugate_gradient(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
29 const Preconditioner& precond,
Index& iters,
30 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();
51 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
52 RealScalar threshold = numext::maxi(RealScalar(tol*tol*rhsNorm2),considerAsZero);
53 RealScalar residualNorm2 = residual.squaredNorm();
54 if (residualNorm2 < threshold)
57 tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
62 p = precond.solve(residual);
64 VectorType z(n), tmp(n);
65 RealScalar absNew = numext::real(residual.dot(p));
69 tmp.noalias() = mat * p;
71 Scalar alpha = absNew / p.dot(tmp);
73 residual -= alpha * tmp;
75 residualNorm2 = residual.squaredNorm();
76 if(residualNorm2 < threshold)
79 z = precond.solve(residual);
81 RealScalar absOld = absNew;
82 absNew = numext::real(residual.dot(z));
83 RealScalar beta = absNew / absOld;
87 tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
93template<
typename _MatrixType,
int _UpLo=
Lower,
99template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
100struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
102 typedef _MatrixType MatrixType;
103 typedef _Preconditioner Preconditioner;
155template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
161 using Base::m_iterations;
163 using Base::m_isInitialized;
165 typedef _MatrixType MatrixType;
166 typedef typename MatrixType::Scalar Scalar;
167 typedef typename MatrixType::RealScalar RealScalar;
168 typedef _Preconditioner Preconditioner;
189 template<
typename MatrixDerived>
195 template<
typename Rhs,
typename Dest>
196 void _solve_vector_with_guess_impl(
const Rhs& b, Dest& x)
const
198 typedef typename Base::MatrixWrapper MatrixWrapper;
199 typedef typename Base::ActualMatrixType ActualMatrixType;
201 TransposeInput = (!MatrixWrapper::MatrixFree)
203 && (!MatrixType::IsRowMajor)
204 && (!NumTraits<Scalar>::IsComplex)
206 typedef typename internal::conditional<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType
const&>::type RowMajorWrapper;
207 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(MatrixWrapper::MatrixFree,UpLo==(
Lower|
Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
208 typedef typename internal::conditional<UpLo==(
Lower|
Upper),
210 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
211 >::type SelfAdjointWrapper;
214 m_error = Base::m_tolerance;
216 RowMajorWrapper row_mat(matrix());
217 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:157
ConjugateGradient()
Definition ConjugateGradient.h:177
ConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition ConjugateGradient.h:190
A preconditioner based on the digonal entries.
Definition BasicPreconditioners.h:37
IterativeSolverBase()
Definition IterativeSolverBase.h:166
Index maxIterations() const
Definition IterativeSolverBase.h:281
@ Lower
Definition Constants.h:209
@ Upper
Definition Constants.h:211
@ Success
Definition Constants.h:442
@ NoConvergence
Definition Constants.h:446
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Definition EigenBase.h:30