LeastSquareConjugateGradient.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) 2015 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_LEAST_SQUARE_CONJUGATE_GRADIENT_H
11 #define EIGEN_LEAST_SQUARE_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 least_square_conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
31  const Preconditioner& precond, Index& iters,
32  typename Dest::RealScalar& tol_error)
33 {
34  using std::sqrt;
35  using std::abs;
36  typedef typename Dest::RealScalar RealScalar;
37  typedef typename Dest::Scalar Scalar;
38  typedef Matrix<Scalar,Dynamic,1> VectorType;
39 
40  RealScalar tol = tol_error;
41  Index maxIters = iters;
42 
43  Index m = mat.rows(), n = mat.cols();
44 
45  VectorType residual = rhs - mat * x;
46  VectorType normal_residual = mat.adjoint() * residual;
47 
48  RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
49  if(rhsNorm2 == 0)
50  {
51  x.setZero();
52  iters = 0;
53  tol_error = 0;
54  return;
55  }
56  RealScalar threshold = tol*tol*rhsNorm2;
57  RealScalar residualNorm2 = normal_residual.squaredNorm();
58  if (residualNorm2 < threshold)
59  {
60  iters = 0;
61  tol_error = sqrt(residualNorm2 / rhsNorm2);
62  return;
63  }
64 
65  VectorType p(n);
66  p = precond.solve(normal_residual); // initial search direction
67 
68  VectorType z(n), tmp(m);
69  RealScalar absNew = numext::real(normal_residual.dot(p)); // the square of the absolute value of r scaled by invM
70  Index i = 0;
71  while(i < maxIters)
72  {
73  tmp.noalias() = mat * p;
74 
75  Scalar alpha = absNew / tmp.squaredNorm(); // the amount we travel on dir
76  x += alpha * p; // update solution
77  residual -= alpha * tmp; // update residual
78  normal_residual.noalias() = mat.adjoint() * residual; // update residual of the normal equation
79 
80  residualNorm2 = normal_residual.squaredNorm();
81  if(residualNorm2 < threshold)
82  break;
83 
84  z = precond.solve(normal_residual); // approximately solve for "A'A z = normal_residual"
85 
86  RealScalar absOld = absNew;
87  absNew = numext::real(normal_residual.dot(z)); // update the absolute value of r
88  RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
89  p = z + beta * p; // update search direction
90  i++;
91  }
92  tol_error = sqrt(residualNorm2 / rhsNorm2);
93  iters = i;
94 }
95 
96 }
97 
98 template< typename MatrixType_,
99  typename Preconditioner_ = LeastSquareDiagonalPreconditioner<typename MatrixType_::Scalar> >
100 class LeastSquaresConjugateGradient;
101 
102 namespace internal {
103 
104 template< typename MatrixType_, typename Preconditioner_>
105 struct traits<LeastSquaresConjugateGradient<MatrixType_,Preconditioner_> >
106 {
107  typedef MatrixType_ MatrixType;
108  typedef Preconditioner_ Preconditioner;
109 };
110 
111 }
112 
150 template< typename MatrixType_, typename Preconditioner_>
151 class LeastSquaresConjugateGradient : public IterativeSolverBase<LeastSquaresConjugateGradient<MatrixType_,Preconditioner_> >
152 {
154  using Base::matrix;
155  using Base::m_error;
156  using Base::m_iterations;
157  using Base::m_info;
158  using Base::m_isInitialized;
159 public:
160  typedef MatrixType_ MatrixType;
161  typedef typename MatrixType::Scalar Scalar;
163  typedef Preconditioner_ Preconditioner;
164 
165 public:
166 
169 
180  template<typename MatrixDerived>
182 
184 
186  template<typename Rhs,typename Dest>
187  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
188  {
191 
194  }
195 
196 };
197 
198 } // end namespace Eigen
199 
200 #endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
Matrix3f m
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
Array< int, 3, 1 > b
int n
RealReturnType real() const
MatrixXcf A
#define EIGEN_DONT_INLINE
Definition: Macros.h:844
float * p
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
Base class for linear iterative solvers.
LeastSquaresConjugateGradient< MatrixType_, Preconditioner_ > & derived()
const ActualMatrixType & matrix() const
A conjugate gradient solver for sparse (or dense) least-square problems.
IterativeSolverBase< LeastSquaresConjugateGradient > Base
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
const ActualMatrixType & matrix() const
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
EIGEN_DONT_INLINE void least_square_conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)