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 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, 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 
157 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
158 class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
159 {
161  using Base::mp_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;
169  typedef typename MatrixType::Index Index;
170  typedef typename MatrixType::RealScalar RealScalar;
171  typedef _Preconditioner Preconditioner;
172 
173  enum {
174  UpLo = _UpLo
175  };
176 
177 public:
178 
180  ConjugateGradient() : Base() {}
181 
192  ConjugateGradient(const MatrixType& A) : Base(A) {}
193 
195 
201  template<typename Rhs,typename Guess>
203  solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
204  {
205  eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
206  eigen_assert(Base::rows()==b.rows()
207  && "ConjugateGradient::solve(): invalid number of rows of the right hand side matrix b");
209  <ConjugateGradient, Rhs, Guess>(*this, b.derived(), x0);
210  }
211 
213  template<typename Rhs,typename Dest>
214  void _solveWithGuess(const Rhs& b, Dest& x) const
215  {
216  m_iterations = Base::maxIterations();
217  m_error = Base::m_tolerance;
218 
219  for(int j=0; j<b.cols(); ++j)
220  {
221  m_iterations = Base::maxIterations();
222  m_error = Base::m_tolerance;
223 
224  typename Dest::ColXpr xj(x,j);
225  internal::conjugate_gradient(mp_matrix->template selfadjointView<UpLo>(), b.col(j), xj,
226  Base::m_preconditioner, m_iterations, m_error);
227  }
228 
229  m_isInitialized = true;
230  m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
231  }
232 
234  template<typename Rhs,typename Dest>
235  void _solve(const Rhs& b, Dest& x) const
236  {
237  x.setOnes();
238  _solveWithGuess(b,x);
239  }
240 
241 protected:
242 
243 };
244 
245 
246 namespace internal {
247 
248 template<typename _MatrixType, int _UpLo, typename _Preconditioner, typename Rhs>
249 struct solve_retval<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner>, Rhs>
250  : solve_retval_base<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner>, Rhs>
251 {
253  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
254 
255  template<typename Dest> void evalTo(Dest& dst) const
256  {
257  dec()._solve(rhs(),dst);
258  }
259 };
260 
261 } // end namespace internal
262 
263 } // end namespace Eigen
264 
265 #endif // EIGEN_CONJUGATE_GRADIENT_H
A preconditioner based on the digonal entries.
IntermediateState sqrt(const Expression &arg)
EIGEN_DONT_INLINE void conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, int &iters, typename Dest::RealScalar &tol_error)
ConjugateGradient(const MatrixType &A)
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ColXpr
Definition: BlockMethods.h:15
A conjugate gradient solver for sparse self-adjoint problems.
const internal::solve_retval_with_guess< ConjugateGradient, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
RealReturnType real() const
void _solveWithGuess(const Rhs &b, Dest &x) const
void _solve(const Rhs &b, Dest &x) const
void rhs(const real_t *x, real_t *f)
MatrixType::RealScalar RealScalar
MatrixType::Scalar Scalar
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
Definition: Solve.h:61
IterativeSolverBase< ConjugateGradient > Base
#define EIGEN_DONT_INLINE
EIGEN_STRONG_INLINE Index cols() const
#define eigen_assert(x)
Base class for linear iterative solvers.
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
_Preconditioner Preconditioner


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:31