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

Standard Cholesky decomposition (LL^T) of a matrix and associated features. More...

#include <LLT.h>

Public Types

enum  { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
enum  { PacketSize = internal::packet_traits<Scalar>::size, AlignmentMask = int(PacketSize)-1, UpLo = _UpLo }
 
typedef Eigen::Index Index
 
typedef _MatrixType MatrixType
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef internal::LLT_Traits< MatrixType, UpLoTraits
 

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 LLTadjoint () const
 
Index cols () const
 
template<typename InputType >
LLTcompute (const EigenBase< InputType > &matrix)
 
template<typename InputType >
LLT< MatrixType, _UpLo > & compute (const EigenBase< InputType > &a)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
 LLT ()
 Default Constructor. More...
 
 LLT (Index size)
 Default Constructor with memory preallocation. More...
 
template<typename InputType >
 LLT (const EigenBase< InputType > &matrix)
 
template<typename InputType >
 LLT (EigenBase< InputType > &matrix)
 Constructs a LDLT factorization from a given matrix. More...
 
Traits::MatrixL matrixL () const
 
const MatrixTypematrixLLT () const
 
Traits::MatrixU matrixU () const
 
template<typename VectorType >
LLT rankUpdate (const VectorType &vec, const RealScalar &sigma=1)
 
template<typename VectorType >
LLT< _MatrixType, _UpLo > rankUpdate (const VectorType &v, const RealScalar &sigma)
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
Index rows () const
 
template<typename Rhs >
const Solve< LLT, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Derived >
void solveInPlace (MatrixBase< Derived > &bAndX) 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
 

Detailed Description

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

Standard Cholesky decomposition (LL^T) of a matrix and associated features.

Template Parameters
_MatrixTypethe type of the matrix of which we are computing the LL^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.

This class performs a LL^T Cholesky decomposition of a symmetric, positive definite matrix A such that A = LL^* = U^*U, where L is lower triangular.

While the Cholesky decomposition is particularly useful to solve selfadjoint problems like D^*D x = b, for that purpose, we recommend the Cholesky decomposition without square root which is more stable and even faster. Nevertheless, this standard Cholesky decomposition remains useful in many other situations like generalised eigen problems with hermitian matrices.

Remember that Cholesky decompositions are not rank-revealing. This LLT decomposition is only stable on positive definite matrices, use LDLT instead for the semidefinite case. Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.

Example:

Output:

This class supports the inplace decomposition mechanism.

See also
MatrixBase::llt(), SelfAdjointView::llt(), class LDLT

Definition at line 52 of file LLT.h.

Member Typedef Documentation

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

Definition at line 63 of file LLT.h.

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

Definition at line 55 of file LLT.h.

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

Definition at line 62 of file LLT.h.

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

Definition at line 61 of file LLT.h.

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

Definition at line 64 of file LLT.h.

template<typename _MatrixType, int _UpLo>
typedef internal::LLT_Traits<MatrixType,UpLo> Eigen::LLT< _MatrixType, _UpLo >::Traits

Definition at line 72 of file LLT.h.

Member Enumeration Documentation

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

Definition at line 56 of file LLT.h.

template<typename _MatrixType, int _UpLo>
anonymous enum
Enumerator
PacketSize 
AlignmentMask 
UpLo 

Definition at line 66 of file LLT.h.

Constructor & Destructor Documentation

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

Default Constructor.

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

Definition at line 80 of file LLT.h.

template<typename _MatrixType, int _UpLo>
Eigen::LLT< _MatrixType, _UpLo >::LLT ( 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
LLT()

Definition at line 88 of file LLT.h.

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

Definition at line 92 of file LLT.h.

template<typename _MatrixType, int _UpLo>
template<typename InputType >
Eigen::LLT< _MatrixType, _UpLo >::LLT ( 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
LLT(const EigenBase&)

Definition at line 107 of file LLT.h.

Member Function Documentation

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

Definition at line 473 of file LLT.h.

template<typename _MatrixType, int _UpLo>
const LLT& Eigen::LLT< _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 193 of file LLT.h.

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

Definition at line 209 of file LLT.h.

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

Definition at line 196 of file LLT.h.

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

Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of matrix

Returns
a reference to *this

Example:

Output:

 

Definition at line 421 of file LLT.h.

template<typename _MatrixType, int _UpLo>
ComputationInfo Eigen::LLT< _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 182 of file LLT.h.

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

Definition at line 122 of file LLT.h.

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

TODO: document the storage layout

Definition at line 168 of file LLT.h.

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

Definition at line 115 of file LLT.h.

template<typename _MatrixType, int _UpLo>
template<typename VectorType >
LLT Eigen::LLT< _MatrixType, _UpLo >::rankUpdate ( const VectorType vec,
const RealScalar sigma = 1 
)
template<typename _MatrixType, int _UpLo>
template<typename VectorType >
LLT<_MatrixType,_UpLo> Eigen::LLT< _MatrixType, _UpLo >::rankUpdate ( const VectorType v,
const RealScalar sigma 
)

Performs a rank one update (or dowdate) of the current decomposition. If A = LL^* before the rank one update, then after it we have LL^* = A + sigma * v v^* where v must be a vector of same dimension.

Definition at line 457 of file LLT.h.

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

Definition at line 157 of file LLT.h.

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

Definition at line 504 of file LLT.h.

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

Definition at line 195 of file LLT.h.

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

Since this LLT class assumes anyway that the matrix A is invertible, the solution theoretically exists and is unique regardless of b.

Example:

Output:

See also
solveInPlace(), MatrixBase::llt(), SelfAdjointView::llt()

Definition at line 140 of file LLT.h.

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

Definition at line 492 of file LLT.h.

Member Data Documentation

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

Definition at line 221 of file LLT.h.

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

Definition at line 220 of file LLT.h.

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

Definition at line 219 of file LLT.h.

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

Definition at line 218 of file LLT.h.


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


hebiros
Author(s): Xavier Artache , Matthew Tesch
autogenerated on Thu Sep 3 2020 04:10:07