Eigen::DiagonalPreconditioner< Scalar_ > Class Template Reference

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

#include <BasicPreconditioners.h>

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

Public Types

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

Public Member Functions

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

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
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
list x
Definition: plotDoE.py:28
Template Parameters
Scalar_the type of the scalar.

\implsparsesolverconcept

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

Member Typedef Documentation

◆ Scalar

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

◆ StorageIndex

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

◆ Vector

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

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar_ >
anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 
@ MaxColsAtCompileTime
Definition: BasicPreconditioners.h:45
@ ColsAtCompileTime
Definition: BasicPreconditioners.h:45
const int Dynamic
Definition: Constants.h:25

Constructor & Destructor Documentation

◆ DiagonalPreconditioner() [1/2]

template<typename Scalar_ >
Eigen::DiagonalPreconditioner< Scalar_ >::DiagonalPreconditioner ( )
inline
47 : m_isInitialized(false) {}
bool m_isInitialized
Definition: BasicPreconditioners.h:100

◆ DiagonalPreconditioner() [2/2]

template<typename Scalar_ >
template<typename MatType >
Eigen::DiagonalPreconditioner< Scalar_ >::DiagonalPreconditioner ( const MatType &  mat)
inlineexplicit
50  : m_invdiag(mat.cols()) {
51  compute(mat);
52  }
Vector m_invdiag
Definition: BasicPreconditioners.h:99
DiagonalPreconditioner & compute(const MatType &mat)
Definition: BasicPreconditioners.h:78
Index cols() const
Definition: SparseMatrix.h:161

References Eigen::DiagonalPreconditioner< Scalar_ >::compute().

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
84  {
85  x = m_invdiag.array() * b.array();
86  }

References b, Eigen::DiagonalPreconditioner< Scalar_ >::m_invdiag, and plotDoE::x.

◆ analyzePattern()

template<typename Scalar_ >
template<typename MatType >
DiagonalPreconditioner& Eigen::DiagonalPreconditioner< Scalar_ >::analyzePattern ( const MatType &  )
inline
58  {
59  return *this;
60  }

◆ cols()

◆ compute()

template<typename Scalar_ >
template<typename MatType >
DiagonalPreconditioner& Eigen::DiagonalPreconditioner< Scalar_ >::compute ( const MatType &  mat)
inline
78  {
79  return factorize(mat);
80  }
DiagonalPreconditioner & factorize(const MatType &mat)
Definition: BasicPreconditioners.h:63

References Eigen::DiagonalPreconditioner< Scalar_ >::factorize().

Referenced by Eigen::DiagonalPreconditioner< Scalar_ >::DiagonalPreconditioner().

◆ factorize()

template<typename Scalar_ >
template<typename MatType >
DiagonalPreconditioner& Eigen::DiagonalPreconditioner< Scalar_ >::factorize ( const MatType &  mat)
inline
63  {
65  for (int j = 0; j < mat.outerSize(); ++j) {
66  typename MatType::InnerIterator it(mat, j);
67  while (it && it.index() != j) ++it;
68  if (it && it.index() == j && it.value() != Scalar(0))
69  m_invdiag(j) = Scalar(1) / it.value();
70  else
71  m_invdiag(j) = Scalar(1);
72  }
73  m_isInitialized = true;
74  return *this;
75  }
SCALAR Scalar
Definition: bench_gemm.cpp:45
Scalar_ Scalar
Definition: BasicPreconditioners.h:40
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
Index outerSize() const
Definition: SparseMatrix.h:166
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), j, Eigen::DiagonalPreconditioner< Scalar_ >::m_invdiag, Eigen::DiagonalPreconditioner< Scalar_ >::m_isInitialized, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::outerSize(), and Eigen::PlainObjectBase< Derived >::resize().

Referenced by Eigen::DiagonalPreconditioner< Scalar_ >::compute().

◆ info()

template<typename Scalar_ >
ComputationInfo Eigen::DiagonalPreconditioner< Scalar_ >::info ( )
inline
96 { return Success; }
@ Success
Definition: Constants.h:440

References Eigen::Success.

◆ rows()

◆ solve()

template<typename Scalar_ >
template<typename Rhs >
const Solve<DiagonalPreconditioner, Rhs> Eigen::DiagonalPreconditioner< Scalar_ >::solve ( const MatrixBase< Rhs > &  b) const
inline
89  {
90  eigen_assert(m_isInitialized && "DiagonalPreconditioner is not initialized.");
91  eigen_assert(m_invdiag.size() == b.rows() &&
92  "DiagonalPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
93  return Solve<DiagonalPreconditioner, Rhs>(*this, b.derived());
94  }
#define eigen_assert(x)
Definition: Macros.h:910

References b, eigen_assert, Eigen::DiagonalPreconditioner< Scalar_ >::m_invdiag, and Eigen::DiagonalPreconditioner< Scalar_ >::m_isInitialized.

Member Data Documentation

◆ m_invdiag

◆ m_isInitialized


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