latest version v1.9 - last update 10 Apr 2010 |
This functor calculates the Cholesky decomposition of a symmetric, positive-definite matrix A. More...
#include <ltiCholeskyDecomposition.h>
Classes | |
class | parameters |
the parameters for the class choleskyDecomposition More... | |
Public Member Functions | |
choleskyDecomposition () | |
choleskyDecomposition (const parameters &par) | |
choleskyDecomposition (const choleskyDecomposition &other) | |
virtual | ~choleskyDecomposition () |
virtual const char * | getTypeName () const |
bool | apply (matrix< T > &srcdest) const |
bool | apply (const matrix< T > &src, matrix< T > &l) const |
bool | apply (matrix< T > &srcdest, const typename parameters::eTriangularType &tType) const |
bool | apply (const matrix< T > &src, matrix< T > &l, const typename parameters::eTriangularType &tType) const |
choleskyDecomposition & | copy (const choleskyDecomposition &other) |
choleskyDecomposition & | operator= (const choleskyDecomposition &other) |
virtual functor * | clone () const |
const parameters & | getParameters () const |
This functor calculates the Cholesky decomposition of a symmetric, positive-definite matrix A.
The Cholesky decomposition is defined as A = L * L' (with L' the transpose of L) or A = U' * U. The cholesky factor L or U are lower or upper triangular matrices, respectively, where L=U'. The nature of the Cholesky factor can be determined by setting parameters::triangularType or using the appropriate apply method. The default is Upper due to shorter calculation time.
The Cholesky decomposition can be used to solve a linear equation system A*x=b by first solving L*y=b and then L'*x=y using the forwardSubstitution and backSubstitution functors, respectively.
lti::choleskyDecomposition< T >::choleskyDecomposition | ( | ) |
default constructor
lti::choleskyDecomposition< T >::choleskyDecomposition | ( | const parameters & | par | ) |
Construct a functor using the given parameters.
lti::choleskyDecomposition< T >::choleskyDecomposition | ( | const choleskyDecomposition< T > & | other | ) |
copy constructor
other | the object to be copied |
virtual lti::choleskyDecomposition< T >::~choleskyDecomposition | ( | ) | [virtual] |
destructor
bool lti::choleskyDecomposition< T >::apply | ( | const matrix< T > & | src, | |
matrix< T > & | l, | |||
const typename parameters::eTriangularType & | tType | |||
) | const |
The matrix given must be symmetric and positive-definite.
The apply method returns the upper or lower triangular matrix version of the Cholosky decomposition, depending on the value of tType.
src | matrix<T> with the source data. | |
l | matrix<T> where the result will be left. | |
tType | type of triangular matrix |
bool lti::choleskyDecomposition< T >::apply | ( | matrix< T > & | srcdest, | |
const typename parameters::eTriangularType & | tType | |||
) | const |
The matrix given must be symmetric and positive-definite.
The apply method returns the upper or lower triangular matrix version of the Cholosky decomposition, depending on the value of tType.
srcdest | matrix<T> with the source data. The result will be left here too. | |
tType | type of triangular matrix |
bool lti::choleskyDecomposition< T >::apply | ( | const matrix< T > & | src, | |
matrix< T > & | l | |||
) | const |
The matrix given must be symmetric and positive-definite.
The apply method returns the upper or lower triangular matrix version of the Cholosky decomposition, depending on the value of parameters::triangularType. The default Upper is faster.
src | matrix<T> with the source data. | |
l | matrix<T> where the result will be left. |
bool lti::choleskyDecomposition< T >::apply | ( | matrix< T > & | srcdest | ) | const |
The matrix given must be symmetric and positive-definite.
The apply method returns the upper or lower triangular matrix version of the Cholosky decomposition, depending on the value of parameters::triangularType. The default Upper is faster.
srcdest | matrix<T> with the source data. The result will be left here too. |
virtual functor* lti::choleskyDecomposition< T >::clone | ( | ) | const [virtual] |
returns a pointer to a clone of this functor.
Implements lti::functor.
choleskyDecomposition& lti::choleskyDecomposition< T >::copy | ( | const choleskyDecomposition< T > & | other | ) |
copy data of "other" functor.
other | the functor to be copied |
Reimplemented from lti::functor.
const parameters& lti::choleskyDecomposition< T >::getParameters | ( | ) | const |
returns used parameters
Reimplemented from lti::functor.
virtual const char* lti::choleskyDecomposition< T >::getTypeName | ( | ) | const [virtual] |
returns the name of this type ("choleskyDecomposition")
Reimplemented from lti::linearAlgebraFunctor.
choleskyDecomposition& lti::choleskyDecomposition< T >::operator= | ( | const choleskyDecomposition< T > & | other | ) |