10 #ifndef EIGEN_CONJUGATE_GRADIENT_H
11 #define EIGEN_CONJUGATE_GRADIENT_H
26 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
28 void conjugate_gradient(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
29 const Preconditioner& precond, Index& iters,
30 typename Dest::RealScalar& tol_error)
34 typedef typename Dest::RealScalar RealScalar;
35 typedef typename Dest::Scalar Scalar;
36 typedef Matrix<Scalar,Dynamic,1> VectorType;
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);
66 RealScalar absNew = numext::real(residual.dot(p));
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;
83 absNew = numext::real(residual.dot(z));
84 RealScalar beta = absNew / absOld;
88 tol_error = sqrt(residualNorm2 / rhsNorm2);
94 template<
typename _MatrixType,
int _UpLo=
Lower,
95 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
100 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
103 typedef _MatrixType MatrixType;
104 typedef _Preconditioner Preconditioner;
152 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
153 class ConjugateGradient :
public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
155 typedef IterativeSolverBase<ConjugateGradient> Base;
156 using Base::mp_matrix;
158 using Base::m_iterations;
160 using Base::m_isInitialized;
162 typedef _MatrixType MatrixType;
163 typedef typename MatrixType::Scalar Scalar;
164 typedef typename MatrixType::RealScalar RealScalar;
165 typedef _Preconditioner Preconditioner;
191 template<
typename Rhs,
typename Dest>
192 void _solve_with_guess_impl(
const Rhs& b, Dest& x)
const
194 typedef Ref<const MatrixType> MatRef;
195 typedef typename internal::conditional<UpLo==(Lower|Upper) && (!MatrixType::IsRowMajor) && (!NumTraits<Scalar>::IsComplex),
196 Transpose<const MatRef>, MatRef
const&>::type RowMajorWrapper;
197 typedef typename internal::conditional<UpLo==(
Lower|
Upper),
199 typename MatRef::template ConstSelfAdjointViewReturnType<UpLo>::Type
200 >::type SelfAdjointWrapper;
202 m_error = Base::m_tolerance;
204 for(Index j=0; j<b.cols(); ++j)
207 m_error = Base::m_tolerance;
209 typename Dest::ColXpr xj(x,j);
210 RowMajorWrapper row_mat(mp_matrix);
211 internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
214 m_isInitialized =
true;
219 using Base::_solve_impl;
220 template<
typename Rhs,
typename Dest>
221 void _solve_impl(
const MatrixBase<Rhs>& b, Dest& x)
const
224 _solve_with_guess_impl(b.derived(),x);
233 #endif // EIGEN_CONJUGATE_GRADIENT_H
ConjugateGradient(const MatrixType &A)
Definition: ConjugateGradient.h:186
Definition: Constants.h:196
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition: ConjugateGradient.h:96
ConjugateGradient()
Definition: ConjugateGradient.h:174
Index maxIterations() const
Definition: IterativeSolverBase.h:155
Definition: Constants.h:198
Definition: Constants.h:424
Definition: Eigen_Colamd.h:54
Definition: Constants.h:428