Public Types | Public Member Functions | Static Protected Member Functions | Protected Attributes | Friends | List of all members
Eigen::LDLT< _MatrixType, _UpLo > Class Template Reference

Robust Cholesky decomposition of a matrix with pivoting. More...

#include <LDLT.h>

Inheritance diagram for Eigen::LDLT< _MatrixType, _UpLo >:
Inheritance graph
[legend]

Public Types

enum  { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, UpLo = _UpLo }
 
typedef SolverBase< LDLTBase
 
typedef _MatrixType MatrixType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTimePermutationType
 
typedef Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
 
typedef internal::LDLT_Traits< MatrixType, UpLoTraits
 
typedef Transpositions< RowsAtCompileTime, MaxRowsAtCompileTimeTranspositionType
 
- Public Types inherited from Eigen::SolverBase< LDLT< _MatrixType, _UpLo > >
enum  
 
typedef internal::conditional< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, ConstTransposeReturnType >, ConstTransposeReturnType >::type AdjointReturnType
 
typedef EigenBase< LDLT< _MatrixType, _UpLo > > Base
 
typedef Scalar CoeffReturnType
 
typedef internal::add_const< Transpose< const LDLT< _MatrixType, _UpLo > > >::type ConstTransposeReturnType
 
typedef internal::traits< LDLT< _MatrixType, _UpLo > >::Scalar Scalar
 
- Public Types inherited from Eigen::EigenBase< LDLT< _MatrixType, _UpLo > >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< LDLT< _MatrixType, _UpLo > >::StorageKind StorageKind
 

Public Member Functions

template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
const LDLTadjoint () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename InputType >
LDLT< MatrixType, _UpLo > & compute (const EigenBase< InputType > &a)
 
template<typename InputType >
LDLTcompute (const EigenBase< InputType > &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
bool isNegative (void) const
 
bool isPositive () const
 
 LDLT ()
 Default Constructor. More...
 
template<typename InputType >
 LDLT (const EigenBase< InputType > &matrix)
 Constructor with decomposition. More...
 
template<typename InputType >
 LDLT (EigenBase< InputType > &matrix)
 Constructs a LDLT factorization from a given matrix. More...
 
 LDLT (Index size)
 Default Constructor with memory preallocation. More...
 
Traits::MatrixL matrixL () const
 
const MatrixTypematrixLDLT () const
 
Traits::MatrixU matrixU () const
 
template<typename Derived >
LDLTrankUpdate (const MatrixBase< Derived > &w, const RealScalar &alpha=1)
 
template<typename Derived >
LDLT< MatrixType, _UpLo > & rankUpdate (const MatrixBase< Derived > &w, const typename LDLT< MatrixType, _UpLo >::RealScalar &sigma)
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
void setZero ()
 
template<typename Derived >
bool solveInPlace (MatrixBase< Derived > &bAndX) const
 
const TranspositionTypetranspositionsP () const
 
Diagonal< const MatrixTypevectorD () const
 
- Public Member Functions inherited from Eigen::SolverBase< LDLT< _MatrixType, _UpLo > >
AdjointReturnType adjoint () const
 
EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & derived ()
 
const EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & derived () const
 
const Solve< LDLT< _MatrixType, _UpLo >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
ConstTransposeReturnType transpose () const
 
 ~SolverBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< LDLT< _MatrixType, _UpLo > >
EIGEN_DEVICE_FUNC void addTo (Dest &dst) const
 
EIGEN_DEVICE_FUNC void applyThisOnTheLeft (Dest &dst) const
 
EIGEN_DEVICE_FUNC void applyThisOnTheRight (Dest &dst) const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & const_cast_derived () const
 
const EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & const_derived () const
 
EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & derived ()
 
const EIGEN_DEVICE_FUNC LDLT< _MatrixType, _UpLo > & derived () const
 
EIGEN_DEVICE_FUNC void evalTo (Dest &dst) const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC void subTo (Dest &dst) const
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

ComputationInfo m_info
 
bool m_isInitialized
 
RealScalar m_l1_norm
 
MatrixType m_matrix
 
internal::SignMatrix m_sign
 
TmpMatrixType m_temporary
 
TranspositionType m_transpositions
 

Friends

class SolverBase< LDLT >
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SolverBase< LDLT< _MatrixType, _UpLo > >
void _check_solve_assertion (const Rhs &b) const
 

Detailed Description

template<typename _MatrixType, int _UpLo>
class Eigen::LDLT< _MatrixType, _UpLo >

Robust Cholesky decomposition of a matrix with pivoting.

Template Parameters
_MatrixTypethe type of the matrix of which to compute the LDL^T Cholesky decomposition
_UpLothe triangular part that will be used for the decompositon: Lower (default) or Upper. The other triangular part won't be read.

Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite matrix $ A $ such that $ A = P^TLDL^*P $, where P is a permutation matrix, L is lower triangular with a unit diagonal and D is a diagonal matrix.

The decomposition uses pivoting to ensure stability, so that D will have zeros in the bottom right rank(A) - n submatrix. Avoiding the square root on D also stabilizes the computation.

Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::ldlt(), SelfAdjointView::ldlt(), class LLT

Definition at line 59 of file LDLT.h.

Member Typedef Documentation

◆ Base

template<typename _MatrixType , int _UpLo>
typedef SolverBase<LDLT> Eigen::LDLT< _MatrixType, _UpLo >::Base

Definition at line 64 of file LDLT.h.

◆ MatrixType

template<typename _MatrixType , int _UpLo>
typedef _MatrixType Eigen::LDLT< _MatrixType, _UpLo >::MatrixType

Definition at line 63 of file LDLT.h.

◆ PermutationType

template<typename _MatrixType , int _UpLo>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::LDLT< _MatrixType, _UpLo >::PermutationType

Definition at line 76 of file LDLT.h.

◆ TmpMatrixType

template<typename _MatrixType , int _UpLo>
typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> Eigen::LDLT< _MatrixType, _UpLo >::TmpMatrixType

Definition at line 73 of file LDLT.h.

◆ Traits

template<typename _MatrixType , int _UpLo>
typedef internal::LDLT_Traits<MatrixType,UpLo> Eigen::LDLT< _MatrixType, _UpLo >::Traits

Definition at line 78 of file LDLT.h.

◆ TranspositionType

template<typename _MatrixType , int _UpLo>
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::LDLT< _MatrixType, _UpLo >::TranspositionType

Definition at line 75 of file LDLT.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType , int _UpLo>
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
UpLo 

Definition at line 68 of file LDLT.h.

Constructor & Destructor Documentation

◆ LDLT() [1/4]

template<typename _MatrixType , int _UpLo>
Eigen::LDLT< _MatrixType, _UpLo >::LDLT ( )
inline

Default Constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via LDLT::compute(const MatrixType&).

Definition at line 85 of file LDLT.h.

◆ LDLT() [2/4]

template<typename _MatrixType , int _UpLo>
Eigen::LDLT< _MatrixType, _UpLo >::LDLT ( Index  size)
inlineexplicit

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also
LDLT()

Definition at line 98 of file LDLT.h.

◆ LDLT() [3/4]

template<typename _MatrixType , int _UpLo>
template<typename InputType >
Eigen::LDLT< _MatrixType, _UpLo >::LDLT ( const EigenBase< InputType > &  matrix)
inlineexplicit

Constructor with decomposition.

This calculates the decomposition for the input matrix.

See also
LDLT(Index size)

Definition at line 113 of file LDLT.h.

◆ LDLT() [4/4]

template<typename _MatrixType , int _UpLo>
template<typename InputType >
Eigen::LDLT< _MatrixType, _UpLo >::LDLT ( EigenBase< InputType > &  matrix)
inlineexplicit

Constructs a LDLT factorization from a given matrix.

This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::Ref.

See also
LDLT(const EigenBase&)

Definition at line 130 of file LDLT.h.

Member Function Documentation

◆ _solve_impl()

template<typename _MatrixType , int _UpLo>
template<typename RhsType , typename DstType >
void Eigen::LDLT< _MatrixType, _UpLo >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const

Definition at line 567 of file LDLT.h.

◆ _solve_impl_transposed()

template<typename _MatrixType , int _UpLo>
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::LDLT< _MatrixType, _UpLo >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const

Definition at line 574 of file LDLT.h.

◆ adjoint()

template<typename _MatrixType , int _UpLo>
const LDLT& Eigen::LDLT< _MatrixType, _UpLo >::adjoint ( ) const
inline
Returns
the adjoint of *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint.

This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as:

x = decomposition.adjoint().solve(b)

Definition at line 247 of file LDLT.h.

◆ check_template_parameters()

template<typename _MatrixType , int _UpLo>
static void Eigen::LDLT< _MatrixType, _UpLo >::check_template_parameters ( )
inlinestaticprotected

Definition at line 273 of file LDLT.h.

◆ cols()

template<typename _MatrixType , int _UpLo>
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index Eigen::LDLT< _MatrixType, _UpLo >::cols ( ) const
inline

Definition at line 250 of file LDLT.h.

◆ compute() [1/2]

template<typename _MatrixType , int _UpLo>
template<typename InputType >
LDLT<MatrixType,_UpLo>& Eigen::LDLT< _MatrixType, _UpLo >::compute ( const EigenBase< InputType > &  a)

Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of matrix

Definition at line 499 of file LDLT.h.

◆ compute() [2/2]

template<typename _MatrixType , int _UpLo>
template<typename InputType >
LDLT& Eigen::LDLT< _MatrixType, _UpLo >::compute ( const EigenBase< InputType > &  matrix)

◆ info()

template<typename _MatrixType , int _UpLo>
ComputationInfo Eigen::LDLT< _MatrixType, _UpLo >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the factorization failed because of a zero pivot.

Definition at line 257 of file LDLT.h.

◆ isNegative()

template<typename _MatrixType , int _UpLo>
bool Eigen::LDLT< _MatrixType, _UpLo >::isNegative ( void  ) const
inline
Returns
true if the matrix is negative (semidefinite)

Definition at line 185 of file LDLT.h.

◆ isPositive()

template<typename _MatrixType , int _UpLo>
bool Eigen::LDLT< _MatrixType, _UpLo >::isPositive ( ) const
inline
Returns
true if the matrix is positive (semidefinite)

Definition at line 178 of file LDLT.h.

◆ matrixL()

template<typename _MatrixType , int _UpLo>
Traits::MatrixL Eigen::LDLT< _MatrixType, _UpLo >::matrixL ( ) const
inline
Returns
a view of the lower triangular matrix L

Definition at line 156 of file LDLT.h.

◆ matrixLDLT()

template<typename _MatrixType , int _UpLo>
const MatrixType& Eigen::LDLT< _MatrixType, _UpLo >::matrixLDLT ( ) const
inline
Returns
the internal LDLT decomposition matrix

TODO: document the storage layout

Definition at line 234 of file LDLT.h.

◆ matrixU()

template<typename _MatrixType , int _UpLo>
Traits::MatrixU Eigen::LDLT< _MatrixType, _UpLo >::matrixU ( ) const
inline
Returns
a view of the upper triangular matrix U

Definition at line 149 of file LDLT.h.

◆ rankUpdate() [1/2]

template<typename _MatrixType , int _UpLo>
template<typename Derived >
LDLT& Eigen::LDLT< _MatrixType, _UpLo >::rankUpdate ( const MatrixBase< Derived > &  w,
const RealScalar alpha = 1 
)

◆ rankUpdate() [2/2]

template<typename _MatrixType , int _UpLo>
template<typename Derived >
LDLT<MatrixType,_UpLo>& Eigen::LDLT< _MatrixType, _UpLo >::rankUpdate ( const MatrixBase< Derived > &  w,
const typename LDLT< MatrixType, _UpLo >::RealScalar sigma 
)

Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.

Parameters
wa vector to be incorporated into the decomposition.
sigmaa scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
See also
setZero()

Definition at line 539 of file LDLT.h.

◆ rcond()

template<typename _MatrixType , int _UpLo>
RealScalar Eigen::LDLT< _MatrixType, _UpLo >::rcond ( ) const
inline
Returns
an estimate of the reciprocal condition number of the matrix of which *this is the LDLT decomposition.

Definition at line 221 of file LDLT.h.

◆ reconstructedMatrix()

template<typename MatrixType , int _UpLo>
MatrixType Eigen::LDLT< MatrixType, _UpLo >::reconstructedMatrix
Returns
the matrix represented by the decomposition, i.e., it returns the product: P^T L D L^* P. This function is provided for debug purpose.

Definition at line 643 of file LDLT.h.

◆ rows()

template<typename _MatrixType , int _UpLo>
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index Eigen::LDLT< _MatrixType, _UpLo >::rows ( ) const
inline

Definition at line 249 of file LDLT.h.

◆ setZero()

template<typename _MatrixType , int _UpLo>
void Eigen::LDLT< _MatrixType, _UpLo >::setZero ( )
inline

Clear any existing decomposition

See also
rankUpdate(w,sigma)

Definition at line 143 of file LDLT.h.

◆ solveInPlace()

template<typename MatrixType , int _UpLo>
template<typename Derived >
bool Eigen::LDLT< MatrixType, _UpLo >::solveInPlace ( MatrixBase< Derived > &  bAndX) const

Definition at line 629 of file LDLT.h.

◆ transpositionsP()

template<typename _MatrixType , int _UpLo>
const TranspositionType& Eigen::LDLT< _MatrixType, _UpLo >::transpositionsP ( ) const
inline
Returns
the permutation matrix P as a transposition sequence.

Definition at line 164 of file LDLT.h.

◆ vectorD()

template<typename _MatrixType , int _UpLo>
Diagonal<const MatrixType> Eigen::LDLT< _MatrixType, _UpLo >::vectorD ( ) const
inline
Returns
the coefficients of the diagonal matrix D

Definition at line 171 of file LDLT.h.

Friends And Related Function Documentation

◆ SolverBase< LDLT >

template<typename _MatrixType , int _UpLo>
friend class SolverBase< LDLT >
friend

Definition at line 65 of file LDLT.h.

Member Data Documentation

◆ m_info

template<typename _MatrixType , int _UpLo>
ComputationInfo Eigen::LDLT< _MatrixType, _UpLo >::m_info
protected

Definition at line 290 of file LDLT.h.

◆ m_isInitialized

template<typename _MatrixType , int _UpLo>
bool Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized
protected

Definition at line 289 of file LDLT.h.

◆ m_l1_norm

template<typename _MatrixType , int _UpLo>
RealScalar Eigen::LDLT< _MatrixType, _UpLo >::m_l1_norm
protected

Definition at line 285 of file LDLT.h.

◆ m_matrix

template<typename _MatrixType , int _UpLo>
MatrixType Eigen::LDLT< _MatrixType, _UpLo >::m_matrix
protected

Definition at line 284 of file LDLT.h.

◆ m_sign

template<typename _MatrixType , int _UpLo>
internal::SignMatrix Eigen::LDLT< _MatrixType, _UpLo >::m_sign
protected

Definition at line 288 of file LDLT.h.

◆ m_temporary

template<typename _MatrixType , int _UpLo>
TmpMatrixType Eigen::LDLT< _MatrixType, _UpLo >::m_temporary
protected

Definition at line 287 of file LDLT.h.

◆ m_transpositions

template<typename _MatrixType , int _UpLo>
TranspositionType Eigen::LDLT< _MatrixType, _UpLo >::m_transpositions
protected

Definition at line 286 of file LDLT.h.


The documentation for this class was generated from the following file:
b
Scalar * b
Definition: benchVecAdd.cpp:17
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12


gtsam
Author(s):
autogenerated on Thu Dec 19 2024 04:09:46