Logo ROOT   6.16/01
Reference Guide
List of all members | Public Member Functions | Private Attributes | List of all members
ROOT::Math::CholeskyDecompGenDim< F > Class Template Reference

template<class F>
class ROOT::Math::CholeskyDecompGenDim< F >

class to compute the Cholesky decomposition of a matrix

class to compute the Cholesky decomposition of a symmetric positive definite matrix when the dimensionality of the problem is not known at compile time

provides routines to check if the decomposition succeeded (i.e. if matrix is positive definite and non-singular), to solve a linear system for the given matrix and to obtain its inverse

the actual functionality is implemented in templated helper classes which have specializations for dimensions N = 1 to 6 to achieve a gain in speed for common matrix sizes

usage example:

// let m be a symmetric positive definite SMatrix (use type float
// for internal computations, matrix size is 4x4)
CholeskyDecompGenDim<float> decomp(4, m);
// check if the decomposition succeeded
if (!decomp) {
std::cerr << "decomposition failed!" << std::endl;
} else {
// let rhs be a vector; we seek a vector x such that m * x = rhs
decomp.Solve(rhs);
// rhs now contains the solution we are looking for
// obtain the inverse of m, put it into m itself
decomp.Invert(m);
}
auto * m
Definition: textangle.C:8

Definition at line 310 of file CholeskyDecomp.h.

Public Member Functions

template<class M >
 CholeskyDecompGenDim (unsigned N, const M &m)
 perform a Cholesky decomposition More...
 
template<typename G >
 CholeskyDecompGenDim (unsigned N, G *m)
 perform a Cholesky decomposition More...
 
 ~CholeskyDecompGenDim ()
 destructor More...
 
template<typename G >
bool getL (G *m) const
 obtain the decomposed matrix L More...
 
template<class M >
bool getL (M &m) const
 obtain the decomposed matrix L More...
 
template<typename G >
bool getLi (G *m) const
 obtain the inverse of the decomposed matrix L More...
 
template<class M >
bool getLi (M &m) const
 obtain the inverse of the decomposed matrix L More...
 
template<typename G >
bool Invert (G *m) const
 place the inverse into m More...
 
template<class M >
bool Invert (M &m) const
 place the inverse into m More...
 
bool ok () const
 returns true if decomposition was successful More...
 
 operator bool () const
 returns true if decomposition was successful More...
 
template<class V >
bool Solve (V &rhs) const
 solves a linear system for the given right hand side More...
 

Private Attributes

FfL
 lower triangular matrix L More...
 
unsigned fN
 dimensionality dimensionality of the problem More...
 
bool fOk
 flag indicating a successful decomposition More...
 

#include <Math/CholeskyDecomp.h>

Constructor & Destructor Documentation

◆ CholeskyDecompGenDim() [1/2]

template<class F >
template<class M >
ROOT::Math::CholeskyDecompGenDim< F >::CholeskyDecompGenDim ( unsigned  N,
const M &  m 
)
inline

perform a Cholesky decomposition

perfrom a Cholesky decomposition of a symmetric positive definite matrix m

this is the constructor to uses with an SMatrix (and objects that behave like an SMatrix in terms of using operator()(int i, int j) for access to elements)

Definition at line 331 of file CholeskyDecomp.h.

◆ CholeskyDecompGenDim() [2/2]

template<class F >
template<typename G >
ROOT::Math::CholeskyDecompGenDim< F >::CholeskyDecompGenDim ( unsigned  N,
G m 
)
inline

perform a Cholesky decomposition

perfrom a Cholesky decomposition of a symmetric positive definite matrix m

this is the constructor to use in special applications where plain arrays are used

NOTE: the matrix is given in packed representation, matrix element m(i,j) (j <= i) is supposed to be in array element (i * (i + 1)) / 2 + j

Definition at line 349 of file CholeskyDecomp.h.

◆ ~CholeskyDecompGenDim()

template<class F >
ROOT::Math::CholeskyDecompGenDim< F >::~CholeskyDecompGenDim ( )
inline

destructor

