Medial Code Documentation
Loading...
Searching...
No Matches
ConjugateGradient.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2011-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_CONJUGATE_GRADIENT_H
11#define EIGEN_CONJUGATE_GRADIENT_H
12
13namespace Eigen {
14
15namespace internal {
16
26template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
27EIGEN_DONT_INLINE
28void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
29 const Preconditioner& precond, Index& iters,
30 typename Dest::RealScalar& tol_error)
31{
32 typedef typename Dest::RealScalar RealScalar;
33 typedef typename Dest::Scalar Scalar;
34 typedef Matrix<Scalar,Dynamic,1> VectorType;
35
36 RealScalar tol = tol_error;
37 Index maxIters = iters;
38
39 Index n = mat.cols();
40
41 VectorType residual = rhs - mat * x; //initial residual
42
43 RealScalar rhsNorm2 = rhs.squaredNorm();
44 if(rhsNorm2 == 0)
45 {
46 x.setZero();
47 iters = 0;
48 tol_error = 0;
49 return;
50 }
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)
55 {
56 iters = 0;
57 tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
58 return;
59 }
60
61 VectorType p(n);
62 p = precond.solve(residual); // initial search direction
63
64 VectorType z(n), tmp(n);
65 RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
66 Index i = 0;
67 while(i < maxIters)
68 {
69 tmp.noalias() = mat * p; // the bottleneck of the algorithm
70
71 Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
72 x += alpha * p; // update solution
73 residual -= alpha * tmp; // update residual
74
75 residualNorm2 = residual.squaredNorm();
76 if(residualNorm2 < threshold)
77 break;
78
79 z = precond.solve(residual); // approximately solve for "A z = residual"
80
81 RealScalar absOld = absNew;
82 absNew = numext::real(residual.dot(z)); // update the absolute value of r
83 RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
84 p = z + beta * p; // update search direction
85 i++;
86 }
87 tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
88 iters = i;
89}
90
91}
92
93template< typename _MatrixType, int _UpLo=Lower,
94 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
95class ConjugateGradient;
96
97namespace internal {
98
99template< typename _MatrixType, int _UpLo, typename _Preconditioner>
101{
102 typedef _MatrixType MatrixType;
104};
105
106}
107
155template< typename _MatrixType, int _UpLo, typename _Preconditioner>
156class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
157{
159 using Base::matrix;
160 using Base::m_error;
161 using Base::m_iterations;
162 using Base::m_info;
163 using Base::m_isInitialized;
164public:
165 typedef _MatrixType MatrixType;
166 typedef typename MatrixType::Scalar Scalar;
167 typedef typename MatrixType::RealScalar RealScalar;
169
170 enum {
171 UpLo = _UpLo
172 };
173
174public:
175
178
189 template<typename MatrixDerived>
190 explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
191
193
195 template<typename Rhs,typename Dest>
196 void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
197 {
198 typedef typename Base::MatrixWrapper MatrixWrapper;
199 typedef typename Base::ActualMatrixType ActualMatrixType;
200 enum {
201 TransposeInput = (!MatrixWrapper::MatrixFree)
202 && (UpLo==(Lower|Upper))
203 && (!MatrixType::IsRowMajor)
204 && (!NumTraits<Scalar>::IsComplex)
205 };
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),
209 RowMajorWrapper,
210 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
211 >::type SelfAdjointWrapper;
212
213 m_iterations = Base::maxIterations();
214 m_error = Base::m_tolerance;
215
216 RowMajorWrapper row_mat(matrix());
217 internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b, x, Base::m_preconditioner, m_iterations, m_error);
218 m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
219 }
220
221protected:
222
223};
224
225} // end namespace Eigen
226
227#endif // EIGEN_CONJUGATE_GRADIENT_H
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition ConjugateGradient.h:157
ConjugateGradient()
Default constructor.
Definition ConjugateGradient.h:177
ConjugateGradient(const EigenBase< MatrixDerived > &A)
Initialize the solver with matrix A for further Ax=b solving.
Definition ConjugateGradient.h:190
Base class for linear iterative solvers.
Definition IterativeSolverBase.h:144
Index maxIterations() const
Definition IterativeSolverBase.h:281
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
@ 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
Namespace containing all symbols from the Eigen library.
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Definition ForwardDeclarations.h:17
Definition inference.c:32