Eigen::IterativeSolverBase< Derived > Class Template Reference

Base class for linear iterative solvers. More...

+ Inheritance diagram for Eigen::IterativeSolverBase< Derived >:

Public Types

enum  {
  ColsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef internal::traits< Derived >::MatrixType MatrixType
 
typedef internal::traits< Derived >::Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_impl (const Rhs &b, Dest &x) const
 
template<typename Rhs , typename DestDerived >
std::enable_if_t< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 > _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const
 
template<typename Rhs , typename DestDerived >
std::enable_if_t< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 > _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const
 
template<typename Rhs , typename DestDerived >
void _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
 
template<typename MatrixDerived >
Derived & analyzePattern (const EigenBase< MatrixDerived > &A)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename MatrixDerived >
Derived & compute (const EigenBase< MatrixDerived > &A)
 
Derived & derived ()
 
const Derived & derived () const
 
RealScalar error () const
 
template<typename MatrixDerived >
Derived & factorize (const EigenBase< MatrixDerived > &A)
 
ComputationInfo info () const
 
Index iterations () const
 
 IterativeSolverBase ()
 
template<typename MatrixDerived >
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 
 IterativeSolverBase (IterativeSolverBase &&)=default
 
Index maxIterations () const
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
Derived & setMaxIterations (Index maxIters)
 
Derived & setTolerance (const RealScalar &tolerance)
 
template<typename Rhs , typename Guess >
const SolveWithGuess< Derived, Rhs, Guess > solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Protected Types

typedef MatrixWrapper::ActualMatrixType ActualMatrixType
 
typedef SparseSolverBase< Derived > Base
 
typedef internal::generic_matrix_wrapper< MatrixTypeMatrixWrapper
 

Protected Member Functions

template<typename InputType >
void grab (const InputType &A)
 
void init ()
 
const ActualMatrixTypematrix () const
 

Protected Attributes

bool m_analysisIsOk
 
RealScalar m_error
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
Index m_iterations
 
MatrixWrapper m_matrixWrapper
 
Index m_maxIterations
 
Preconditioner m_preconditioner
 
RealScalar m_tolerance
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

Detailed Description

template<typename Derived>
class Eigen::IterativeSolverBase< Derived >

Base class for linear iterative solvers.

See also
class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner

Definition at line 145 of file IterativeSolverBase.h.

Member Typedef Documentation

◆ ActualMatrixType

template<typename Derived >
typedef MatrixWrapper::ActualMatrixType Eigen::IterativeSolverBase< Derived >::ActualMatrixType
protected

Definition at line 422 of file IterativeSolverBase.h.

◆ Base

template<typename Derived >
typedef SparseSolverBase<Derived> Eigen::IterativeSolverBase< Derived >::Base
protected

Definition at line 148 of file IterativeSolverBase.h.

◆ MatrixType

template<typename Derived >
typedef internal::traits<Derived>::MatrixType Eigen::IterativeSolverBase< Derived >::MatrixType

Definition at line 152 of file IterativeSolverBase.h.

◆ MatrixWrapper

template<typename Derived >
typedef internal::generic_matrix_wrapper<MatrixType> Eigen::IterativeSolverBase< Derived >::MatrixWrapper
protected

Definition at line 421 of file IterativeSolverBase.h.

◆ Preconditioner

template<typename Derived >
typedef internal::traits<Derived>::Preconditioner Eigen::IterativeSolverBase< Derived >::Preconditioner

Definition at line 153 of file IterativeSolverBase.h.

◆ RealScalar

template<typename Derived >
typedef MatrixType::RealScalar Eigen::IterativeSolverBase< Derived >::RealScalar

Definition at line 156 of file IterativeSolverBase.h.

◆ Scalar

template<typename Derived >
typedef MatrixType::Scalar Eigen::IterativeSolverBase< Derived >::Scalar

Definition at line 154 of file IterativeSolverBase.h.

◆ StorageIndex

template<typename Derived >
typedef MatrixType::StorageIndex Eigen::IterativeSolverBase< Derived >::StorageIndex

Definition at line 155 of file IterativeSolverBase.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Derived >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 158 of file IterativeSolverBase.h.

158  {
159  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
160  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
161  };

Constructor & Destructor Documentation

◆ IterativeSolverBase() [1/3]

template<typename Derived >
Eigen::IterativeSolverBase< Derived >::IterativeSolverBase ( )
inline

Default constructor.

Definition at line 168 of file IterativeSolverBase.h.

169  {
170  init();
171  }

◆ IterativeSolverBase() [2/3]

template<typename Derived >
template<typename MatrixDerived >
Eigen::IterativeSolverBase< Derived >::IterativeSolverBase ( const EigenBase< MatrixDerived > &  A)
inlineexplicit

Initialize the solver with matrix A for further Ax=b solving.

This constructor is a shortcut for the default constructor followed by a call to compute().

Warning
this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.

Definition at line 184 of file IterativeSolverBase.h.

185  : m_matrixWrapper(A.derived())
186  {
187  init();
188  compute(matrix());
189  }
MatrixXcf A
Derived & compute(const EigenBase< MatrixDerived > &A)
const ActualMatrixType & matrix() const

◆ IterativeSolverBase() [3/3]

template<typename Derived >
Eigen::IterativeSolverBase< Derived >::IterativeSolverBase ( IterativeSolverBase< Derived > &&  )
default

◆ ~IterativeSolverBase()

template<typename Derived >
Eigen::IterativeSolverBase< Derived >::~IterativeSolverBase ( )
inline

Definition at line 194 of file IterativeSolverBase.h.

194 {}

Member Function Documentation

◆ _solve_impl()

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::IterativeSolverBase< Derived >::_solve_impl ( const Rhs &  b,
Dest &  x 
) const
inline

Definition at line 405 of file IterativeSolverBase.h.

406  {
407  x.setZero();
408  derived()._solve_with_guess_impl(b,x);
409  }
Array< int, 3, 1 > b

◆ _solve_with_guess_impl() [1/3]

template<typename Derived >
template<typename Rhs , typename DestDerived >
std::enable_if_t<Rhs::ColsAtCompileTime!=1 && DestDerived::ColsAtCompileTime!=1> Eigen::IterativeSolverBase< Derived >::_solve_with_guess_impl ( const Rhs &  b,
MatrixBase< DestDerived > &  aDest 
) const
inline

Definition at line 373 of file IterativeSolverBase.h.

374  {
375  eigen_assert(rows()==b.rows());
376 
377  Index rhsCols = b.cols();
378  DestDerived& dest(aDest.derived());
379  ComputationInfo global_info = Success;
380  for(Index k=0; k<rhsCols; ++k)
381  {
382  typename DestDerived::ColXpr xk(dest,k);
383  typename Rhs::ConstColXpr bk(b,k);
384  derived()._solve_vector_with_guess_impl(bk,xk);
385 
386  // The call to _solve_vector_with_guess updates m_info, so if it failed for a previous column
387  // we need to restore it to the worst value.
389  global_info = NumericalIssue;
390  else if(m_info==NoConvergence)
391  global_info = NoConvergence;
392  }
393  m_info = global_info;
394  }
#define eigen_assert(x)
Definition: Macros.h:902
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
ComputationInfo
Definition: Constants.h:444
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82

◆ _solve_with_guess_impl() [2/3]

template<typename Derived >
template<typename Rhs , typename DestDerived >
std::enable_if_t<Rhs::ColsAtCompileTime==1 || DestDerived::ColsAtCompileTime==1> Eigen::IterativeSolverBase< Derived >::_solve_with_guess_impl ( const Rhs &  b,
MatrixBase< DestDerived > &  dest 
) const
inline

Definition at line 398 of file IterativeSolverBase.h.

399  {
400  derived()._solve_vector_with_guess_impl(b,dest.derived());
401  }

◆ _solve_with_guess_impl() [3/3]

template<typename Derived >
template<typename Rhs , typename DestDerived >
void Eigen::IterativeSolverBase< Derived >::_solve_with_guess_impl ( const Rhs &  b,
SparseMatrixBase< DestDerived > &  aDest 
) const
inline

Definition at line 339 of file IterativeSolverBase.h.

340  {
341  eigen_assert(rows()==b.rows());
342 
343  Index rhsCols = b.cols();
344  Index size = b.rows();
345  DestDerived& dest(aDest.derived());
346  typedef typename DestDerived::Scalar DestScalar;
349  // We do not directly fill dest because sparse expressions have to be free of aliasing issue.
350  // For non square least-square problems, b and dest might not have the same size whereas they might alias each-other.
351  typename DestDerived::PlainObject tmp(cols(),rhsCols);
352  ComputationInfo global_info = Success;
353  for(Index k=0; k<rhsCols; ++k)
354  {
355  tb = b.col(k);
356  tx = dest.col(k);
357  derived()._solve_vector_with_guess_impl(tb,tx);
358  tmp.col(k) = tx.sparseView(0);
359 
360  // The call to _solve_vector_with_guess_impl updates m_info, so if it failed for a previous column
361  // we need to restore it to the worst value.
363  global_info = NumericalIssue;
364  else if(m_info==NoConvergence)
365  global_info = NoConvergence;
366  }
367  m_info = global_info;
368  dest.swap(tmp);
369  }
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182

◆ analyzePattern()

template<typename Derived >
template<typename MatrixDerived >
Derived& Eigen::IterativeSolverBase< Derived >::analyzePattern ( const EigenBase< MatrixDerived > &  A)
inline

Initializes the iterative solver for the sparsity pattern of the matrix A for further solving Ax=b problems.

Currently, this function mostly calls analyzePattern on the preconditioner. In the future we might, for instance, implement column reordering for faster matrix vector products.

Definition at line 202 of file IterativeSolverBase.h.

203  {
204  grab(A.derived());
205  m_preconditioner.analyzePattern(matrix());
206  m_isInitialized = true;
207  m_analysisIsOk = true;
208  m_info = m_preconditioner.info();
209  return derived();
210  }
void grab(const InputType &A)

◆ cols()

template<typename Derived >
EIGEN_CONSTEXPR Index Eigen::IterativeSolverBase< Derived >::cols ( void  ) const
inline

Definition at line 258 of file IterativeSolverBase.h.

258 { return matrix().cols(); }

◆ compute()

template<typename Derived >
template<typename MatrixDerived >
Derived& Eigen::IterativeSolverBase< Derived >::compute ( const EigenBase< MatrixDerived > &  A)
inline

Initializes the iterative solver with the matrix A for further solving Ax=b problems.

Currently, this function mostly initializes/computes the preconditioner. In the future we might, for instance, implement column reordering for faster matrix vector products.

Warning
this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.

Definition at line 243 of file IterativeSolverBase.h.

244  {
245  grab(A.derived());
246  m_preconditioner.compute(matrix());
247  m_isInitialized = true;
248  m_analysisIsOk = true;
249  m_factorizationIsOk = true;
250  m_info = m_preconditioner.info();
251  return derived();
252  }

◆ derived() [1/2]

template<typename Derived >
Derived& Eigen::SparseSolverBase< Derived >::derived
inline

Definition at line 83 of file SparseSolverBase.h.

83 { return *static_cast<Derived*>(this); }

◆ derived() [2/2]

template<typename Derived >
const Derived& Eigen::SparseSolverBase< Derived >::derived
inline

Definition at line 84 of file SparseSolverBase.h.

84 { return *static_cast<const Derived*>(this); }

◆ error()

template<typename Derived >
RealScalar Eigen::IterativeSolverBase< Derived >::error ( ) const
inline
Returns
the tolerance error reached during the last solve. It is a close approximation of the true relative residual error |Ax-b|/|b|.

Definition at line 310 of file IterativeSolverBase.h.

311  {
312  eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
313  return m_error;
314  }

◆ factorize()

template<typename Derived >
template<typename MatrixDerived >
Derived& Eigen::IterativeSolverBase< Derived >::factorize ( const EigenBase< MatrixDerived > &  A)
inline

Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems.

Currently, this function mostly calls factorize on the preconditioner.

Warning
this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.

Definition at line 222 of file IterativeSolverBase.h.

223  {
224  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
225  grab(A.derived());
226  m_preconditioner.factorize(matrix());
227  m_factorizationIsOk = true;
228  m_info = m_preconditioner.info();
229  return derived();
230  }

◆ grab()

template<typename Derived >
template<typename InputType >
void Eigen::IterativeSolverBase< Derived >::grab ( const InputType &  A)
inlineprotected

Definition at line 430 of file IterativeSolverBase.h.

431  {
432  m_matrixWrapper.grab(A);
433  }

◆ info()

template<typename Derived >
ComputationInfo Eigen::IterativeSolverBase< Derived >::info ( ) const
inline
Returns
Success if the iterations converged, and NoConvergence otherwise.

Definition at line 331 of file IterativeSolverBase.h.

332  {
333  eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
334  return m_info;
335  }

◆ init()

template<typename Derived >
void Eigen::IterativeSolverBase< Derived >::init ( )
inlineprotected

Definition at line 412 of file IterativeSolverBase.h.

413  {
414  m_isInitialized = false;
415  m_analysisIsOk = false;
416  m_factorizationIsOk = false;
417  m_maxIterations = -1;
418  m_tolerance = NumTraits<Scalar>::epsilon();
419  }

◆ iterations()

template<typename Derived >
Index Eigen::IterativeSolverBase< Derived >::iterations ( ) const
inline
Returns
the number of iterations performed during the last solve

Definition at line 301 of file IterativeSolverBase.h.

302  {
303  eigen_assert(m_isInitialized && "IterativeSolverBase is not initialized.");
304  return m_iterations;
305  }

◆ matrix()

template<typename Derived >
const ActualMatrixType& Eigen::IterativeSolverBase< Derived >::matrix ( ) const
inlineprotected

Definition at line 424 of file IterativeSolverBase.h.

425  {
426  return m_matrixWrapper.matrix();
427  }

◆ maxIterations()

template<typename Derived >
Index Eigen::IterativeSolverBase< Derived >::maxIterations ( ) const
inline
Returns
the max number of iterations. It is either the value set by setMaxIterations or, by default, twice the number of columns of the matrix.

Definition at line 286 of file IterativeSolverBase.h.

287  {
288  return (m_maxIterations<0) ? 2*matrix().cols() : m_maxIterations;
289  }

◆ preconditioner() [1/2]

template<typename Derived >
Preconditioner& Eigen::IterativeSolverBase< Derived >::preconditioner ( )
inline
Returns
a read-write reference to the preconditioner for custom configuration.

Definition at line 277 of file IterativeSolverBase.h.

277 { return m_preconditioner; }

◆ preconditioner() [2/2]

template<typename Derived >
const Preconditioner& Eigen::IterativeSolverBase< Derived >::preconditioner ( ) const
inline
Returns
a read-only reference to the preconditioner.

Definition at line 280 of file IterativeSolverBase.h.

280 { return m_preconditioner; }

◆ rows()

template<typename Derived >
EIGEN_CONSTEXPR Index Eigen::IterativeSolverBase< Derived >::rows ( void  ) const
inline

Definition at line 255 of file IterativeSolverBase.h.

255 { return matrix().rows(); }

◆ setMaxIterations()

template<typename Derived >
Derived& Eigen::IterativeSolverBase< Derived >::setMaxIterations ( Index  maxIters)
inline

Sets the max number of iterations. Default is twice the number of columns of the matrix.

Definition at line 294 of file IterativeSolverBase.h.

295  {
296  m_maxIterations = maxIters;
297  return derived();
298  }

◆ setTolerance()

template<typename Derived >
Derived& Eigen::IterativeSolverBase< Derived >::setTolerance ( const RealScalar tolerance)
inline

Sets the tolerance threshold used by the stopping criteria.

This value is used as an upper bound to the relative residual error: |Ax-b|/|b|. The default value is the machine precision given by NumTraits<Scalar>::epsilon()

Definition at line 270 of file IterativeSolverBase.h.

271  {
273  return derived();
274  }

◆ solveWithGuess()

template<typename Derived >
template<typename Rhs , typename Guess >
const SolveWithGuess<Derived, Rhs, Guess> Eigen::IterativeSolverBase< Derived >::solveWithGuess ( const MatrixBase< Rhs > &  b,
const Guess &  x0 
) const
inline
Returns
the solution x of \( A x = b \) using the current decomposition of A and x0 as an initial solution.
See also
solve(), compute()

Definition at line 323 of file IterativeSolverBase.h.

324  {
325  eigen_assert(m_isInitialized && "Solver is not initialized.");
326  eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
327  return SolveWithGuess<Derived, Rhs, Guess>(derived(), b.derived(), x0);
328  }

◆ tolerance()

template<typename Derived >
RealScalar Eigen::IterativeSolverBase< Derived >::tolerance ( ) const
inline
Returns
the tolerance threshold used by the stopping criteria.
See also
setTolerance()

Definition at line 263 of file IterativeSolverBase.h.

263 { return m_tolerance; }

Member Data Documentation

◆ m_analysisIsOk

template<typename Derived >
bool Eigen::IterativeSolverBase< Derived >::m_analysisIsOk
mutableprotected

Definition at line 444 of file IterativeSolverBase.h.

◆ m_error

template<typename Derived >
RealScalar Eigen::IterativeSolverBase< Derived >::m_error
mutableprotected

Definition at line 441 of file IterativeSolverBase.h.

◆ m_factorizationIsOk

template<typename Derived >
bool Eigen::IterativeSolverBase< Derived >::m_factorizationIsOk
protected

Definition at line 444 of file IterativeSolverBase.h.

◆ m_info

template<typename Derived >
ComputationInfo Eigen::IterativeSolverBase< Derived >::m_info
mutableprotected

Definition at line 443 of file IterativeSolverBase.h.

◆ m_isInitialized

template<typename Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

Definition at line 123 of file SparseSolverBase.h.

◆ m_iterations

template<typename Derived >
Index Eigen::IterativeSolverBase< Derived >::m_iterations
mutableprotected

Definition at line 442 of file IterativeSolverBase.h.

◆ m_matrixWrapper

template<typename Derived >
MatrixWrapper Eigen::IterativeSolverBase< Derived >::m_matrixWrapper
protected

Definition at line 435 of file IterativeSolverBase.h.

◆ m_maxIterations

template<typename Derived >
Index Eigen::IterativeSolverBase< Derived >::m_maxIterations
protected

Definition at line 438 of file IterativeSolverBase.h.

◆ m_preconditioner

template<typename Derived >
Preconditioner Eigen::IterativeSolverBase< Derived >::m_preconditioner
protected

Definition at line 436 of file IterativeSolverBase.h.

◆ m_tolerance

template<typename Derived >
RealScalar Eigen::IterativeSolverBase< Derived >::m_tolerance
protected

Definition at line 439 of file IterativeSolverBase.h.


The documentation for this class was generated from the following file: