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

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

#include <LDLT.h>

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
  UpLo = _UpLo
}
 
typedef Eigen::Index Index
 
typedef _MatrixType MatrixType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTimePermutationType
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
 
typedef internal::LDLT_Traits< MatrixType, UpLoTraits
 
typedef Transpositions< RowsAtCompileTime, MaxRowsAtCompileTimeTranspositionType
 

Public Member Functions

template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
const LDLTadjoint () const
 
Index cols () const
 
template<typename InputType >
LDLTcompute (const EigenBase< InputType > &matrix)
 
template<typename InputType >
LDLT< MatrixType, _UpLo > & compute (const EigenBase< InputType > &a)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
bool isNegative (void) const
 
bool isPositive () const
 
 LDLT ()
 Default Constructor. More...
 
 LDLT (Index size)
 Default Constructor with memory preallocation. 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...
 
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
 
Index rows () const
 
void setZero ()
 
template<typename Rhs >
const Solve< LDLT, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Derived >
bool solveInPlace (MatrixBase< Derived > &bAndX) const
 
const TranspositionTypetranspositionsP () const
 
Diagonal< const MatrixTypevectorD () 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
 

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 L 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 50 of file LDLT.h.

Member Typedef Documentation

template<typename _MatrixType, int _UpLo>
typedef Eigen::Index Eigen::LDLT< _MatrixType, _UpLo >::Index
Deprecated:
since Eigen 3.3

Definition at line 63 of file LDLT.h.

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

Definition at line 53 of file LDLT.h.

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

Definition at line 68 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::LDLT< _MatrixType, _UpLo >::RealScalar

Definition at line 62 of file LDLT.h.

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

Definition at line 61 of file LDLT.h.

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

Definition at line 64 of file LDLT.h.

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

Definition at line 65 of file LDLT.h.

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

Definition at line 70 of file LDLT.h.

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

Definition at line 67 of file LDLT.h.

Member Enumeration Documentation

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

Definition at line 54 of file LDLT.h.

Constructor & Destructor Documentation

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 77 of file LDLT.h.

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 90 of file LDLT.h.

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 105 of file LDLT.h.

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 122 of file LDLT.h.

Member Function Documentation

template<typename _MatrixType, int _UpLo>
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void Eigen::LDLT< _MatrixType, _UpLo >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
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 561 of file LDLT.h.

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 243 of file LDLT.h.

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

Definition at line 267 of file LDLT.h.

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

Definition at line 246 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
template<typename InputType >
LDLT& Eigen::LDLT< _MatrixType, _UpLo >::compute ( const EigenBase< InputType > &  matrix)
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 493 of file LDLT.h.

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

Reports whether previous computation was successful.

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

Definition at line 253 of file LDLT.h.

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 177 of file LDLT.h.

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

Definition at line 170 of file LDLT.h.

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 148 of file LDLT.h.

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 230 of file LDLT.h.

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 141 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
template<typename Derived >
LDLT& Eigen::LDLT< _MatrixType, _UpLo >::rankUpdate ( const MatrixBase< Derived > &  w,
const RealScalar alpha = 1 
)
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 533 of file LDLT.h.

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 217 of file LDLT.h.

template<typename MatrixType , int _UpLo>
MatrixType Eigen::LDLT< MatrixType, _UpLo >::reconstructedMatrix ( ) const
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 628 of file LDLT.h.

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

Definition at line 245 of file LDLT.h.

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

Clear any existing decomposition

See also
rankUpdate(w,sigma)

Definition at line 135 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
template<typename Rhs >
const Solve<LDLT, Rhs> Eigen::LDLT< _MatrixType, _UpLo >::solve ( const MatrixBase< Rhs > &  b) const
inline
Returns
a solution x of $ A x = b $ using the current decomposition of A.

This function also supports in-place solves using the syntax x = decompositionObject.solve(x) .

More precisely, this method solves $ A x = b $ using the decomposition $ A = P^T L D L^* P $ by solving the systems $ P^T y_1 = b $, $ L y_2 = y_1 $, $ D y_3 = y_2 $, $ L^* y_4 = y_3 $ and $ P x = y_4 $ in succession. If the matrix $ A $ is singular, then $ D $ will also be singular (all the other matrices are invertible). In that case, the least-square solution of $ D y_3 = y_2 $ is computed. This does not mean that this function computes the least-square solution of $ A x = b $ is $ A $ is singular.

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

Definition at line 200 of file LDLT.h.

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

Definition at line 614 of file LDLT.h.

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 156 of file LDLT.h.

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 163 of file LDLT.h.

Member Data Documentation

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

Definition at line 284 of file LDLT.h.

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

Definition at line 283 of file LDLT.h.

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

Definition at line 279 of file LDLT.h.

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

Definition at line 278 of file LDLT.h.

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

Definition at line 282 of file LDLT.h.

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

Definition at line 281 of file LDLT.h.

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

Definition at line 280 of file LDLT.h.


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


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:52:50