Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
IterativeSolverBase.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011 Gael Guennebaud <[email protected]>
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_ITERATIVE_SOLVER_BASE_H
11 #define EIGEN_ITERATIVE_SOLVER_BASE_H
12 
13 namespace Eigen {
14 
20 template< typename Derived>
21 class IterativeSolverBase : internal::noncopyable
22 {
23 public:
24  typedef typename internal::traits<Derived>::MatrixType MatrixType;
25  typedef typename internal::traits<Derived>::Preconditioner Preconditioner;
26  typedef typename MatrixType::Scalar Scalar;
27  typedef typename MatrixType::Index Index;
28  typedef typename MatrixType::RealScalar RealScalar;
29 
30 public:
31 
32  Derived& derived() { return *static_cast<Derived*>(this); }
33  const Derived& derived() const { return *static_cast<const Derived*>(this); }
34 
37  : mp_matrix(0)
38  {
39  init();
40  }
41 
52  template<typename InputDerived>
54  {
55  init();
56  compute(A.derived());
57  }
58 
60 
66  template<typename InputDerived>
68  {
69  grabInput(A.derived());
70  m_preconditioner.analyzePattern(*mp_matrix);
71  m_isInitialized = true;
72  m_analysisIsOk = true;
73  m_info = Success;
74  return derived();
75  }
76 
86  template<typename InputDerived>
88  {
89  grabInput(A.derived());
90  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
91  m_preconditioner.factorize(*mp_matrix);
92  m_factorizationIsOk = true;
93  m_info = Success;
94  return derived();
95  }
96 
107  template<typename InputDerived>
108  Derived& compute(const EigenBase<InputDerived>& A)
109  {
110  grabInput(A.derived());
111  m_preconditioner.compute(*mp_matrix);
112  m_isInitialized = true;
113  m_analysisIsOk = true;
114  m_factorizationIsOk = true;
115  m_info = Success;
116  return derived();
117  }
118 
120  Index rows() const { return mp_matrix ? mp_matrix->rows() : 0; }
122  Index cols() const { return mp_matrix ? mp_matrix->cols() : 0; }
123 
125  RealScalar tolerance() const { return m_tolerance; }
126 
128  Derived& setTolerance(const RealScalar& tolerance)
129  {
130  m_tolerance = tolerance;
131  return derived();
132  }
133 
135  Preconditioner& preconditioner() { return m_preconditioner; }
136 
138  const Preconditioner& preconditioner() const { return m_preconditioner; }
139 
141  int maxIterations() const
142  {
143  return (mp_matrix && m_maxIterations<0) ? mp_matrix->cols() : m_maxIterations;
144  }
145 
147  Derived& setMaxIterations(int maxIters)
148  {
149  m_maxIterations = maxIters;
150  return derived();
151  }
152 
154  int iterations() const
155  {
156  eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
157  return m_iterations;
158  }
159 
161  RealScalar error() const
162  {
163  eigen_assert(m_isInitialized && "ConjugateGradient is not initialized.");
164  return m_error;
165  }
166 
171  template<typename Rhs> inline const internal::solve_retval<Derived, Rhs>
172  solve(const MatrixBase<Rhs>& b) const
173  {
174  eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
175  eigen_assert(rows()==b.rows()
176  && "IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
177  return internal::solve_retval<Derived, Rhs>(derived(), b.derived());
178  }
179 
184  template<typename Rhs>
185  inline const internal::sparse_solve_retval<IterativeSolverBase, Rhs>
186  solve(const SparseMatrixBase<Rhs>& b) const
187  {
188  eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
189  eigen_assert(rows()==b.rows()
190  && "IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
191  return internal::sparse_solve_retval<IterativeSolverBase, Rhs>(*this, b.derived());
192  }
193 
196  {
197  eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
198  return m_info;
199  }
200 
202  template<typename Rhs, typename DestScalar, int DestOptions, typename DestIndex>
203  void _solve_sparse(const Rhs& b, SparseMatrix<DestScalar,DestOptions,DestIndex> &dest) const
204  {
205  eigen_assert(rows()==b.rows());
206 
207  int rhsCols = b.cols();
208  int size = b.rows();
211  for(int k=0; k<rhsCols; ++k)
212  {
213  tb = b.col(k);
214  tx = derived().solve(tb);
215  dest.col(k) = tx.sparseView(0);
216  }
217  }
218 
219 protected:
220 
221  template<typename InputDerived>
222  void grabInput(const EigenBase<InputDerived>& A)
223  {
224  // we const cast to prevent the creation of a MatrixType temporary by the compiler.
225  grabInput_impl(A.const_cast_derived());
226  }
227 
228  template<typename InputDerived>
229  void grabInput_impl(const EigenBase<InputDerived>& A)
230  {
231  m_copyMatrix = A;
232  mp_matrix = &m_copyMatrix;
233  }
234 
235  void grabInput_impl(MatrixType& A)
236  {
237  if(MatrixType::RowsAtCompileTime==Dynamic && MatrixType::ColsAtCompileTime==Dynamic)
238  m_copyMatrix.resize(0,0);
239  mp_matrix = &A;
240  }
241 
242  void init()
243  {
244  m_isInitialized = false;
245  m_analysisIsOk = false;
246  m_factorizationIsOk = false;
247  m_maxIterations = -1;
248  m_tolerance = NumTraits<Scalar>::epsilon();
249  }
250  MatrixType m_copyMatrix;
251  const MatrixType* mp_matrix;
252  Preconditioner m_preconditioner;
253 
254  int m_maxIterations;
255  RealScalar m_tolerance;
256 
257  mutable RealScalar m_error;
258  mutable int m_iterations;
259  mutable ComputationInfo m_info;
260  mutable bool m_isInitialized, m_analysisIsOk, m_factorizationIsOk;
261 };
262 
263 namespace internal {
264 
265 template<typename Derived, typename Rhs>
266 struct sparse_solve_retval<IterativeSolverBase<Derived>, Rhs>
267  : sparse_solve_retval_base<IterativeSolverBase<Derived>, Rhs>
268 {
269  typedef IterativeSolverBase<Derived> Dec;
270  EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
271 
272  template<typename Dest> void evalTo(Dest& dst) const
273  {
274  dec().derived()._solve_sparse(rhs(),dst);
275  }
276 };
277 
278 } // end namespace internal
279 
280 } // end namespace Eigen
281 
282 #endif // EIGEN_ITERATIVE_SOLVER_BASE_H
ComputationInfo info() const
Definition: IterativeSolverBase.h:195
Preconditioner & preconditioner()
Definition: IterativeSolverBase.h:135
A versatible sparse matrix representation.
Definition: SparseMatrix.h:85
Derived & analyzePattern(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:67
const int Dynamic
Definition: Constants.h:21
Definition: EigenBase.h:26
const internal::sparse_solve_retval< IterativeSolverBase, Rhs > solve(const SparseMatrixBase< Rhs > &b) const
Definition: IterativeSolverBase.h:186
IterativeSolverBase()
Definition: IterativeSolverBase.h:36
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:239
Derived & derived()
Definition: EigenBase.h:34
Derived & factorize(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:87
Derived & setTolerance(const RealScalar &tolerance)
Definition: IterativeSolverBase.h:128
RealScalar error() const
Definition: IterativeSolverBase.h:161
IterativeSolverBase(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:53
RealScalar tolerance() const
Definition: IterativeSolverBase.h:125
const Preconditioner & preconditioner() const
Definition: IterativeSolverBase.h:138
Definition: Constants.h:376
const internal::solve_retval< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: IterativeSolverBase.h:172
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
Index rows() const
Definition: SparseMatrixBase.h:159
ComputationInfo
Definition: Constants.h:374
Derived & setMaxIterations(int maxIters)
Definition: IterativeSolverBase.h:147
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:21
int iterations() const
Definition: IterativeSolverBase.h:154
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
ColXpr col(Index i)
Definition: SparseMatrixBase.h:733
int maxIterations() const
Definition: IterativeSolverBase.h:141
Derived & compute(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:108