Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-19 09:06:20

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2015 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
0011 #define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
0012 
0013 namespace RivetEigen { 
0014 
0015 namespace internal {
0016 
0017 /** \internal Low-level conjugate gradient algorithm for least-square problems
0018   * \param mat The matrix A
0019   * \param rhs The right hand side vector b
0020   * \param x On input and initial solution, on output the computed solution.
0021   * \param precond A preconditioner being able to efficiently solve for an
0022   *                approximation of A'Ax=b (regardless of b)
0023   * \param iters On input the max number of iteration, on output the number of performed iterations.
0024   * \param tol_error On input the tolerance error, on output an estimation of the relative error.
0025   */
0026 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
0027 EIGEN_DONT_INLINE
0028 void least_square_conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
0029                                      const Preconditioner& precond, Index& iters,
0030                                      typename Dest::RealScalar& tol_error)
0031 {
0032   using std::sqrt;
0033   using std::abs;
0034   typedef typename Dest::RealScalar RealScalar;
0035   typedef typename Dest::Scalar Scalar;
0036   typedef Matrix<Scalar,Dynamic,1> VectorType;
0037   
0038   RealScalar tol = tol_error;
0039   Index maxIters = iters;
0040   
0041   Index m = mat.rows(), n = mat.cols();
0042 
0043   VectorType residual        = rhs - mat * x;
0044   VectorType normal_residual = mat.adjoint() * residual;
0045 
0046   RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
0047   if(rhsNorm2 == 0) 
0048   {
0049     x.setZero();
0050     iters = 0;
0051     tol_error = 0;
0052     return;
0053   }
0054   RealScalar threshold = tol*tol*rhsNorm2;
0055   RealScalar residualNorm2 = normal_residual.squaredNorm();
0056   if (residualNorm2 < threshold)
0057   {
0058     iters = 0;
0059     tol_error = sqrt(residualNorm2 / rhsNorm2);
0060     return;
0061   }
0062   
0063   VectorType p(n);
0064   p = precond.solve(normal_residual);                         // initial search direction
0065 
0066   VectorType z(n), tmp(m);
0067   RealScalar absNew = numext::real(normal_residual.dot(p));  // the square of the absolute value of r scaled by invM
0068   Index i = 0;
0069   while(i < maxIters)
0070   {
0071     tmp.noalias() = mat * p;
0072 
0073     Scalar alpha = absNew / tmp.squaredNorm();      // the amount we travel on dir
0074     x += alpha * p;                                 // update solution
0075     residual -= alpha * tmp;                        // update residual
0076     normal_residual = mat.adjoint() * residual;     // update residual of the normal equation
0077     
0078     residualNorm2 = normal_residual.squaredNorm();
0079     if(residualNorm2 < threshold)
0080       break;
0081     
0082     z = precond.solve(normal_residual);             // approximately solve for "A'A z = normal_residual"
0083 
0084     RealScalar absOld = absNew;
0085     absNew = numext::real(normal_residual.dot(z));  // update the absolute value of r
0086     RealScalar beta = absNew / absOld;              // calculate the Gram-Schmidt value used to create the new search direction
0087     p = z + beta * p;                               // update search direction
0088     i++;
0089   }
0090   tol_error = sqrt(residualNorm2 / rhsNorm2);
0091   iters = i;
0092 }
0093 
0094 }
0095 
0096 template< typename _MatrixType,
0097           typename _Preconditioner = LeastSquareDiagonalPreconditioner<typename _MatrixType::Scalar> >
0098 class LeastSquaresConjugateGradient;
0099 
0100 namespace internal {
0101 
0102 template< typename _MatrixType, typename _Preconditioner>
0103 struct traits<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
0104 {
0105   typedef _MatrixType MatrixType;
0106   typedef _Preconditioner Preconditioner;
0107 };
0108 
0109 }
0110 
0111 /** \ingroup IterativeLinearSolvers_Module
0112   * \brief A conjugate gradient solver for sparse (or dense) least-square problems
0113   *
0114   * This class allows to solve for A x = b linear problems using an iterative conjugate gradient algorithm.
0115   * The matrix A can be non symmetric and rectangular, but the matrix A' A should be positive-definite to guaranty stability.
0116   * Otherwise, the SparseLU or SparseQR classes might be preferable.
0117   * The matrix A and the vectors x and b can be either dense or sparse.
0118   *
0119   * \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
0120   * \tparam _Preconditioner the type of the preconditioner. Default is LeastSquareDiagonalPreconditioner
0121   *
0122   * \implsparsesolverconcept
0123   * 
0124   * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
0125   * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
0126   * and NumTraits<Scalar>::epsilon() for the tolerance.
0127   * 
0128   * This class can be used as the direct solver classes. Here is a typical usage example:
0129     \code
0130     int m=1000000, n = 10000;
0131     VectorXd x(n), b(m);
0132     SparseMatrix<double> A(m,n);
0133     // fill A and b
0134     LeastSquaresConjugateGradient<SparseMatrix<double> > lscg;
0135     lscg.compute(A);
0136     x = lscg.solve(b);
0137     std::cout << "#iterations:     " << lscg.iterations() << std::endl;
0138     std::cout << "estimated error: " << lscg.error()      << std::endl;
0139     // update b, and solve again
0140     x = lscg.solve(b);
0141     \endcode
0142   * 
0143   * By default the iterations start with x=0 as an initial guess of the solution.
0144   * One can control the start using the solveWithGuess() method.
0145   * 
0146   * \sa class ConjugateGradient, SparseLU, SparseQR
0147   */
0148 template< typename _MatrixType, typename _Preconditioner>
0149 class LeastSquaresConjugateGradient : public IterativeSolverBase<LeastSquaresConjugateGradient<_MatrixType,_Preconditioner> >
0150 {
0151   typedef IterativeSolverBase<LeastSquaresConjugateGradient> Base;
0152   using Base::matrix;
0153   using Base::m_error;
0154   using Base::m_iterations;
0155   using Base::m_info;
0156   using Base::m_isInitialized;
0157 public:
0158   typedef _MatrixType MatrixType;
0159   typedef typename MatrixType::Scalar Scalar;
0160   typedef typename MatrixType::RealScalar RealScalar;
0161   typedef _Preconditioner Preconditioner;
0162 
0163 public:
0164 
0165   /** Default constructor. */
0166   LeastSquaresConjugateGradient() : Base() {}
0167 
0168   /** Initialize the solver with matrix \a A for further \c Ax=b solving.
0169     * 
0170     * This constructor is a shortcut for the default constructor followed
0171     * by a call to compute().
0172     * 
0173     * \warning this class stores a reference to the matrix A as well as some
0174     * precomputed values that depend on it. Therefore, if \a A is changed
0175     * this class becomes invalid. Call compute() to update it with the new
0176     * matrix A, or modify a copy of A.
0177     */
0178   template<typename MatrixDerived>
0179   explicit LeastSquaresConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
0180 
0181   ~LeastSquaresConjugateGradient() {}
0182 
0183   /** \internal */
0184   template<typename Rhs,typename Dest>
0185   void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const
0186   {
0187     m_iterations = Base::maxIterations();
0188     m_error = Base::m_tolerance;
0189 
0190     internal::least_square_conjugate_gradient(matrix(), b, x, Base::m_preconditioner, m_iterations, m_error);
0191     m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
0192   }
0193 
0194 };
0195 
0196 } // end namespace RivetEigen
0197 
0198 #endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H