Eigen::DiagonalPreconditioner< Scalar_ > Class Template Reference

A preconditioner based on the digonal entries. More...

+ Inheritance diagram for Eigen::DiagonalPreconditioner< Scalar_ >:

Public Types

enum  {
  ColsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef Vector::StorageIndex StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_impl (const Rhs &b, Dest &x) const
 
template<typename MatType >
DiagonalPreconditioneranalyzePattern (const MatType &)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename MatType >
DiagonalPreconditionercompute (const MatType &mat)
 
 DiagonalPreconditioner ()
 
template<typename MatType >
 DiagonalPreconditioner (const MatType &mat)
 
template<typename MatType >
DiagonalPreconditionerfactorize (const MatType &mat)
 
ComputationInfo info ()
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
template<typename Rhs >
const Solve< DiagonalPreconditioner, Rhs > solve (const MatrixBase< Rhs > &b) const
 

Protected Attributes

Vector m_invdiag
 
bool m_isInitialized
 

Private Types

typedef Scalar_ Scalar
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 

Detailed Description

template<typename Scalar_>
class Eigen::DiagonalPreconditioner< Scalar_ >

A preconditioner based on the digonal entries.

This class allows to approximately solve for A.x = b problems assuming A is a diagonal matrix. In other words, this preconditioner neglects all off diagonal entries and, in Eigen's language, solves for:

A.diagonal().asDiagonal() . x = b
Array< int, 3, 1 > b
MatrixXcf A
Template Parameters
Scalar_the type of the scalar.

This class follows the sparse solver concept .

This preconditioner is suitable for both selfadjoint and general problems. The diagonal entries are pre-inverted and stored into a dense vector.

Note
A variant that has yet to be implemented would attempt to preserve the norm of each column.
See also
class LeastSquareDiagonalPreconditioner, class ConjugateGradient

Definition at line 38 of file BasicPreconditioners.h.

Member Typedef Documentation

◆ Scalar

template<typename Scalar_ >
typedef Scalar_ Eigen::DiagonalPreconditioner< Scalar_ >::Scalar
private

Definition at line 40 of file BasicPreconditioners.h.

◆ StorageIndex

template<typename Scalar_ >
typedef Vector::StorageIndex Eigen::DiagonalPreconditioner< Scalar_ >::StorageIndex

Definition at line 43 of file BasicPreconditioners.h.

◆ Vector

template<typename Scalar_ >
typedef Matrix<Scalar,Dynamic,1> Eigen::DiagonalPreconditioner< Scalar_ >::Vector
private

Definition at line 41 of file BasicPreconditioners.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar_ >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 44 of file BasicPreconditioners.h.

Constructor & Destructor Documentation

◆ DiagonalPreconditioner() [1/2]

template<typename Scalar_ >
Eigen::DiagonalPreconditioner< Scalar_ >::DiagonalPreconditioner ( )
inline

Definition at line 49 of file BasicPreconditioners.h.

◆ DiagonalPreconditioner() [2/2]

template<typename Scalar_ >
template<typename MatType >
Eigen::DiagonalPreconditioner< Scalar_ >::DiagonalPreconditioner ( const MatType &  mat)
inlineexplicit

Definition at line 52 of file BasicPreconditioners.h.

52  : m_invdiag(mat.cols())
53  {
54  compute(mat);
55  }
DiagonalPreconditioner & compute(const MatType &mat)

Member Function Documentation

◆ _solve_impl()

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

Definition at line 91 of file BasicPreconditioners.h.

92  {
93  x = m_invdiag.array() * b.array() ;
94  }

◆ analyzePattern()

template<typename Scalar_ >
template<typename MatType >
DiagonalPreconditioner& Eigen::DiagonalPreconditioner< Scalar_ >::analyzePattern ( const MatType &  )
inline

Definition at line 61 of file BasicPreconditioners.h.

62  {
63  return *this;
64  }

◆ cols()

template<typename Scalar_ >
EIGEN_CONSTEXPR Index Eigen::DiagonalPreconditioner< Scalar_ >::cols ( void  ) const
inline

Definition at line 58 of file BasicPreconditioners.h.

58 { return m_invdiag.size(); }

◆ compute()

template<typename Scalar_ >
template<typename MatType >
DiagonalPreconditioner& Eigen::DiagonalPreconditioner< Scalar_ >::compute ( const MatType &  mat)
inline

Definition at line 84 of file BasicPreconditioners.h.

85  {
86  return factorize(mat);
87  }
DiagonalPreconditioner & factorize(const MatType &mat)

◆ factorize()

template<typename Scalar_ >
template<typename MatType >
DiagonalPreconditioner& Eigen::DiagonalPreconditioner< Scalar_ >::factorize ( const MatType &  mat)
inline

Definition at line 67 of file BasicPreconditioners.h.

68  {
69  m_invdiag.resize(mat.cols());
70  for(int j=0; j<mat.outerSize(); ++j)
71  {
72  typename MatType::InnerIterator it(mat,j);
73  while(it && it.index()!=j) ++it;
74  if(it && it.index()==j && it.value()!=Scalar(0))
75  m_invdiag(j) = Scalar(1)/it.value();
76  else
77  m_invdiag(j) = Scalar(1);
78  }
79  m_isInitialized = true;
80  return *this;
81  }
constexpr void resize(Index rows, Index cols)
std::ptrdiff_t j

◆ info()

template<typename Scalar_ >
ComputationInfo Eigen::DiagonalPreconditioner< Scalar_ >::info ( )
inline

Definition at line 105 of file BasicPreconditioners.h.

105 { return Success; }
@ Success
Definition: Constants.h:446

◆ rows()

template<typename Scalar_ >
EIGEN_CONSTEXPR Index Eigen::DiagonalPreconditioner< Scalar_ >::rows ( void  ) const
inline

Definition at line 57 of file BasicPreconditioners.h.

57 { return m_invdiag.size(); }

◆ solve()

template<typename Scalar_ >
template<typename Rhs >
const Solve<DiagonalPreconditioner, Rhs> Eigen::DiagonalPreconditioner< Scalar_ >::solve ( const MatrixBase< Rhs > &  b) const
inline

Definition at line 97 of file BasicPreconditioners.h.

98  {
99  eigen_assert(m_isInitialized && "DiagonalPreconditioner is not initialized.");
100  eigen_assert(m_invdiag.size()==b.rows()
101  && "DiagonalPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
102  return Solve<DiagonalPreconditioner, Rhs>(*this, b.derived());
103  }
#define eigen_assert(x)
Definition: Macros.h:902

Member Data Documentation

◆ m_invdiag

template<typename Scalar_ >
Vector Eigen::DiagonalPreconditioner< Scalar_ >::m_invdiag
protected

Definition at line 108 of file BasicPreconditioners.h.

◆ m_isInitialized

template<typename Scalar_ >
bool Eigen::DiagonalPreconditioner< Scalar_ >::m_isInitialized
protected

Definition at line 109 of file BasicPreconditioners.h.


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