Public Types | Public Member Functions | Protected Attributes
Eigen::LDLT< _MatrixType, _UpLo > Class Template Reference

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

#include <LDLT.h>

List of all members.

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options & ~RowMajorBit, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, UpLo = _UpLo
}
typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef PermutationMatrix
< RowsAtCompileTime,
MaxRowsAtCompileTime
PermutationType
typedef NumTraits< typename
MatrixType::Scalar >::Real 
RealScalar
typedef MatrixType::Scalar Scalar
typedef Matrix< Scalar,
RowsAtCompileTime, 1, Options,
MaxRowsAtCompileTime, 1 > 
TmpMatrixType
typedef internal::LDLT_Traits
< MatrixType, UpLo
Traits
typedef Transpositions
< RowsAtCompileTime,
MaxRowsAtCompileTime
TranspositionType

Public Member Functions

Index cols () const
LDLTcompute (const MatrixType &matrix)
ComputationInfo info () const
 Reports whether previous computation was successful.
bool isNegative (void) const
bool isPositive () const
 LDLT ()
 Default Constructor.
 LDLT (Index size)
 Default Constructor with memory preallocation.
 LDLT (const MatrixType &matrix)
 Constructor with decomposition.
Traits::MatrixL matrixL () const
const MatrixTypematrixLDLT () const
Traits::MatrixU matrixU () const
template<typename Derived >
LDLTrankUpdate (const MatrixBase< Derived > &w, RealScalar alpha=1)
template<typename Derived >
LDLT< MatrixType, _UpLo > & rankUpdate (const MatrixBase< Derived > &w, typename NumTraits< typename MatrixType::Scalar >::Real sigma)
MatrixType reconstructedMatrix () const
Index rows () const
void setZero ()
template<typename Rhs >
const internal::solve_retval
< LDLT, Rhs > 
solve (const MatrixBase< Rhs > &b) const
template<typename Derived >
bool solveInPlace (MatrixBase< Derived > &bAndX) const
const TranspositionTypetranspositionsP () const
Diagonal< const MatrixTypevectorD () const

Protected Attributes

bool m_isInitialized
MatrixType m_matrix
int 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.

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.

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

Definition at line 45 of file LDLT.h.


Member Typedef Documentation

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

Definition at line 59 of file LDLT.h.

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

Definition at line 48 of file LDLT.h.

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

Definition at line 63 of file LDLT.h.

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

Definition at line 58 of file LDLT.h.

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

Definition at line 57 of file LDLT.h.

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

Definition at line 60 of file LDLT.h.

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

Definition at line 65 of file LDLT.h.

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

Definition at line 62 of file LDLT.h.


Member Enumeration Documentation

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

Definition at line 49 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 72 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
Eigen::LDLT< _MatrixType, _UpLo >::LDLT ( Index  size) [inline]

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

template<typename _MatrixType, int _UpLo>
Eigen::LDLT< _MatrixType, _UpLo >::LDLT ( const MatrixType matrix) [inline]

Constructor with decomposition.

This calculates the decomposition for the input matrix.

See also:
LDLT(Index size)

Definition at line 92 of file LDLT.h.


Member Function Documentation

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

Definition at line 214 of file LDLT.h.

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

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

Definition at line 428 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 matrix.appears to be negative.

Definition at line 221 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 153 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 139 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 117 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 205 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 110 of file LDLT.h.

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

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

Definition at line 213 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 104 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
template<typename Rhs >
const internal::solve_retval<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()

Definition at line 176 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 534 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 125 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 132 of file LDLT.h.


Member Data Documentation

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

Definition at line 239 of file LDLT.h.

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

Definition at line 235 of file LDLT.h.

template<typename _MatrixType, int _UpLo>
int Eigen::LDLT< _MatrixType, _UpLo >::m_sign [protected]

Definition at line 238 of file LDLT.h.

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

Definition at line 237 of file LDLT.h.

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

Definition at line 236 of file LDLT.h.


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


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:12:54