Eigen  3.2.8
IterativeSolverBase.h
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_ITERATIVE_SOLVER_BASE_H
00011 #define EIGEN_ITERATIVE_SOLVER_BASE_H
00012 
00013 namespace Eigen { 
00014 
00020 template< typename Derived>
00021 class IterativeSolverBase : internal::noncopyable
00022 {
00023 public:
00024   typedef typename internal::traits<Derived>::MatrixType MatrixType;
00025   typedef typename internal::traits<Derived>::Preconditioner Preconditioner;
00026   typedef typename MatrixType::Scalar Scalar;
00027   typedef typename MatrixType::Index Index;
00028   typedef typename MatrixType::RealScalar RealScalar;
00029 
00030 public:
00031 
00032   Derived& derived() { return *static_cast<Derived*>(this); }
00033   const Derived& derived() const { return *static_cast<const Derived*>(this); }
00034 
00036   IterativeSolverBase()
00037     : mp_matrix(0)
00038   {
00039     init();
00040   }
00041 
00052   template<typename InputDerived>
00053   IterativeSolverBase(const EigenBase<InputDerived>& A)
00054   {
00055     init();
00056     compute(A.derived());
00057   }
00058 
00059   ~IterativeSolverBase() {}
00060   
00066   template<typename InputDerived>
00067   Derived& analyzePattern(const EigenBase<InputDerived>& A)
00068   {
00069     grabInput(A.derived());
00070     m_preconditioner.analyzePattern(*mp_matrix);
00071     m_isInitialized = true;
00072     m_analysisIsOk = true;
00073     m_info = Success;
00074     return derived();
00075   }
00076   
00086   template<typename InputDerived>
00087   Derived& factorize(const EigenBase<InputDerived>& A)
00088   {
00089     grabInput(A.derived());
00090     eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); 
00091     m_preconditioner.factorize(*mp_matrix);
00092     m_factorizationIsOk = true;
00093     m_info = Success;
00094     return derived();
00095   }
00096 
00107   template<typename InputDerived>
00108   Derived& compute(const EigenBase<InputDerived>& A)
00109   {
00110     grabInput(A.derived());
00111     m_preconditioner.compute(*mp_matrix);
00112     m_isInitialized = true;
00113     m_analysisIsOk = true;
00114     m_factorizationIsOk = true;
00115     m_info = Success;
00116     return derived();
00117   }
00118 
00120   Index rows() const { return mp_matrix ? mp_matrix->rows() : 0; }
00122   Index cols() const { return mp_matrix ? mp_matrix->cols() : 0; }
00123 
00125   RealScalar tolerance() const { return m_tolerance; }
00126   
00128   Derived& setTolerance(const RealScalar& tolerance)
00129   {
00130     m_tolerance = tolerance;
00131     return derived();
00132   }
00133 
00135   Preconditioner& preconditioner() { return m_preconditioner; }
00136   
00138   const Preconditioner& preconditioner() const { return m_preconditioner; }
00139 
00141   int maxIterations() const
00142   {
00143     return (mp_matrix && m_maxIterations<0) ? mp_matrix->cols() : m_maxIterations;
00144   }
00145   
00147   Derived& setMaxIterations(int maxIters)
00148   {
00149     m_maxIterations = maxIters;
00150     return derived();
00151   }
00152 
00154   int iterations() const
00155   {
00156     eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
00157     return m_iterations;
00158   }
00159 
00161   RealScalar error() const
00162   {
00163     eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
00164     return m_error;
00165   }
00166 
00171   template<typename Rhs> inline const internal::solve_retval<Derived, Rhs>
00172   solve(const MatrixBase<Rhs>& b) const
00173   {
00174     eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
00175     eigen_assert(rows()==b.rows()
00176               && "IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
00177     return internal::solve_retval<Derived, Rhs>(derived(), b.derived());
00178   }
00179   
00184   template<typename Rhs>
00185   inline const internal::sparse_solve_retval<IterativeSolverBase, Rhs>
00186   solve(const SparseMatrixBase<Rhs>& b) const
00187   {
00188     eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
00189     eigen_assert(rows()==b.rows()
00190               && "IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
00191     return internal::sparse_solve_retval<IterativeSolverBase, Rhs>(*this, b.derived());
00192   }
00193 
00195   ComputationInfo info() const
00196   {
00197     eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
00198     return m_info;
00199   }
00200   
00202   template<typename Rhs, typename DestScalar, int DestOptions, typename DestIndex>
00203   void _solve_sparse(const Rhs& b, SparseMatrix<DestScalar,DestOptions,DestIndex> &dest) const
00204   {
00205     eigen_assert(rows()==b.rows());
00206     
00207     int rhsCols = b.cols();
00208     int size = b.rows();
00209     Eigen::Matrix<DestScalar,Dynamic,1> tb(size);
00210     Eigen::Matrix<DestScalar,Dynamic,1> tx(size);
00211     for(int k=0; k<rhsCols; ++k)
00212     {
00213       tb = b.col(k);
00214       tx = derived().solve(tb);
00215       dest.col(k) = tx.sparseView(0);
00216     }
00217   }
00218 
00219 protected:
00220 
00221   template<typename InputDerived>
00222   void grabInput(const EigenBase<InputDerived>& A)
00223   {
00224     // we const cast to prevent the creation of a MatrixType temporary by the compiler.
00225     grabInput_impl(A.const_cast_derived());
00226   }
00227 
00228   template<typename InputDerived>
00229   void grabInput_impl(const EigenBase<InputDerived>& A)
00230   {
00231     m_copyMatrix = A;
00232     mp_matrix = &m_copyMatrix;
00233   }
00234 
00235   void grabInput_impl(MatrixType& A)
00236   {
00237     if(MatrixType::RowsAtCompileTime==Dynamic && MatrixType::ColsAtCompileTime==Dynamic)
00238       m_copyMatrix.resize(0,0);
00239     mp_matrix = &A;
00240   }
00241 
00242   void init()
00243   {
00244     m_isInitialized = false;
00245     m_analysisIsOk = false;
00246     m_factorizationIsOk = false;
00247     m_maxIterations = -1;
00248     m_tolerance = NumTraits<Scalar>::epsilon();
00249   }
00250   MatrixType m_copyMatrix;
00251   const MatrixType* mp_matrix;
00252   Preconditioner m_preconditioner;
00253 
00254   int m_maxIterations;
00255   RealScalar m_tolerance;
00256   
00257   mutable RealScalar m_error;
00258   mutable int m_iterations;
00259   mutable ComputationInfo m_info;
00260   mutable bool m_isInitialized, m_analysisIsOk, m_factorizationIsOk;
00261 };
00262 
00263 namespace internal {
00264  
00265 template<typename Derived, typename Rhs>
00266 struct sparse_solve_retval<IterativeSolverBase<Derived>, Rhs>
00267   : sparse_solve_retval_base<IterativeSolverBase<Derived>, Rhs>
00268 {
00269   typedef IterativeSolverBase<Derived> Dec;
00270   EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
00271 
00272   template<typename Dest> void evalTo(Dest& dst) const
00273   {
00274     dec().derived()._solve_sparse(rhs(),dst);
00275   }
00276 };
00277 
00278 } // end namespace internal
00279 
00280 } // end namespace Eigen
00281 
00282 #endif // EIGEN_ITERATIVE_SOLVER_BASE_H
 All Classes Functions Variables Typedefs Enumerations Enumerator Friends