Eigen  3.2.7
ConjugateGradient.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011 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>
27 EIGEN_DONT_INLINE
28 void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
29  const Preconditioner& precond, int& 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;
36  typedef Matrix<Scalar,Dynamic,1> VectorType;
37 
38  RealScalar tol = tol_error;
39  int maxIters = iters;
40 
41  int 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  int 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 residue
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 
144 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
145 class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
146 {
147  typedef IterativeSolverBase<ConjugateGradient> Base;
148  using Base::mp_matrix;
149  using Base::m_error;
150  using Base::m_iterations;
151  using Base::m_info;
152  using Base::m_isInitialized;
153 public:
154  typedef _MatrixType MatrixType;
155  typedef typename MatrixType::Scalar Scalar;
156  typedef typename MatrixType::Index Index;
157  typedef typename MatrixType::RealScalar RealScalar;
158  typedef _Preconditioner Preconditioner;
159 
160  enum {
161  UpLo = _UpLo
162  };
163 
164 public:
165 
167  ConjugateGradient() : Base() {}
168 
179  template<typename MatrixDerived>
180  explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
181 
182  ~ConjugateGradient() {}
183 
189  template<typename Rhs,typename Guess>
190  inline const internal::solve_retval_with_guess<ConjugateGradient, Rhs, Guess>
191  solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
192  {
193  eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
194  eigen_assert(Base::rows()==b.rows()
195  && "ConjugateGradient::solve(): invalid number of rows of the right hand side matrix b");
196  return internal::solve_retval_with_guess
197  <ConjugateGradient, Rhs, Guess>(*this, b.derived(), x0);
198  }
199 
201  template<typename Rhs,typename Dest>
202  void _solveWithGuess(const Rhs& b, Dest& x) const
203  {
204  typedef typename internal::conditional<UpLo==(Lower|Upper),
205  const MatrixType&,
207  >::type MatrixWrapperType;
208  m_iterations = Base::maxIterations();
209  m_error = Base::m_tolerance;
210 
211  for(int j=0; j<b.cols(); ++j)
212  {
213  m_iterations = Base::maxIterations();
214  m_error = Base::m_tolerance;
215 
216  typename Dest::ColXpr xj(x,j);
217  internal::conjugate_gradient(MatrixWrapperType(*mp_matrix), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
218  }
219 
220  m_isInitialized = true;
221  m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
222  }
223 
225  template<typename Rhs,typename Dest>
226  void _solve(const Rhs& b, Dest& x) const
227  {
228  x.setZero();
229  _solveWithGuess(b,x);
230  }
231 
232 protected:
233 
234 };
235 
236 
237 namespace internal {
238 
239 template<typename _MatrixType, int _UpLo, typename _Preconditioner, typename Rhs>
240 struct solve_retval<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner>, Rhs>
241  : solve_retval_base<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner>, Rhs>
242 {
243  typedef ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> Dec;
244  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
245 
246  template<typename Dest> void evalTo(Dest& dst) const
247  {
248  dec()._solve(rhs(),dst);
249  }
250 };
251 
252 } // end namespace internal
253 
254 } // end namespace Eigen
255 
256 #endif // EIGEN_CONJUGATE_GRADIENT_H
ConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition: ConjugateGradient.h:180
Definition: Constants.h:167
ConjugateGradient()
Definition: ConjugateGradient.h:167
const internal::solve_retval_with_guess< ConjugateGradient, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: ConjugateGradient.h:191
Definition: LDLT.h:16
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
Definition: SparseSelfAdjointView.h:49
A conjugate gradient solver for sparse self-adjoint problems.
Definition: ConjugateGradient.h:96
Definition: Constants.h:169
Definition: EigenBase.h:26
Definition: Eigen_Colamd.h:54
Definition: Constants.h:380
Definition: Constants.h:376
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
int maxIterations() const
Definition: IterativeSolverBase.h:141