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 namespace Eigen {
14 
15 namespace internal {
16 
26 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
28 void least_square_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 m = mat.rows(), n = mat.cols();
42 
43  VectorType residual = rhs - mat * x;
44  VectorType normal_residual = mat.adjoint() * residual;
45 
46  RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
47  if(rhsNorm2 == 0)
48  {
49  x.setZero();
50  iters = 0;
51  tol_error = 0;
52  return;
53  }
54  RealScalar threshold = tol*tol*rhsNorm2;
55  RealScalar residualNorm2 = normal_residual.squaredNorm();
56  if (residualNorm2 < threshold)
57  {
58  iters = 0;
59  tol_error = sqrt(residualNorm2 / rhsNorm2);
60  return;
61  }
62 
63  VectorType p(n);
64  p = precond.solve(normal_residual); // initial search direction
65 
66  VectorType z(n), tmp(m);
67  RealScalar absNew = numext::real(normal_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;
72 
73  Scalar alpha = absNew / tmp.squaredNorm(); // the amount we travel on dir
74  x += alpha * p; // update solution
75  residual -= alpha * tmp; // update residual
76  normal_residual = mat.adjoint() * residual; // update residual of the normal equation
77 
78  residualNorm2 = normal_residual.squaredNorm();
79  if(residualNorm2 < threshold)
80  break;
81 
82  z = precond.solve(normal_residual); // approximately solve for "A'A z = normal_residual"
83 
84  RealScalar absOld = absNew;
85  absNew = numext::real(normal_residual.dot(z)); // update the absolute value of r
86  RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
87  p = z + beta * p; // update search direction
88  i++;
89  }
90  tol_error = sqrt(residualNorm2 / rhsNorm2);
91  iters = i;
92 }
93 
94 }
95 
96 template< typename _MatrixType,
97  typename _Preconditioner = LeastSquareDiagonalPreconditioner<typename _MatrixType::Scalar> >
99 
100 namespace internal {
101 
102 template< typename _MatrixType, typename _Preconditioner>
103 struct traits<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
104 {
105  typedef _MatrixType MatrixType;
106  typedef _Preconditioner Preconditioner;
107 };
108 
109 }
110 
148 template< typename _MatrixType, typename _Preconditioner>
149 class LeastSquaresConjugateGradient : public IterativeSolverBase<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
150 {
152  using Base::matrix;
153  using Base::m_error;
154  using Base::m_iterations;
155  using Base::m_info;
156  using Base::m_isInitialized;
157 public:
158  typedef _MatrixType MatrixType;
159  typedef typename MatrixType::Scalar Scalar;
161  typedef _Preconditioner Preconditioner;
162 
163 public:
164 
167 
178  template<typename MatrixDerived>
180 
182 
184  template<typename Rhs,typename Dest>
185  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
186  {
189 
192  }
193 
194 };
195 
196 } // end namespace Eigen
197 
198 #endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
Eigen
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
Eigen::LeastSquaresConjugateGradient::Scalar
MatrixType::Scalar Scalar
Definition: LeastSquareConjugateGradient.h:159
Eigen::internal::traits< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >::MatrixType
_MatrixType MatrixType
Definition: LeastSquareConjugateGradient.h:105
alpha
RealScalar alpha
Definition: level1_cplx_impl.h:147
MatrixType
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
b
Scalar * b
Definition: benchVecAdd.cpp:17
Eigen::EigenBase
Definition: EigenBase.h:29
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
Eigen::LeastSquaresConjugateGradient::RealScalar
MatrixType::RealScalar RealScalar
Definition: LeastSquareConjugateGradient.h:160
Eigen::Success
@ Success
Definition: Constants.h:442
real
float real
Definition: datatypes.h:10
Eigen::LeastSquaresConjugateGradient::Base
IterativeSolverBase< LeastSquaresConjugateGradient > Base
Definition: LeastSquareConjugateGradient.h:151
mat
MatrixXf mat
Definition: Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
beta
double beta(double a, double b)
Definition: beta.c:61
Eigen::internal::traits< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >::Preconditioner
_Preconditioner Preconditioner
Definition: LeastSquareConjugateGradient.h:106
n
int n
Definition: BiCGSTAB_simple.cpp:1
Eigen::NoConvergence
@ NoConvergence
Definition: Constants.h:446
A
Definition: test_numpy_dtypes.cpp:298
Eigen::LeastSquaresConjugateGradient::LeastSquaresConjugateGradient
LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:166
Eigen::LeastSquaresConjugateGradient::_solve_vector_with_guess_impl
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: LeastSquareConjugateGradient.h:185
Eigen::IterativeSolverBase::maxIterations
Index maxIterations() const
Definition: IterativeSolverBase.h:281
Eigen::IterativeSolverBase::m_info
ComputationInfo m_info
Definition: IterativeSolverBase.h:438
Eigen::LeastSquaresConjugateGradient::m_iterations
Index m_iterations
Definition: IterativeSolverBase.h:437
pybind_wrapper_test_script.z
z
Definition: pybind_wrapper_test_script.py:61
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
Eigen::IterativeSolverBase
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:143
RealScalar
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:47
Eigen::IterativeSolverBase::m_tolerance
RealScalar m_tolerance
Definition: IterativeSolverBase.h:434
Eigen::LeastSquaresConjugateGradient
A conjugate gradient solver for sparse (or dense) least-square problems.
Definition: LeastSquareConjugateGradient.h:98
Eigen::internal::traits
Definition: ForwardDeclarations.h:17
Eigen::LeastSquaresConjugateGradient::Preconditioner
_Preconditioner Preconditioner
Definition: LeastSquareConjugateGradient.h:161
Eigen::internal::Rhs
@ Rhs
Definition: TensorContractionMapper.h:18
Eigen::LeastSquaresConjugateGradient::MatrixType
_MatrixType MatrixType
Definition: LeastSquareConjugateGradient.h:158
Eigen::IterativeSolverBase::m_iterations
Index m_iterations
Definition: IterativeSolverBase.h:437
p
float * p
Definition: Tutorial_Map_using.cpp:9
Eigen::LeastSquaresConjugateGradient::LeastSquaresConjugateGradient
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition: LeastSquareConjugateGradient.h:179
Eigen::IterativeSolverBase::matrix
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:419
Eigen::LeastSquaresConjugateGradient::m_info
ComputationInfo m_info
Definition: IterativeSolverBase.h:438
gtsam::tol
const G double tol
Definition: Group.h:79
Eigen::LeastSquaresConjugateGradient::m_error
RealScalar m_error
Definition: IterativeSolverBase.h:436
Eigen::Matrix< Scalar, Dynamic, 1 >
Eigen::LeastSquaresConjugateGradient::~LeastSquaresConjugateGradient
~LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:181
abs
#define abs(x)
Definition: datatypes.h:17
internal
Definition: BandTriangularSolver.h:13
VectorType
Definition: FFTW.cpp:65
Eigen::IterativeSolverBase::m_preconditioner
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:431
ceres::sqrt
Jet< T, N > sqrt(const Jet< T, N > &f)
Definition: jet.h:418
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Eigen::internal::least_square_conjugate_gradient
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)
Definition: LeastSquareConjugateGradient.h:28
Eigen::IterativeSolverBase::m_error
RealScalar m_error
Definition: IterativeSolverBase.h:436
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
EIGEN_DONT_INLINE
#define EIGEN_DONT_INLINE
Definition: Macros.h:940
Eigen::SparseSolverBase::m_isInitialized
bool m_isInitialized
Definition: SparseSolverBase.h:119


gtsam
Author(s):
autogenerated on Thu Jun 13 2024 03:02:59