Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::DiagonalPreconditioner< _Scalar > Class Template Reference

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

#include <src/eigen/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h>

+ Inheritance diagram for Eigen::DiagonalPreconditioner< _Scalar >:
+ Collaboration 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)
 
Index rows () const
 
Index cols () const
 
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
Template Parameters
_Scalarthe 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 
42 {
45 };
@ MaxColsAtCompileTime
Definition BasicPreconditioners.h:44
@ ColsAtCompileTime
Definition BasicPreconditioners.h:43
const int Dynamic
Definition Constants.h:21

Constructor & Destructor Documentation

◆ DiagonalPreconditioner() [1/2]

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

◆ DiagonalPreconditioner() [2/2]

template<typename _Scalar >
template<typename MatType >
Eigen::DiagonalPreconditioner< _Scalar >::DiagonalPreconditioner ( const MatType &  mat)
inlineexplicit
50 : m_invdiag(mat.cols())
51 {
52 compute(mat);
53 }
DiagonalPreconditioner & compute(const MatType &mat)
Definition BasicPreconditioners.h:82
Vector m_invdiag
Definition BasicPreconditioners.h:106

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

+ Here is the call graph for this function:

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
90 {
91 x = m_invdiag.array() * b.array() ;
92 }
TCoord< P > x(const P &p)
Definition geometry_traits.hpp:297

References Eigen::DiagonalPreconditioner< _Scalar >::m_invdiag.

◆ analyzePattern()

template<typename _Scalar >
template<typename MatType >
DiagonalPreconditioner & Eigen::DiagonalPreconditioner< _Scalar >::analyzePattern ( const MatType &  )
inline
60 {
61 return *this;
62 }

◆ cols()

template<typename _Scalar >
Index Eigen::DiagonalPreconditioner< _Scalar >::cols ( ) const
inline

◆ compute()

template<typename _Scalar >
template<typename MatType >
DiagonalPreconditioner & Eigen::DiagonalPreconditioner< _Scalar >::compute ( const MatType &  mat)
inline
83 {
84 return factorize(mat);
85 }
DiagonalPreconditioner & factorize(const MatType &mat)
Definition BasicPreconditioners.h:65

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

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

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ factorize()

template<typename _Scalar >
template<typename MatType >
DiagonalPreconditioner & Eigen::DiagonalPreconditioner< _Scalar >::factorize ( const MatType &  mat)
inline
66 {
67 m_invdiag.resize(mat.cols());
68 for(int j=0; j<mat.outerSize(); ++j)
69 {
70 typename MatType::InnerIterator it(mat,j);
71 while(it && it.index()!=j) ++it;
72 if(it && it.index()==j && it.value()!=Scalar(0))
73 m_invdiag(j) = Scalar(1)/it.value();
74 else
75 m_invdiag(j) = Scalar(1);
76 }
77 m_isInitialized = true;
78 return *this;
79 }
_Scalar Scalar
Definition BasicPreconditioners.h:38
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279

References Eigen::DiagonalPreconditioner< _Scalar >::m_invdiag, Eigen::DiagonalPreconditioner< _Scalar >::m_isInitialized, and Eigen::PlainObjectBase< Derived >::resize().

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

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ info()

template<typename _Scalar >
ComputationInfo Eigen::DiagonalPreconditioner< _Scalar >::info ( )
inline
103{ return Success; }
@ Success
Definition Constants.h:432

References Eigen::Success.

◆ rows()

template<typename _Scalar >
Index Eigen::DiagonalPreconditioner< _Scalar >::rows ( ) const
inline

◆ solve()

template<typename _Scalar >
template<typename Rhs >
const Solve< DiagonalPreconditioner, Rhs > Eigen::DiagonalPreconditioner< _Scalar >::solve ( const MatrixBase< Rhs > &  b) const
inline
96 {
97 eigen_assert(m_isInitialized && "DiagonalPreconditioner is not initialized.");
98 eigen_assert(m_invdiag.size()==b.rows()
99 && "DiagonalPreconditioner::solve(): invalid number of rows of the right hand side matrix b");
100 return Solve<DiagonalPreconditioner, Rhs>(*this, b.derived());
101 }
#define eigen_assert(x)
Definition Macros.h:579

References 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: