ROOT 6.18/05 Reference Guide |
class to compute the Cholesky decomposition of a matrix
class to compute the Cholesky decomposition of a symmetric positive definite matrix
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:
Definition at line 76 of file CholeskyDecomp.h.
Public Member Functions | |
template<class M > | |
CholeskyDecomp (const M &m) | |
perform a Cholesky decomposition More... | |
template<typename G > | |
CholeskyDecomp (G *m) | |
perform a Cholesky decomposition 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 | |
F | fL [N *(N+1)/2] |
lower triangular matrix L More... | |
bool | fOk |
flag indicating a successful decomposition More... | |
#include <Math/CholeskyDecomp.h>
|
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 94 of file CholeskyDecomp.h.
|
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 112 of file CholeskyDecomp.h.
|
inline |
obtain the decomposed matrix L
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 208 of file CholeskyDecomp.h.
|
inline |
obtain the decomposed matrix L
This is the method to use with a plain array.
Definition at line 183 of file CholeskyDecomp.h.
|
inline |
obtain the inverse of the decomposed matrix L
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 258 of file CholeskyDecomp.h.
|
inline |
obtain the inverse of the decomposed matrix L
This is the method to use with a plain array.
Definition at line 227 of file CholeskyDecomp.h.
|
inline |
place the inverse into m
This is the method to use with a plain array.
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 166 of file CholeskyDecomp.h.
|
inline |
place the inverse into m
This is the method to use with an SMatrix.
Definition at line 149 of file CholeskyDecomp.h.
|
inline |
returns true if decomposition was successful
Definition at line 123 of file CholeskyDecomp.h.
|
inline |
returns true if decomposition was successful
Definition at line 126 of file CholeskyDecomp.h.
|
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]
Definition at line 136 of file CholeskyDecomp.h.
lower triangular matrix L
lower triangular matrix L, packed storage, with diagonal elements pre-inverted
Definition at line 82 of file CholeskyDecomp.h.
|
private |
flag indicating a successful decomposition
Definition at line 84 of file CholeskyDecomp.h.