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) 2011-2014 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_ITERATIVE_SOLVER_BASE_H 11*bf2c3715SXin Li #define EIGEN_ITERATIVE_SOLVER_BASE_H 12*bf2c3715SXin Li 13*bf2c3715SXin Li namespace Eigen { 14*bf2c3715SXin Li 15*bf2c3715SXin Li namespace internal { 16*bf2c3715SXin Li 17*bf2c3715SXin Li template<typename MatrixType> 18*bf2c3715SXin Li struct is_ref_compatible_impl 19*bf2c3715SXin Li { 20*bf2c3715SXin Li private: 21*bf2c3715SXin Li template <typename T0> 22*bf2c3715SXin Li struct any_conversion 23*bf2c3715SXin Li { 24*bf2c3715SXin Li template <typename T> any_conversion(const volatile T&); 25*bf2c3715SXin Li template <typename T> any_conversion(T&); 26*bf2c3715SXin Li }; 27*bf2c3715SXin Li struct yes {int a[1];}; 28*bf2c3715SXin Li struct no {int a[2];}; 29*bf2c3715SXin Li 30*bf2c3715SXin Li template<typename T> 31*bf2c3715SXin Li static yes test(const Ref<const T>&, int); 32*bf2c3715SXin Li template<typename T> 33*bf2c3715SXin Li static no test(any_conversion<T>, ...); 34*bf2c3715SXin Li 35*bf2c3715SXin Li public: 36*bf2c3715SXin Li static MatrixType ms_from; 37*bf2c3715SXin Li enum { value = sizeof(test<MatrixType>(ms_from, 0))==sizeof(yes) }; 38*bf2c3715SXin Li }; 39*bf2c3715SXin Li 40*bf2c3715SXin Li template<typename MatrixType> 41*bf2c3715SXin Li struct is_ref_compatible 42*bf2c3715SXin Li { 43*bf2c3715SXin Li enum { value = is_ref_compatible_impl<typename remove_all<MatrixType>::type>::value }; 44*bf2c3715SXin Li }; 45*bf2c3715SXin Li 46*bf2c3715SXin Li template<typename MatrixType, bool MatrixFree = !internal::is_ref_compatible<MatrixType>::value> 47*bf2c3715SXin Li class generic_matrix_wrapper; 48*bf2c3715SXin Li 49*bf2c3715SXin Li // We have an explicit matrix at hand, compatible with Ref<> 50*bf2c3715SXin Li template<typename MatrixType> 51*bf2c3715SXin Li class generic_matrix_wrapper<MatrixType,false> 52*bf2c3715SXin Li { 53*bf2c3715SXin Li public: 54*bf2c3715SXin Li typedef Ref<const MatrixType> ActualMatrixType; 55*bf2c3715SXin Li template<int UpLo> struct ConstSelfAdjointViewReturnType { 56*bf2c3715SXin Li typedef typename ActualMatrixType::template ConstSelfAdjointViewReturnType<UpLo>::Type Type; 57*bf2c3715SXin Li }; 58*bf2c3715SXin Li 59*bf2c3715SXin Li enum { 60*bf2c3715SXin Li MatrixFree = false 61*bf2c3715SXin Li }; 62*bf2c3715SXin Li generic_matrix_wrapper()63*bf2c3715SXin Li generic_matrix_wrapper() 64*bf2c3715SXin Li : m_dummy(0,0), m_matrix(m_dummy) 65*bf2c3715SXin Li {} 66*bf2c3715SXin Li 67*bf2c3715SXin Li template<typename InputType> generic_matrix_wrapper(const InputType & mat)68*bf2c3715SXin Li generic_matrix_wrapper(const InputType &mat) 69*bf2c3715SXin Li : m_matrix(mat) 70*bf2c3715SXin Li {} 71*bf2c3715SXin Li matrix()72*bf2c3715SXin Li const ActualMatrixType& matrix() const 73*bf2c3715SXin Li { 74*bf2c3715SXin Li return m_matrix; 75*bf2c3715SXin Li } 76*bf2c3715SXin Li 77*bf2c3715SXin Li template<typename MatrixDerived> grab(const EigenBase<MatrixDerived> & mat)78*bf2c3715SXin Li void grab(const EigenBase<MatrixDerived> &mat) 79*bf2c3715SXin Li { 80*bf2c3715SXin Li m_matrix.~Ref<const MatrixType>(); 81*bf2c3715SXin Li ::new (&m_matrix) Ref<const MatrixType>(mat.derived()); 82*bf2c3715SXin Li } 83*bf2c3715SXin Li grab(const Ref<const MatrixType> & mat)84*bf2c3715SXin Li void grab(const Ref<const MatrixType> &mat) 85*bf2c3715SXin Li { 86*bf2c3715SXin Li if(&(mat.derived()) != &m_matrix) 87*bf2c3715SXin Li { 88*bf2c3715SXin Li m_matrix.~Ref<const MatrixType>(); 89*bf2c3715SXin Li ::new (&m_matrix) Ref<const MatrixType>(mat); 90*bf2c3715SXin Li } 91*bf2c3715SXin Li } 92*bf2c3715SXin Li 93*bf2c3715SXin Li protected: 94*bf2c3715SXin Li MatrixType m_dummy; // used to default initialize the Ref<> object 95*bf2c3715SXin Li ActualMatrixType m_matrix; 96*bf2c3715SXin Li }; 97*bf2c3715SXin Li 98*bf2c3715SXin Li // MatrixType is not compatible with Ref<> -> matrix-free wrapper 99*bf2c3715SXin Li template<typename MatrixType> 100*bf2c3715SXin Li class generic_matrix_wrapper<MatrixType,true> 101*bf2c3715SXin Li { 102*bf2c3715SXin Li public: 103*bf2c3715SXin Li typedef MatrixType ActualMatrixType; 104*bf2c3715SXin Li template<int UpLo> struct ConstSelfAdjointViewReturnType 105*bf2c3715SXin Li { 106*bf2c3715SXin Li typedef ActualMatrixType Type; 107*bf2c3715SXin Li }; 108*bf2c3715SXin Li 109*bf2c3715SXin Li enum { 110*bf2c3715SXin Li MatrixFree = true 111*bf2c3715SXin Li }; 112*bf2c3715SXin Li generic_matrix_wrapper()113*bf2c3715SXin Li generic_matrix_wrapper() 114*bf2c3715SXin Li : mp_matrix(0) 115*bf2c3715SXin Li {} 116*bf2c3715SXin Li generic_matrix_wrapper(const MatrixType & mat)117*bf2c3715SXin Li generic_matrix_wrapper(const MatrixType &mat) 118*bf2c3715SXin Li : mp_matrix(&mat) 119*bf2c3715SXin Li {} 120*bf2c3715SXin Li matrix()121*bf2c3715SXin Li const ActualMatrixType& matrix() const 122*bf2c3715SXin Li { 123*bf2c3715SXin Li return *mp_matrix; 124*bf2c3715SXin Li } 125*bf2c3715SXin Li grab(const MatrixType & mat)126*bf2c3715SXin Li void grab(const MatrixType &mat) 127*bf2c3715SXin Li { 128*bf2c3715SXin Li mp_matrix = &mat; 129*bf2c3715SXin Li } 130*bf2c3715SXin Li 131*bf2c3715SXin Li protected: 132*bf2c3715SXin Li const ActualMatrixType *mp_matrix; 133*bf2c3715SXin Li }; 134*bf2c3715SXin Li 135*bf2c3715SXin Li } 136*bf2c3715SXin Li 137*bf2c3715SXin Li /** \ingroup IterativeLinearSolvers_Module 138*bf2c3715SXin Li * \brief Base class for linear iterative solvers 139*bf2c3715SXin Li * 140*bf2c3715SXin Li * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner 141*bf2c3715SXin Li */ 142*bf2c3715SXin Li template< typename Derived> 143*bf2c3715SXin Li class IterativeSolverBase : public SparseSolverBase<Derived> 144*bf2c3715SXin Li { 145*bf2c3715SXin Li protected: 146*bf2c3715SXin Li typedef SparseSolverBase<Derived> Base; 147*bf2c3715SXin Li using Base::m_isInitialized; 148*bf2c3715SXin Li 149*bf2c3715SXin Li public: 150*bf2c3715SXin Li typedef typename internal::traits<Derived>::MatrixType MatrixType; 151*bf2c3715SXin Li typedef typename internal::traits<Derived>::Preconditioner Preconditioner; 152*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar; 153*bf2c3715SXin Li typedef typename MatrixType::StorageIndex StorageIndex; 154*bf2c3715SXin Li typedef typename MatrixType::RealScalar RealScalar; 155*bf2c3715SXin Li 156*bf2c3715SXin Li enum { 157*bf2c3715SXin Li ColsAtCompileTime = MatrixType::ColsAtCompileTime, 158*bf2c3715SXin Li MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime 159*bf2c3715SXin Li }; 160*bf2c3715SXin Li 161*bf2c3715SXin Li public: 162*bf2c3715SXin Li 163*bf2c3715SXin Li using Base::derived; 164*bf2c3715SXin Li 165*bf2c3715SXin Li /** Default constructor. */ IterativeSolverBase()166*bf2c3715SXin Li IterativeSolverBase() 167*bf2c3715SXin Li { 168*bf2c3715SXin Li init(); 169*bf2c3715SXin Li } 170*bf2c3715SXin Li 171*bf2c3715SXin Li /** Initialize the solver with matrix \a A for further \c Ax=b solving. 172*bf2c3715SXin Li * 173*bf2c3715SXin Li * This constructor is a shortcut for the default constructor followed 174*bf2c3715SXin Li * by a call to compute(). 175*bf2c3715SXin Li * 176*bf2c3715SXin Li * \warning this class stores a reference to the matrix A as well as some 177*bf2c3715SXin Li * precomputed values that depend on it. Therefore, if \a A is changed 178*bf2c3715SXin Li * this class becomes invalid. Call compute() to update it with the new 179*bf2c3715SXin Li * matrix A, or modify a copy of A. 180*bf2c3715SXin Li */ 181*bf2c3715SXin Li template<typename MatrixDerived> IterativeSolverBase(const EigenBase<MatrixDerived> & A)182*bf2c3715SXin Li explicit IterativeSolverBase(const EigenBase<MatrixDerived>& A) 183*bf2c3715SXin Li : m_matrixWrapper(A.derived()) 184*bf2c3715SXin Li { 185*bf2c3715SXin Li init(); 186*bf2c3715SXin Li compute(matrix()); 187*bf2c3715SXin Li } 188*bf2c3715SXin Li ~IterativeSolverBase()189*bf2c3715SXin Li ~IterativeSolverBase() {} 190*bf2c3715SXin Li 191*bf2c3715SXin Li /** Initializes the iterative solver for the sparsity pattern of the matrix \a A for further solving \c Ax=b problems. 192*bf2c3715SXin Li * 193*bf2c3715SXin Li * Currently, this function mostly calls analyzePattern on the preconditioner. In the future 194*bf2c3715SXin Li * we might, for instance, implement column reordering for faster matrix vector products. 195*bf2c3715SXin Li */ 196*bf2c3715SXin Li template<typename MatrixDerived> analyzePattern(const EigenBase<MatrixDerived> & A)197*bf2c3715SXin Li Derived& analyzePattern(const EigenBase<MatrixDerived>& A) 198*bf2c3715SXin Li { 199*bf2c3715SXin Li grab(A.derived()); 200*bf2c3715SXin Li m_preconditioner.analyzePattern(matrix()); 201*bf2c3715SXin Li m_isInitialized = true; 202*bf2c3715SXin Li m_analysisIsOk = true; 203*bf2c3715SXin Li m_info = m_preconditioner.info(); 204*bf2c3715SXin Li return derived(); 205*bf2c3715SXin Li } 206*bf2c3715SXin Li 207*bf2c3715SXin Li /** Initializes the iterative solver with the numerical values of the matrix \a A for further solving \c Ax=b problems. 208*bf2c3715SXin Li * 209*bf2c3715SXin Li * Currently, this function mostly calls factorize on the preconditioner. 210*bf2c3715SXin Li * 211*bf2c3715SXin Li * \warning this class stores a reference to the matrix A as well as some 212*bf2c3715SXin Li * precomputed values that depend on it. Therefore, if \a A is changed 213*bf2c3715SXin Li * this class becomes invalid. Call compute() to update it with the new 214*bf2c3715SXin Li * matrix A, or modify a copy of A. 215*bf2c3715SXin Li */ 216*bf2c3715SXin Li template<typename MatrixDerived> factorize(const EigenBase<MatrixDerived> & A)217*bf2c3715SXin Li Derived& factorize(const EigenBase<MatrixDerived>& A) 218*bf2c3715SXin Li { 219*bf2c3715SXin Li eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); 220*bf2c3715SXin Li grab(A.derived()); 221*bf2c3715SXin Li m_preconditioner.factorize(matrix()); 222*bf2c3715SXin Li m_factorizationIsOk = true; 223*bf2c3715SXin Li m_info = m_preconditioner.info(); 224*bf2c3715SXin Li return derived(); 225*bf2c3715SXin Li } 226*bf2c3715SXin Li 227*bf2c3715SXin Li /** Initializes the iterative solver with the matrix \a A for further solving \c Ax=b problems. 228*bf2c3715SXin Li * 229*bf2c3715SXin Li * Currently, this function mostly initializes/computes the preconditioner. In the future 230*bf2c3715SXin Li * we might, for instance, implement column reordering for faster matrix vector products. 231*bf2c3715SXin Li * 232*bf2c3715SXin Li * \warning this class stores a reference to the matrix A as well as some 233*bf2c3715SXin Li * precomputed values that depend on it. Therefore, if \a A is changed 234*bf2c3715SXin Li * this class becomes invalid. Call compute() to update it with the new 235*bf2c3715SXin Li * matrix A, or modify a copy of A. 236*bf2c3715SXin Li */ 237*bf2c3715SXin Li template<typename MatrixDerived> compute(const EigenBase<MatrixDerived> & A)238*bf2c3715SXin Li Derived& compute(const EigenBase<MatrixDerived>& A) 239*bf2c3715SXin Li { 240*bf2c3715SXin Li grab(A.derived()); 241*bf2c3715SXin Li m_preconditioner.compute(matrix()); 242*bf2c3715SXin Li m_isInitialized = true; 243*bf2c3715SXin Li m_analysisIsOk = true; 244*bf2c3715SXin Li m_factorizationIsOk = true; 245*bf2c3715SXin Li m_info = m_preconditioner.info(); 246*bf2c3715SXin Li return derived(); 247*bf2c3715SXin Li } 248*bf2c3715SXin Li 249*bf2c3715SXin Li /** \internal */ rows()250*bf2c3715SXin Li EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return matrix().rows(); } 251*bf2c3715SXin Li 252*bf2c3715SXin Li /** \internal */ cols()253*bf2c3715SXin Li EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return matrix().cols(); } 254*bf2c3715SXin Li 255*bf2c3715SXin Li /** \returns the tolerance threshold used by the stopping criteria. 256*bf2c3715SXin Li * \sa setTolerance() 257*bf2c3715SXin Li */ tolerance()258*bf2c3715SXin Li RealScalar tolerance() const { return m_tolerance; } 259*bf2c3715SXin Li 260*bf2c3715SXin Li /** Sets the tolerance threshold used by the stopping criteria. 261*bf2c3715SXin Li * 262*bf2c3715SXin Li * This value is used as an upper bound to the relative residual error: |Ax-b|/|b|. 263*bf2c3715SXin Li * The default value is the machine precision given by NumTraits<Scalar>::epsilon() 264*bf2c3715SXin Li */ setTolerance(const RealScalar & tolerance)265*bf2c3715SXin Li Derived& setTolerance(const RealScalar& tolerance) 266*bf2c3715SXin Li { 267*bf2c3715SXin Li m_tolerance = tolerance; 268*bf2c3715SXin Li return derived(); 269*bf2c3715SXin Li } 270*bf2c3715SXin Li 271*bf2c3715SXin Li /** \returns a read-write reference to the preconditioner for custom configuration. */ preconditioner()272*bf2c3715SXin Li Preconditioner& preconditioner() { return m_preconditioner; } 273*bf2c3715SXin Li 274*bf2c3715SXin Li /** \returns a read-only reference to the preconditioner. */ preconditioner()275*bf2c3715SXin Li const Preconditioner& preconditioner() const { return m_preconditioner; } 276*bf2c3715SXin Li 277*bf2c3715SXin Li /** \returns the max number of iterations. 278*bf2c3715SXin Li * It is either the value set by setMaxIterations or, by default, 279*bf2c3715SXin Li * twice the number of columns of the matrix. 280*bf2c3715SXin Li */ maxIterations()281*bf2c3715SXin Li Index maxIterations() const 282*bf2c3715SXin Li { 283*bf2c3715SXin Li return (m_maxIterations<0) ? 2*matrix().cols() : m_maxIterations; 284*bf2c3715SXin Li } 285*bf2c3715SXin Li 286*bf2c3715SXin Li /** Sets the max number of iterations. 287*bf2c3715SXin Li * Default is twice the number of columns of the matrix. 288*bf2c3715SXin Li */ setMaxIterations(Index maxIters)289*bf2c3715SXin Li Derived& setMaxIterations(Index maxIters) 290*bf2c3715SXin Li { 291*bf2c3715SXin Li m_maxIterations = maxIters; 292*bf2c3715SXin Li return derived(); 293*bf2c3715SXin Li } 294*bf2c3715SXin Li 295*bf2c3715SXin Li /** \returns the number of iterations performed during the last solve */ iterations()296*bf2c3715SXin Li Index iterations() const 297*bf2c3715SXin Li { 298*bf2c3715SXin Li eigen_assert(m_isInitialized && "ConjugateGradient is not initialized."); 299*bf2c3715SXin Li return m_iterations; 300*bf2c3715SXin Li } 301*bf2c3715SXin Li 302*bf2c3715SXin Li /** \returns the tolerance error reached during the last solve. 303*bf2c3715SXin Li * It is a close approximation of the true relative residual error |Ax-b|/|b|. 304*bf2c3715SXin Li */ error()305*bf2c3715SXin Li RealScalar error() const 306*bf2c3715SXin Li { 307*bf2c3715SXin Li eigen_assert(m_isInitialized && "ConjugateGradient is not initialized."); 308*bf2c3715SXin Li return m_error; 309*bf2c3715SXin Li } 310*bf2c3715SXin Li 311*bf2c3715SXin Li /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A 312*bf2c3715SXin Li * and \a x0 as an initial solution. 313*bf2c3715SXin Li * 314*bf2c3715SXin Li * \sa solve(), compute() 315*bf2c3715SXin Li */ 316*bf2c3715SXin Li template<typename Rhs,typename Guess> 317*bf2c3715SXin Li inline const SolveWithGuess<Derived, Rhs, Guess> solveWithGuess(const MatrixBase<Rhs> & b,const Guess & x0)318*bf2c3715SXin Li solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const 319*bf2c3715SXin Li { 320*bf2c3715SXin Li eigen_assert(m_isInitialized && "Solver is not initialized."); 321*bf2c3715SXin Li eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b"); 322*bf2c3715SXin Li return SolveWithGuess<Derived, Rhs, Guess>(derived(), b.derived(), x0); 323*bf2c3715SXin Li } 324*bf2c3715SXin Li 325*bf2c3715SXin Li /** \returns Success if the iterations converged, and NoConvergence otherwise. */ info()326*bf2c3715SXin Li ComputationInfo info() const 327*bf2c3715SXin Li { 328*bf2c3715SXin Li eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized."); 329*bf2c3715SXin Li return m_info; 330*bf2c3715SXin Li } 331*bf2c3715SXin Li 332*bf2c3715SXin Li /** \internal */ 333*bf2c3715SXin Li template<typename Rhs, typename DestDerived> _solve_with_guess_impl(const Rhs & b,SparseMatrixBase<DestDerived> & aDest)334*bf2c3715SXin Li void _solve_with_guess_impl(const Rhs& b, SparseMatrixBase<DestDerived> &aDest) const 335*bf2c3715SXin Li { 336*bf2c3715SXin Li eigen_assert(rows()==b.rows()); 337*bf2c3715SXin Li 338*bf2c3715SXin Li Index rhsCols = b.cols(); 339*bf2c3715SXin Li Index size = b.rows(); 340*bf2c3715SXin Li DestDerived& dest(aDest.derived()); 341*bf2c3715SXin Li typedef typename DestDerived::Scalar DestScalar; 342*bf2c3715SXin Li Eigen::Matrix<DestScalar,Dynamic,1> tb(size); 343*bf2c3715SXin Li Eigen::Matrix<DestScalar,Dynamic,1> tx(cols()); 344*bf2c3715SXin Li // We do not directly fill dest because sparse expressions have to be free of aliasing issue. 345*bf2c3715SXin Li // For non square least-square problems, b and dest might not have the same size whereas they might alias each-other. 346*bf2c3715SXin Li typename DestDerived::PlainObject tmp(cols(),rhsCols); 347*bf2c3715SXin Li ComputationInfo global_info = Success; 348*bf2c3715SXin Li for(Index k=0; k<rhsCols; ++k) 349*bf2c3715SXin Li { 350*bf2c3715SXin Li tb = b.col(k); 351*bf2c3715SXin Li tx = dest.col(k); 352*bf2c3715SXin Li derived()._solve_vector_with_guess_impl(tb,tx); 353*bf2c3715SXin Li tmp.col(k) = tx.sparseView(0); 354*bf2c3715SXin Li 355*bf2c3715SXin Li // The call to _solve_vector_with_guess_impl updates m_info, so if it failed for a previous column 356*bf2c3715SXin Li // we need to restore it to the worst value. 357*bf2c3715SXin Li if(m_info==NumericalIssue) 358*bf2c3715SXin Li global_info = NumericalIssue; 359*bf2c3715SXin Li else if(m_info==NoConvergence) 360*bf2c3715SXin Li global_info = NoConvergence; 361*bf2c3715SXin Li } 362*bf2c3715SXin Li m_info = global_info; 363*bf2c3715SXin Li dest.swap(tmp); 364*bf2c3715SXin Li } 365*bf2c3715SXin Li 366*bf2c3715SXin Li template<typename Rhs, typename DestDerived> 367*bf2c3715SXin Li typename internal::enable_if<Rhs::ColsAtCompileTime!=1 && DestDerived::ColsAtCompileTime!=1>::type _solve_with_guess_impl(const Rhs & b,MatrixBase<DestDerived> & aDest)368*bf2c3715SXin Li _solve_with_guess_impl(const Rhs& b, MatrixBase<DestDerived> &aDest) const 369*bf2c3715SXin Li { 370*bf2c3715SXin Li eigen_assert(rows()==b.rows()); 371*bf2c3715SXin Li 372*bf2c3715SXin Li Index rhsCols = b.cols(); 373*bf2c3715SXin Li DestDerived& dest(aDest.derived()); 374*bf2c3715SXin Li ComputationInfo global_info = Success; 375*bf2c3715SXin Li for(Index k=0; k<rhsCols; ++k) 376*bf2c3715SXin Li { 377*bf2c3715SXin Li typename DestDerived::ColXpr xk(dest,k); 378*bf2c3715SXin Li typename Rhs::ConstColXpr bk(b,k); 379*bf2c3715SXin Li derived()._solve_vector_with_guess_impl(bk,xk); 380*bf2c3715SXin Li 381*bf2c3715SXin Li // The call to _solve_vector_with_guess updates m_info, so if it failed for a previous column 382*bf2c3715SXin Li // we need to restore it to the worst value. 383*bf2c3715SXin Li if(m_info==NumericalIssue) 384*bf2c3715SXin Li global_info = NumericalIssue; 385*bf2c3715SXin Li else if(m_info==NoConvergence) 386*bf2c3715SXin Li global_info = NoConvergence; 387*bf2c3715SXin Li } 388*bf2c3715SXin Li m_info = global_info; 389*bf2c3715SXin Li } 390*bf2c3715SXin Li 391*bf2c3715SXin Li template<typename Rhs, typename DestDerived> 392*bf2c3715SXin Li typename internal::enable_if<Rhs::ColsAtCompileTime==1 || DestDerived::ColsAtCompileTime==1>::type _solve_with_guess_impl(const Rhs & b,MatrixBase<DestDerived> & dest)393*bf2c3715SXin Li _solve_with_guess_impl(const Rhs& b, MatrixBase<DestDerived> &dest) const 394*bf2c3715SXin Li { 395*bf2c3715SXin Li derived()._solve_vector_with_guess_impl(b,dest.derived()); 396*bf2c3715SXin Li } 397*bf2c3715SXin Li 398*bf2c3715SXin Li /** \internal default initial guess = 0 */ 399*bf2c3715SXin Li template<typename Rhs,typename Dest> _solve_impl(const Rhs & b,Dest & x)400*bf2c3715SXin Li void _solve_impl(const Rhs& b, Dest& x) const 401*bf2c3715SXin Li { 402*bf2c3715SXin Li x.setZero(); 403*bf2c3715SXin Li derived()._solve_with_guess_impl(b,x); 404*bf2c3715SXin Li } 405*bf2c3715SXin Li 406*bf2c3715SXin Li protected: init()407*bf2c3715SXin Li void init() 408*bf2c3715SXin Li { 409*bf2c3715SXin Li m_isInitialized = false; 410*bf2c3715SXin Li m_analysisIsOk = false; 411*bf2c3715SXin Li m_factorizationIsOk = false; 412*bf2c3715SXin Li m_maxIterations = -1; 413*bf2c3715SXin Li m_tolerance = NumTraits<Scalar>::epsilon(); 414*bf2c3715SXin Li } 415*bf2c3715SXin Li 416*bf2c3715SXin Li typedef internal::generic_matrix_wrapper<MatrixType> MatrixWrapper; 417*bf2c3715SXin Li typedef typename MatrixWrapper::ActualMatrixType ActualMatrixType; 418*bf2c3715SXin Li matrix()419*bf2c3715SXin Li const ActualMatrixType& matrix() const 420*bf2c3715SXin Li { 421*bf2c3715SXin Li return m_matrixWrapper.matrix(); 422*bf2c3715SXin Li } 423*bf2c3715SXin Li 424*bf2c3715SXin Li template<typename InputType> grab(const InputType & A)425*bf2c3715SXin Li void grab(const InputType &A) 426*bf2c3715SXin Li { 427*bf2c3715SXin Li m_matrixWrapper.grab(A); 428*bf2c3715SXin Li } 429*bf2c3715SXin Li 430*bf2c3715SXin Li MatrixWrapper m_matrixWrapper; 431*bf2c3715SXin Li Preconditioner m_preconditioner; 432*bf2c3715SXin Li 433*bf2c3715SXin Li Index m_maxIterations; 434*bf2c3715SXin Li RealScalar m_tolerance; 435*bf2c3715SXin Li 436*bf2c3715SXin Li mutable RealScalar m_error; 437*bf2c3715SXin Li mutable Index m_iterations; 438*bf2c3715SXin Li mutable ComputationInfo m_info; 439*bf2c3715SXin Li mutable bool m_analysisIsOk, m_factorizationIsOk; 440*bf2c3715SXin Li }; 441*bf2c3715SXin Li 442*bf2c3715SXin Li } // end namespace Eigen 443*bf2c3715SXin Li 444*bf2c3715SXin Li #endif // EIGEN_ITERATIVE_SOLVER_BASE_H 445