TooN
Public Member Functions | List of all members
TooN::Cholesky< Size, Precision > Class Template Reference

Decomposes a positive-semidefinite symmetric matrix A (such as a covariance) into L*D*L^T, where L is lower-triangular and D is diagonal. More...

#include <Cholesky.h>

Public Member Functions

template<class P2 , class B2 >
 Cholesky (const Matrix< Size, Size, P2, B2 > &m)
 Construct the Cholesky decomposition of a matrix. More...
 
 Cholesky (int size)
 Constructor for Size=Dynamic.
 
template<class P2 , class B2 >
void compute (const Matrix< Size, Size, P2, B2 > &m)
 Compute the LDL^T decomposition of another matrix. More...
 
template<int Size2, class P2 , class B2 >
Vector< Size, Precision > backsub (const Vector< Size2, P2, B2 > &v) const
 Compute x = A^-1*v Run time is O(N^2)
 
template<int Size2, int C2, class P2 , class B2 >
Matrix< Size, C2, Precision > backsub (const Matrix< Size2, C2, P2, B2 > &m) const
 overload
 
Matrix< Size, Size, Precision > get_inverse ()
 Compute A^-1 and store in M Run time is O(N^3)
 
Precision determinant ()
 Compute the determinant.
 
template<int Size2, typename P2 , typename B2 >
Precision mahalanobis (const Vector< Size2, P2, B2 > &v) const
 
Matrix< Size, Size, Precision > get_unscaled_L () const
 
Matrix< Size, Size, Precision > get_D () const
 
Matrix< Size, Size, Precision > get_L () const
 
int rank () const
 

Detailed Description

template<int Size = Dynamic, class Precision = DefaultPrecision>
class TooN::Cholesky< Size, Precision >

Decomposes a positive-semidefinite symmetric matrix A (such as a covariance) into L*D*L^T, where L is lower-triangular and D is diagonal.

Also can compute the classic A = L*L^T, with L lower triangular. The LDL^T form is faster to compute than the classical Cholesky decomposition. Use get_unscaled_L() and get_D() to access the individual matrices of L*D*L^T decomposition. Use get_L() to access the lower triangular matrix of the classic Cholesky decomposition L*L^T. The decomposition can be used to compute A^-1*x, A^-1*M, M*A^-1*M^T, and A^-1 itself, though the latter rarely needs to be explicitly represented. Also efficiently computes det(A) and rank(A). It can be used as follows:

// Declare some matrices.
Matrix<3> A = ...; // we'll pretend it is pos-def
Matrix<2,3> M;
Matrix<2> B;
Vector<3> y = make_Vector(2,3,4);
// create the Cholesky decomposition of A
Cholesky<3> chol(A);
// compute x = A^-1 * y
x = cholA.backsub(y);
//compute A^-1
Matrix<3> Ainv = cholA.get_inverse();

Cholesky decomposition of a symmetric matrix. Only the lower half of the matrix is considered This uses the non-sqrt version of the decomposition giving symmetric M = L*D*L.T() where the diagonal of L contains ones

Parameters
Sizethe size of the matrix
Precisionthe precision of the entries in the matrix and its decomposition

Constructor & Destructor Documentation

◆ Cholesky()

template<int Size = Dynamic, class Precision = DefaultPrecision>
template<class P2 , class B2 >
TooN::Cholesky< Size, Precision >::Cholesky ( const Matrix< Size, Size, P2, B2 > &  m)
inline

Construct the Cholesky decomposition of a matrix.

This initialises the class, and performs the decomposition immediately. Run time is O(N^3)

Member Function Documentation

◆ compute()

template<int Size = Dynamic, class Precision = DefaultPrecision>
template<class P2 , class B2 >
void TooN::Cholesky< Size, Precision >::compute ( const Matrix< Size, Size, P2, B2 > &  m)
inline

Compute the LDL^T decomposition of another matrix.

Run time is O(N^3)


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