xref: /aosp_15_r20/external/eigen/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2015 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li #ifndef EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
11*bf2c3715SXin Li #define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
12*bf2c3715SXin Li 
13*bf2c3715SXin Li namespace Eigen {
14*bf2c3715SXin Li 
15*bf2c3715SXin Li namespace internal {
16*bf2c3715SXin Li 
17*bf2c3715SXin Li /** \internal Low-level conjugate gradient algorithm for least-square problems
18*bf2c3715SXin Li   * \param mat The matrix A
19*bf2c3715SXin Li   * \param rhs The right hand side vector b
20*bf2c3715SXin Li   * \param x On input and initial solution, on output the computed solution.
21*bf2c3715SXin Li   * \param precond A preconditioner being able to efficiently solve for an
22*bf2c3715SXin Li   *                approximation of A'Ax=b (regardless of b)
23*bf2c3715SXin Li   * \param iters On input the max number of iteration, on output the number of performed iterations.
24*bf2c3715SXin Li   * \param tol_error On input the tolerance error, on output an estimation of the relative error.
25*bf2c3715SXin Li   */
26*bf2c3715SXin Li template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
27*bf2c3715SXin Li EIGEN_DONT_INLINE
least_square_conjugate_gradient(const MatrixType & mat,const Rhs & rhs,Dest & x,const Preconditioner & precond,Index & iters,typename Dest::RealScalar & tol_error)28*bf2c3715SXin Li void least_square_conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
29*bf2c3715SXin Li                                      const Preconditioner& precond, Index& iters,
30*bf2c3715SXin Li                                      typename Dest::RealScalar& tol_error)
31*bf2c3715SXin Li {
32*bf2c3715SXin Li   using std::sqrt;
33*bf2c3715SXin Li   using std::abs;
34*bf2c3715SXin Li   typedef typename Dest::RealScalar RealScalar;
35*bf2c3715SXin Li   typedef typename Dest::Scalar Scalar;
36*bf2c3715SXin Li   typedef Matrix<Scalar,Dynamic,1> VectorType;
37*bf2c3715SXin Li 
38*bf2c3715SXin Li   RealScalar tol = tol_error;
39*bf2c3715SXin Li   Index maxIters = iters;
40*bf2c3715SXin Li 
41*bf2c3715SXin Li   Index m = mat.rows(), n = mat.cols();
42*bf2c3715SXin Li 
43*bf2c3715SXin Li   VectorType residual        = rhs - mat * x;
44*bf2c3715SXin Li   VectorType normal_residual = mat.adjoint() * residual;
45*bf2c3715SXin Li 
46*bf2c3715SXin Li   RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
47*bf2c3715SXin Li   if(rhsNorm2 == 0)
48*bf2c3715SXin Li   {
49*bf2c3715SXin Li     x.setZero();
50*bf2c3715SXin Li     iters = 0;
51*bf2c3715SXin Li     tol_error = 0;
52*bf2c3715SXin Li     return;
53*bf2c3715SXin Li   }
54*bf2c3715SXin Li   RealScalar threshold = tol*tol*rhsNorm2;
55*bf2c3715SXin Li   RealScalar residualNorm2 = normal_residual.squaredNorm();
56*bf2c3715SXin Li   if (residualNorm2 < threshold)
57*bf2c3715SXin Li   {
58*bf2c3715SXin Li     iters = 0;
59*bf2c3715SXin Li     tol_error = sqrt(residualNorm2 / rhsNorm2);
60*bf2c3715SXin Li     return;
61*bf2c3715SXin Li   }
62*bf2c3715SXin Li 
63*bf2c3715SXin Li   VectorType p(n);
64*bf2c3715SXin Li   p = precond.solve(normal_residual);                         // initial search direction
65*bf2c3715SXin Li 
66*bf2c3715SXin Li   VectorType z(n), tmp(m);
67*bf2c3715SXin Li   RealScalar absNew = numext::real(normal_residual.dot(p));  // the square of the absolute value of r scaled by invM
68*bf2c3715SXin Li   Index i = 0;
69*bf2c3715SXin Li   while(i < maxIters)
70*bf2c3715SXin Li   {
71*bf2c3715SXin Li     tmp.noalias() = mat * p;
72*bf2c3715SXin Li 
73*bf2c3715SXin Li     Scalar alpha = absNew / tmp.squaredNorm();      // the amount we travel on dir
74*bf2c3715SXin Li     x += alpha * p;                                 // update solution
75*bf2c3715SXin Li     residual -= alpha * tmp;                        // update residual
76*bf2c3715SXin Li     normal_residual = mat.adjoint() * residual;     // update residual of the normal equation
77*bf2c3715SXin Li 
78*bf2c3715SXin Li     residualNorm2 = normal_residual.squaredNorm();
79*bf2c3715SXin Li     if(residualNorm2 < threshold)
80*bf2c3715SXin Li       break;
81*bf2c3715SXin Li 
82*bf2c3715SXin Li     z = precond.solve(normal_residual);             // approximately solve for "A'A z = normal_residual"
83*bf2c3715SXin Li 
84*bf2c3715SXin Li     RealScalar absOld = absNew;
85*bf2c3715SXin Li     absNew = numext::real(normal_residual.dot(z));  // update the absolute value of r
86*bf2c3715SXin Li     RealScalar beta = absNew / absOld;              // calculate the Gram-Schmidt value used to create the new search direction
87*bf2c3715SXin Li     p = z + beta * p;                               // update search direction
88*bf2c3715SXin Li     i++;
89*bf2c3715SXin Li   }
90*bf2c3715SXin Li   tol_error = sqrt(residualNorm2 / rhsNorm2);
91*bf2c3715SXin Li   iters = i;
92*bf2c3715SXin Li }
93*bf2c3715SXin Li 
94*bf2c3715SXin Li }
95*bf2c3715SXin Li 
96*bf2c3715SXin Li template< typename _MatrixType,
97*bf2c3715SXin Li           typename _Preconditioner = LeastSquareDiagonalPreconditioner<typename _MatrixType::Scalar> >
98*bf2c3715SXin Li class LeastSquaresConjugateGradient;
99*bf2c3715SXin Li 
100*bf2c3715SXin Li namespace internal {
101*bf2c3715SXin Li 
102*bf2c3715SXin Li template< typename _MatrixType, typename _Preconditioner>
103*bf2c3715SXin Li struct traits<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
104*bf2c3715SXin Li {
105*bf2c3715SXin Li   typedef _MatrixType MatrixType;
106*bf2c3715SXin Li   typedef _Preconditioner Preconditioner;
107*bf2c3715SXin Li };
108*bf2c3715SXin Li 
109*bf2c3715SXin Li }
110*bf2c3715SXin Li 
111*bf2c3715SXin Li /** \ingroup IterativeLinearSolvers_Module
112*bf2c3715SXin Li   * \brief A conjugate gradient solver for sparse (or dense) least-square problems
113*bf2c3715SXin Li   *
114*bf2c3715SXin Li   * This class allows to solve for A x = b linear problems using an iterative conjugate gradient algorithm.
115*bf2c3715SXin Li   * The matrix A can be non symmetric and rectangular, but the matrix A' A should be positive-definite to guaranty stability.
116*bf2c3715SXin Li   * Otherwise, the SparseLU or SparseQR classes might be preferable.
117*bf2c3715SXin Li   * The matrix A and the vectors x and b can be either dense or sparse.
118*bf2c3715SXin Li   *
119*bf2c3715SXin Li   * \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
120*bf2c3715SXin Li   * \tparam _Preconditioner the type of the preconditioner. Default is LeastSquareDiagonalPreconditioner
121*bf2c3715SXin Li   *
122*bf2c3715SXin Li   * \implsparsesolverconcept
123*bf2c3715SXin Li   *
124*bf2c3715SXin Li   * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
125*bf2c3715SXin Li   * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
126*bf2c3715SXin Li   * and NumTraits<Scalar>::epsilon() for the tolerance.
127*bf2c3715SXin Li   *
128*bf2c3715SXin Li   * This class can be used as the direct solver classes. Here is a typical usage example:
129*bf2c3715SXin Li     \code
130*bf2c3715SXin Li     int m=1000000, n = 10000;
131*bf2c3715SXin Li     VectorXd x(n), b(m);
132*bf2c3715SXin Li     SparseMatrix<double> A(m,n);
133*bf2c3715SXin Li     // fill A and b
134*bf2c3715SXin Li     LeastSquaresConjugateGradient<SparseMatrix<double> > lscg;
135*bf2c3715SXin Li     lscg.compute(A);
136*bf2c3715SXin Li     x = lscg.solve(b);
137*bf2c3715SXin Li     std::cout << "#iterations:     " << lscg.iterations() << std::endl;
138*bf2c3715SXin Li     std::cout << "estimated error: " << lscg.error()      << std::endl;
139*bf2c3715SXin Li     // update b, and solve again
140*bf2c3715SXin Li     x = lscg.solve(b);
141*bf2c3715SXin Li     \endcode
142*bf2c3715SXin Li   *
143*bf2c3715SXin Li   * By default the iterations start with x=0 as an initial guess of the solution.
144*bf2c3715SXin Li   * One can control the start using the solveWithGuess() method.
145*bf2c3715SXin Li   *
146*bf2c3715SXin Li   * \sa class ConjugateGradient, SparseLU, SparseQR
147*bf2c3715SXin Li   */
148*bf2c3715SXin Li template< typename _MatrixType, typename _Preconditioner>
149*bf2c3715SXin Li class LeastSquaresConjugateGradient : public IterativeSolverBase<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
150*bf2c3715SXin Li {
151*bf2c3715SXin Li   typedef IterativeSolverBase<LeastSquaresConjugateGradient> Base;
152*bf2c3715SXin Li   using Base::matrix;
153*bf2c3715SXin Li   using Base::m_error;
154*bf2c3715SXin Li   using Base::m_iterations;
155*bf2c3715SXin Li   using Base::m_info;
156*bf2c3715SXin Li   using Base::m_isInitialized;
157*bf2c3715SXin Li public:
158*bf2c3715SXin Li   typedef _MatrixType MatrixType;
159*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
160*bf2c3715SXin Li   typedef typename MatrixType::RealScalar RealScalar;
161*bf2c3715SXin Li   typedef _Preconditioner Preconditioner;
162*bf2c3715SXin Li 
163*bf2c3715SXin Li public:
164*bf2c3715SXin Li 
165*bf2c3715SXin Li   /** Default constructor. */
166*bf2c3715SXin Li   LeastSquaresConjugateGradient() : Base() {}
167*bf2c3715SXin Li 
168*bf2c3715SXin Li   /** Initialize the solver with matrix \a A for further \c Ax=b solving.
169*bf2c3715SXin Li     *
170*bf2c3715SXin Li     * This constructor is a shortcut for the default constructor followed
171*bf2c3715SXin Li     * by a call to compute().
172*bf2c3715SXin Li     *
173*bf2c3715SXin Li     * \warning this class stores a reference to the matrix A as well as some
174*bf2c3715SXin Li     * precomputed values that depend on it. Therefore, if \a A is changed
175*bf2c3715SXin Li     * this class becomes invalid. Call compute() to update it with the new
176*bf2c3715SXin Li     * matrix A, or modify a copy of A.
177*bf2c3715SXin Li     */
178*bf2c3715SXin Li   template<typename MatrixDerived>
179*bf2c3715SXin Li   explicit LeastSquaresConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
180*bf2c3715SXin Li 
181*bf2c3715SXin Li   ~LeastSquaresConjugateGradient() {}
182*bf2c3715SXin Li 
183*bf2c3715SXin Li   /** \internal */
184*bf2c3715SXin Li   template<typename Rhs,typename Dest>
185*bf2c3715SXin Li   void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
186*bf2c3715SXin Li   {
187*bf2c3715SXin Li     m_iterations = Base::maxIterations();
188*bf2c3715SXin Li     m_error = Base::m_tolerance;
189*bf2c3715SXin Li 
190*bf2c3715SXin Li     internal::least_square_conjugate_gradient(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
191*bf2c3715SXin Li     m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
192*bf2c3715SXin Li   }
193*bf2c3715SXin Li 
194*bf2c3715SXin Li };
195*bf2c3715SXin Li 
196*bf2c3715SXin Li } // end namespace Eigen
197*bf2c3715SXin Li 
198*bf2c3715SXin Li #endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
199