ConjugateGradient.h
Go to the documentation of this file.
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 
13 namespace Eigen {
14 
15 namespace internal {
16 
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)
31 {
32  using std::sqrt;
33  using std::abs;
34  typedef typename Dest::RealScalar RealScalar;
35  typedef typename Dest::Scalar Scalar;
37 
38  RealScalar tol = tol_error;
39  Index maxIters = iters;
40 
41  Index n = mat.cols();
42 
43  VectorType residual = rhs - mat * x; //initial residual
44 
45  RealScalar rhsNorm2 = rhs.squaredNorm();
46  if(rhsNorm2 == 0)
47  {
48  x.setZero();
49  iters = 0;
50  tol_error = 0;
51  return;
52  }
53  RealScalar threshold = tol*tol*rhsNorm2;
54  RealScalar residualNorm2 = residual.squaredNorm();
55  if (residualNorm2 < threshold)
56  {
57  iters = 0;
58  tol_error = sqrt(residualNorm2 / rhsNorm2);
59  return;
60  }
61 
62  VectorType p(n);
63  p = precond.solve(residual); // initial search direction
64 
65  VectorType z(n), tmp(n);
66  RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
67  Index i = 0;
68  while(i < maxIters)
69  {
70  tmp.noalias() = mat * p; // the bottleneck of the algorithm
71 
72  Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
73  x += alpha * p; // update solution
74  residual -= alpha * tmp; // update residual
75 
76  residualNorm2 = residual.squaredNorm();
77  if(residualNorm2 < threshold)
78  break;
79 
80  z = precond.solve(residual); // approximately solve for "A z = residual"
81 
82  RealScalar absOld = absNew;
83  absNew = numext::real(residual.dot(z)); // update the absolute value of r
84  RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
85  p = z + beta * p; // update search direction
86  i++;
87  }
88  tol_error = sqrt(residualNorm2 / rhsNorm2);
89  iters = i;
90 }
91 
92 }
93 
94 template< typename _MatrixType, int _UpLo=Lower,
95  typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
97 
98 namespace internal {
99 
100 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
101 struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
102 {
103  typedef _MatrixType MatrixType;
104  typedef _Preconditioner Preconditioner;
105 };
106 
107 }
108 
156 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
157 class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
158 {
160  using Base::matrix;
161  using Base::m_error;
162  using Base::m_iterations;
163  using Base::m_info;
164  using Base::m_isInitialized;
165 public:
166  typedef _MatrixType MatrixType;
167  typedef typename MatrixType::Scalar Scalar;
168  typedef typename MatrixType::RealScalar RealScalar;
169  typedef _Preconditioner Preconditioner;
170 
171  enum {
172  UpLo = _UpLo
173  };
174 
175 public:
176 
178  ConjugateGradient() : Base() {}
179 
190  template<typename MatrixDerived>
191  explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
192 
194 
196  template<typename Rhs,typename Dest>
197  void _solve_with_guess_impl(const Rhs& b, Dest& x) const
198  {
199  typedef typename Base::MatrixWrapper MatrixWrapper;
200  typedef typename Base::ActualMatrixType ActualMatrixType;
201  enum {
202  TransposeInput = (!MatrixWrapper::MatrixFree)
203  && (UpLo==(Lower|Upper))
204  && (!MatrixType::IsRowMajor)
206  };
207  typedef typename internal::conditional<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType const&>::type RowMajorWrapper;
208  EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(MatrixWrapper::MatrixFree,UpLo==(Lower|Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
209  typedef typename internal::conditional<UpLo==(Lower|Upper),
210  RowMajorWrapper,
211  typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
212  >::type SelfAdjointWrapper;
213  m_iterations = Base::maxIterations();
214  m_error = Base::m_tolerance;
215 
216  for(Index j=0; j<b.cols(); ++j)
217  {
218  m_iterations = Base::maxIterations();
219  m_error = Base::m_tolerance;
220 
221  typename Dest::ColXpr xj(x,j);
222  RowMajorWrapper row_mat(matrix());
223  internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
224  }
225 
226  m_isInitialized = true;
227  m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
228  }
229 
231  using Base::_solve_impl;
232  template<typename Rhs,typename Dest>
233  void _solve_impl(const MatrixBase<Rhs>& b, Dest& x) const
234  {
235  x.setZero();
236  _solve_with_guess_impl(b.derived(),x);
237  }
238 
239 protected:
240 
241 };
242 
243 } // end namespace Eigen
244 
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
Definition: LDLT.h:16
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ColXpr
Definition: BlockMethods.h:14
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:122
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
#define EIGEN_IMPLIES(a, b)
Definition: Macros.h:899
#define EIGEN_DONT_INLINE
Definition: Macros.h:515
MatrixWrapper::ActualMatrixType ActualMatrixType
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
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.
Definition: MatrixBase.h:48
_Preconditioner Preconditioner


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:08:04