Definition at line 359 of file CholeskyDecomp.h.

Member Function Documentation

◆ getL() [1/2]

template<class F >
template<typename G >
bool ROOT::Math::CholeskyDecompGenDim< F >::getL ( G m) const
inline

obtain the decomposed matrix L

Returns
if the decomposition was successful

NOTE: the matrix is given in packed representation, matrix element m(i,j) (j <= i) is supposed to be in array element (i * (i + 1)) / 2 + j

Definition at line 448 of file CholeskyDecomp.h.

◆ getL() [2/2]

template<class F >
template<class M >
bool ROOT::Math::CholeskyDecompGenDim< F >::getL ( M &  m) const
inline

obtain the decomposed matrix L

This is the method to use with a plain array.

Returns
if the decomposition was successful

Definition at line 423 of file CholeskyDecomp.h.

◆ getLi() [1/2]

template<class F >
template<typename G >
bool ROOT::Math::CholeskyDecompGenDim< F >::getLi ( G m) const
inline

obtain the inverse of the decomposed matrix L

Returns
if the decomposition was successful

NOTE: the matrix is given in packed representation, matrix element m(j,i) (j <= i) is supposed to be in array element (i * (i + 1)) / 2 + j

Definition at line 498 of file CholeskyDecomp.h.

◆ getLi() [2/2]

template<class F >
template<class M >
bool ROOT::Math::CholeskyDecompGenDim< F >::getLi ( M &  m) const
inline

obtain the inverse of the decomposed matrix L

This is the method to use with a plain array.

Returns
if the decomposition was successful

Definition at line 467 of file CholeskyDecomp.h.

◆ Invert() [1/2]

template<class F >
template<typename G >
bool ROOT::Math::CholeskyDecompGenDim< F >::Invert ( G m) const
inline

place the inverse into m

This is the method to use with a plain array.

Returns
if the decomposition was successful

NOTE: the matrix is given in packed representation, matrix element m(i,j) (j <= i) is supposed to be in array element (i * (i + 1)) / 2 + j

Definition at line 406 of file CholeskyDecomp.h.

◆ Invert() [2/2]

template<class F >
template<class M >
bool ROOT::Math::CholeskyDecompGenDim< F >::Invert ( M &  m) const
inline

place the inverse into m

This is the method to use with an SMatrix.

Returns
if the decomposition was successful

Definition at line 389 of file CholeskyDecomp.h.

◆ ok()

template<class F >
bool ROOT::Math::CholeskyDecompGenDim< F >::ok ( ) const
inline

returns true if decomposition was successful

Returns
true if decomposition was successful

Definition at line 363 of file CholeskyDecomp.h.

◆ operator bool()

template<class F >
ROOT::Math::CholeskyDecompGenDim< F >::operator bool ( ) const
inline

returns true if decomposition was successful

Returns
true if decomposition was successful

Definition at line 366 of file CholeskyDecomp.h.

◆ Solve()

template<class F >
template<class V >
bool ROOT::Math::CholeskyDecompGenDim< F >::Solve ( V &  rhs) const
inline

solves a linear system for the given right hand side

Note that you can use both SVector classes and plain arrays for rhs. (Make sure that the sizes match!). It will work with any vector implementing the operator [i]

Returns
if the decomposition was successful

Definition at line 376 of file CholeskyDecomp.h.

Member Data Documentation

◆ fL

template<class F >
F* ROOT::Math::CholeskyDecompGenDim< F >::fL
private

lower triangular matrix L

lower triangular matrix L, packed storage, with diagonal elements pre-inverted

Definition at line 319 of file CholeskyDecomp.h.

◆ fN

template<class F >
unsigned ROOT::Math::CholeskyDecompGenDim< F >::fN
private

dimensionality dimensionality of the problem

Definition at line 315 of file CholeskyDecomp.h.

◆ fOk

template<class F >
bool ROOT::Math::CholeskyDecompGenDim< F >::fOk
private

flag indicating a successful decomposition

Definition at line 321 of file CholeskyDecomp.h.


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