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 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
29 template <typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
31  const Preconditioner& precond, Index& iters,
32  typename Dest::RealScalar& tol_error) {
33  using std::abs;
34  using std::sqrt;
35  typedef typename Dest::RealScalar RealScalar;
36  typedef typename Dest::Scalar Scalar;
38 
39  RealScalar tol = tol_error;
40  Index maxIters = iters;
41 
42  Index m = mat.rows(), n = mat.cols();
43 
44  VectorType residual = rhs - mat * x;
45  VectorType normal_residual = mat.adjoint() * residual;
46 
47  RealScalar rhsNorm2 = (mat.adjoint() * rhs).squaredNorm();
48  if (rhsNorm2 == 0) {
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  iters = 0;
58  tol_error = sqrt(residualNorm2 / rhsNorm2);
59  return;
60  }
61 
62  VectorType p(n);
63  p = precond.solve(normal_residual); // initial search direction
64 
65  VectorType z(n), tmp(m);
66  RealScalar absNew = numext::real(normal_residual.dot(p)); // the square of the absolute value of r scaled by invM
67  Index i = 0;
68  while (i < maxIters) {
69  tmp.noalias() = mat * p;
70 
71  Scalar alpha = absNew / tmp.squaredNorm(); // the amount we travel on dir
72  x += alpha * p; // update solution
73  residual -= alpha * tmp; // update residual
74  normal_residual.noalias() = mat.adjoint() * residual; // update residual of the normal equation
75 
76  residualNorm2 = normal_residual.squaredNorm();
77  if (residualNorm2 < threshold) break;
78 
79  z = precond.solve(normal_residual); // approximately solve for "A'A z = normal_residual"
80 
81  RealScalar absOld = absNew;
82  absNew = numext::real(normal_residual.dot(z)); // update the absolute value of r
83  RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
84  p = z + beta * p; // update search direction
85  i++;
86  }
87  tol_error = sqrt(residualNorm2 / rhsNorm2);
88  iters = i;
89 }
90 
91 } // namespace internal
92 
93 template <typename MatrixType_,
94  typename Preconditioner_ = LeastSquareDiagonalPreconditioner<typename MatrixType_::Scalar> >
95 class LeastSquaresConjugateGradient;
96 
97 namespace internal {
98 
99 template <typename MatrixType_, typename Preconditioner_>
100 struct traits<LeastSquaresConjugateGradient<MatrixType_, Preconditioner_> > {
101  typedef MatrixType_ MatrixType;
102  typedef Preconditioner_ Preconditioner;
103 };
104 
105 } // namespace internal
106 
145 template <typename MatrixType_, typename Preconditioner_>
147  : public IterativeSolverBase<LeastSquaresConjugateGradient<MatrixType_, Preconditioner_> > {
149  using Base::m_error;
150  using Base::m_info;
151  using Base::m_isInitialized;
152  using Base::m_iterations;
153  using Base::matrix;
154 
155  public:
156  typedef MatrixType_ MatrixType;
157  typedef typename MatrixType::Scalar Scalar;
159  typedef Preconditioner_ Preconditioner;
160 
161  public:
164 
175  template <typename MatrixDerived>
177 
179 
181  template <typename Rhs, typename Dest>
182  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const {
185 
188  }
189 };
190 
191 } // end namespace Eigen
192 
193 #endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define EIGEN_DONT_INLINE
Definition: Macros.h:853
float * p
Definition: Tutorial_Map_using.cpp:9
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:124
Index maxIterations() const
Definition: IterativeSolverBase.h:251
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
RealScalar m_error
Definition: IterativeSolverBase.h:387
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:382
Index m_iterations
Definition: IterativeSolverBase.h:388
bool m_isInitialized
Definition: SparseSolverBase.h:110
RealScalar m_tolerance
Definition: IterativeSolverBase.h:385
LeastSquaresConjugateGradient< MatrixType_, Preconditioner_ > & derived()
Definition: SparseSolverBase.h:76
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
A conjugate gradient solver for sparse (or dense) least-square problems.
Definition: LeastSquareConjugateGradient.h:147
MatrixType_ MatrixType
Definition: LeastSquareConjugateGradient.h:156
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
IterativeSolverBase< LeastSquaresConjugateGradient > Base
Definition: LeastSquareConjugateGradient.h:148
~LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:178
MatrixType::RealScalar RealScalar
Definition: LeastSquareConjugateGradient.h:158
void _solve_vector_with_guess_impl(const Rhs &b, Dest &x) const
Definition: LeastSquareConjugateGradient.h:182
Preconditioner_ Preconditioner
Definition: LeastSquareConjugateGradient.h:159
MatrixType::Scalar Scalar
Definition: LeastSquareConjugateGradient.h:157
RealScalar m_error
Definition: IterativeSolverBase.h:387
Index m_iterations
Definition: IterativeSolverBase.h:388
LeastSquaresConjugateGradient()
Definition: LeastSquareConjugateGradient.h:163
LeastSquaresConjugateGradient(const EigenBase< MatrixDerived > &A)
Definition: LeastSquareConjugateGradient.h:176
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
const AdjointReturnType adjoint() const
Definition: SparseMatrixBase.h:360
Index cols() const
Definition: SparseMatrix.h:161
Index rows() const
Definition: SparseMatrix.h:159
float real
Definition: datatypes.h:10
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
RealScalar alpha
Definition: level1_cplx_impl.h:151
int * m
Definition: level2_cplx_impl.h:294
Scalar beta
Definition: level2_cplx_impl.h:36
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
@ Rhs
Definition: TensorContractionMapper.h:20
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:30
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
unsigned Preconditioner
----------------------—Domain Properties------------------------—
Definition: space_time_oscillating_cylinder.cc:725
Definition: Eigen_Colamd.h:49
list x
Definition: plotDoE.py:28
Definition: EigenBase.h:33
MatrixType_ MatrixType
Definition: LeastSquareConjugateGradient.h:101
Preconditioner_ Preconditioner
Definition: LeastSquareConjugateGradient.h:102
Definition: ForwardDeclarations.h:21
Definition: fft_test_shared.h:66