Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions | Protected Attributes
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, Options = MatrixType::Options, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
enum  { PacketSize = internal::packet_traits<Scalar>::size, AlignmentMask = int(PacketSize)-1, UpLo = _UpLo }
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef MatrixType::Index Index
 
typedef internal::LLT_Traits< MatrixType, UpLo > Traits
 

Public Member Functions

 LLT ()
 Default Constructor. More...
 
 LLT (Index size)
 Default Constructor with memory preallocation. More...
 
 LLT (const MatrixType &matrix)
 
Traits::MatrixU matrixU () const
 
Traits::MatrixL matrixL () const
 
template<typename Rhs >
const internal::solve_retval< LLT, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Derived >
void solveInPlace (MatrixBase< Derived > &bAndX) const
 
LLTcompute (const MatrixType &matrix)
 
const MatrixType & matrixLLT () const
 
MatrixType reconstructedMatrix () const
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
Index rows () const
 
Index cols () 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)
 

Protected Attributes

MatrixType m_matrix
 
bool m_isInitialized
 
ComputationInfo m_info
 

Detailed Description

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

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

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:

See also
MatrixBase::llt(), class LDLT

Definition at line 50 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 78 of file LLT.h.

78 : m_matrix(), m_isInitialized(false) {}
template<typename _MatrixType, int _UpLo>
Eigen::LLT< _MatrixType, _UpLo >::LLT ( 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
LLT()

Definition at line 86 of file LLT.h.

86  : m_matrix(size, size),
87  m_isInitialized(false) {}

Member Function Documentation

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

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

Returns
a reference to *this

Example:

Output:

 

Definition at line 385 of file LLT.h.

386 {
387  eigen_assert(a.rows()==a.cols());
388  const Index size = a.rows();
389  m_matrix.resize(size, size);
390  m_matrix = a;
391 
392  m_isInitialized = true;
393  bool ok = Traits::inplace_decomposition(m_matrix);
394  m_info = ok ? Success : NumericalIssue;
395 
396  return *this;
397 }
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 164 of file LLT.h.

165  {
166  eigen_assert(m_isInitialized && "LLT is not initialized.");
167  return m_info;
168  }
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 104 of file LLT.h.

105  {
106  eigen_assert(m_isInitialized && "LLT is not initialized.");
107  return Traits::getL(m_matrix);
108  }
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 150 of file LLT.h.

151  {
152  eigen_assert(m_isInitialized && "LLT is not initialized.");
153  return m_matrix;
154  }
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 97 of file LLT.h.

98  {
99  eigen_assert(m_isInitialized && "LLT is not initialized.");
100  return Traits::getU(m_matrix);
101  }
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 406 of file LLT.h.

407 {
408  EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
409  eigen_assert(v.size()==m_matrix.cols());
410  eigen_assert(m_isInitialized);
411  if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
412  m_info = NumericalIssue;
413  else
414  m_info = Success;
415 
416  return *this;
417 }
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 462 of file LLT.h.

463 {
464  eigen_assert(m_isInitialized && "LLT is not initialized.");
465  return matrixL() * matrixL().adjoint().toDenseMatrix();
466 }
Traits::MatrixL matrixL() const
Definition: LLT.h:104
template<typename _MatrixType, int _UpLo>
template<typename Rhs >
const internal::solve_retval<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()

Definition at line 122 of file LLT.h.

123  {
124  eigen_assert(m_isInitialized && "LLT is not initialized.");
125  eigen_assert(m_matrix.rows()==b.rows()
126  && "LLT::solve(): invalid number of rows of the right hand side matrix b");
127  return internal::solve_retval<LLT, Rhs>(*this, b.derived());
128  }

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