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 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
28 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
30 void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
31  const Preconditioner& precond, Index& iters,
32  typename Dest::RealScalar& tol_error)
33 {
34  typedef typename Dest::RealScalar RealScalar;
35  typedef typename Dest::Scalar Scalar;
36  typedef Matrix<Scalar,Dynamic,1> VectorType;
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  const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
54  RealScalar threshold = numext::maxi(RealScalar(tol*tol*rhsNorm2),considerAsZero);
55  RealScalar residualNorm2 = residual.squaredNorm();
56  if (residualNorm2 < threshold)
57  {
58  iters = 0;
59  tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
60  return;
61  }
62 
63  VectorType p(n);
64  p = precond.solve(residual); // initial search direction
65 
66  VectorType z(n), tmp(n);
67  RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
68  Index i = 0;
69  while(i < maxIters)
70  {
71  tmp.noalias() = mat * p; // the bottleneck of the algorithm
72 
73  Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
74  x += alpha * p; // update solution
75  residual -= alpha * tmp; // update residual
76 
77  residualNorm2 = residual.squaredNorm();
78  if(residualNorm2 < threshold)
79  break;
80 
81  z = precond.solve(residual); // approximately solve for "A z = residual"
82 
83  RealScalar absOld = absNew;
84  absNew = numext::real(residual.dot(z)); // update the absolute value of r
85  RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
86  p = z + beta * p; // update search direction
87  i++;
88  }
89  tol_error = numext::sqrt(residualNorm2 / rhsNorm2);
90  iters = i;
91 }
92 
93 }
94 
95 template< typename MatrixType_, int UpLo_=Lower,
96  typename Preconditioner_ = DiagonalPreconditioner<typename MatrixType_::Scalar> >
97 class ConjugateGradient;
98 
99 namespace internal {
100 
101 template< typename MatrixType_, int UpLo_, typename Preconditioner_>
102 struct traits<ConjugateGradient<MatrixType_,UpLo_,Preconditioner_> >
103 {
104  typedef MatrixType_ MatrixType;
105  typedef Preconditioner_ Preconditioner;
106 };
107 
108 }
109 
157 template< typename MatrixType_, int UpLo_, typename Preconditioner_>
158 class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<MatrixType_,UpLo_,Preconditioner_> >
159 {
161  using Base::matrix;
162  using Base::m_error;
163  using Base::m_iterations;
164  using Base::m_info;
165  using Base::m_isInitialized;
166 public:
167  typedef MatrixType_ MatrixType;
168  typedef typename MatrixType::Scalar Scalar;
170  typedef Preconditioner_ Preconditioner;
171 
172  enum {
173  UpLo = UpLo_
174  };
175 
176 public:
177 
180 
191  template<typename MatrixDerived>
193 
195 
197  template<typename Rhs,typename Dest>
198  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
199  {
200  typedef typename Base::MatrixWrapper MatrixWrapper;
201  typedef typename Base::ActualMatrixType ActualMatrixType;
202  enum {
203  TransposeInput = (!MatrixWrapper::MatrixFree)
204  && (UpLo==(Lower|Upper))
205  && (!MatrixType::IsRowMajor)
207  };
208  typedef std::conditional_t<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType const&> RowMajorWrapper;
209  EIGEN_STATIC_ASSERT(internal::check_implication(MatrixWrapper::MatrixFree,UpLo==(Lower|Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
210  typedef std::conditional_t<UpLo==(Lower|Upper),
211  RowMajorWrapper,
212  typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
213  > SelfAdjointWrapper;
214 
217 
218  RowMajorWrapper row_mat(matrix());
221  }
222 
223 protected:
224 
225 };
226 
227 } // end namespace Eigen
228 
229 #endif // EIGEN_CONJUGATE_GRADIENT_H
Array< int, 3, 1 > b
int n
RealReturnType real() const
MatrixXcf A
#define EIGEN_DONT_INLINE
Definition: Macros.h:844
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
float * p
Matrix< float, 1, Dynamic > MatrixType
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
IterativeSolverBase< ConjugateGradient > Base
MatrixType::Scalar Scalar
MatrixType::RealScalar RealScalar
ConjugateGradient(const EigenBase< MatrixDerived > &A)
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
const ActualMatrixType & matrix() const
Preconditioner_ Preconditioner
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
Base class for linear iterative solvers.
internal::generic_matrix_wrapper< MatrixType > MatrixWrapper
MatrixWrapper::ActualMatrixType ActualMatrixType
ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ > & derived()
const ActualMatrixType & matrix() const
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
constexpr bool check_implication(bool a, bool b)
Definition: Meta.h:579
EIGEN_DONT_INLINE void conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
EIGEN_ALWAYS_INLINE float sqrt(const float &x)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231