CMS 3D CMS Logo

Classes | Public Member Functions | Private Attributes

ROOT::Math::CholeskyDecomp< F, N > Class Template Reference

class to compute the Cholesky decomposition of a matrix More...

#include <CholeskyDecomp.h>

List of all members.

Classes

class  PackedArrayAdapter
 adapter for packed arrays (to SMatrix indexing conventions) More...

Public Member Functions

template<class M >
 CholeskyDecomp (const M &m)
 perform a Cholesky decomposition
template<typename G >
 CholeskyDecomp (G *m)
 perform a Cholesky decomposition
template<typename G >
bool Invert (G *m) const
 place the inverse into m
template<class M >
bool Invert (M &m) const
 place the inverse into m
bool ok () const
 returns true if decomposition was successful
 operator bool () const
 returns true if decomposition was successful
template<class V >
bool Solve (V &rhs) const
 solves a linear system for the given right hand side

Private Attributes

fL [N *(N+1)/2]
 lower triangular matrix L
bool fOk
 flag indicating a successful decomposition

Detailed Description

template<class F, unsigned N>
class ROOT::Math::CholeskyDecomp< F, N >

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:

 // let m be a symmetric positive definite SMatrix (use type float
 // for internal computations, matrix size is 4x4)
 CholeskyDecomp<float, 4> decomp(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);
 }

Definition at line 63 of file CholeskyDecomp.h.


Constructor & Destructor Documentation

template<class F, unsigned N>
template<class M >
ROOT::Math::CholeskyDecomp< F, N >::CholeskyDecomp ( 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 96 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fL, and ROOT::Math::CholeskyDecomp< F, N >::fOk.

   {
      using CholeskyDecompHelpers::_decomposer;
      fOk = _decomposer<F, N, M>()(fL, m);
   }
template<class F, unsigned N>
template<typename G >
ROOT::Math::CholeskyDecomp< F, N >::CholeskyDecomp ( 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 113 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fL, and ROOT::Math::CholeskyDecomp< F, N >::fOk.

   {
      using CholeskyDecompHelpers::_decomposer;
      fOk = _decomposer<F, N, PackedArrayAdapter<G> >()(
         fL, PackedArrayAdapter<G>(m));
   }

Member Function Documentation

template<class F, unsigned N>
template<class M >
bool ROOT::Math::CholeskyDecomp< F, N >::Invert ( M &  m) const [inline]

place the inverse into m

place the inverse into m

This is the method to use with an SMatrix.

Returns:
if the decomposition was successful

Definition at line 149 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fL, and ROOT::Math::CholeskyDecomp< F, N >::fOk.

Referenced by invertPosDefMatrix().

   {
      using CholeskyDecompHelpers::_inverter;
      if (fOk) _inverter<F,N,M>()(m, fL); return fOk;
   }
template<class F, unsigned N>
template<typename G >
bool ROOT::Math::CholeskyDecomp< F, N >::Invert ( G *  m) const [inline]

place the inverse into m

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 166 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fL, and ROOT::Math::CholeskyDecomp< F, N >::fOk.

   {
      using CholeskyDecompHelpers::_inverter;
      if (fOk) {
         PackedArrayAdapter<G> adapted(m);
         _inverter<F,N,PackedArrayAdapter<G> >()(adapted, fL);
      }
      return fOk;
   }
template<class F, unsigned N>
bool ROOT::Math::CholeskyDecomp< F, N >::ok ( ) const [inline]

returns true if decomposition was successful

Returns:
true if decomposition was successful

Definition at line 122 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fOk.

{ return fOk; }
template<class F, unsigned N>
ROOT::Math::CholeskyDecomp< F, N >::operator bool ( ) const [inline]

returns true if decomposition was successful

Returns:
true if decomposition was successful

Definition at line 125 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fOk.

{ return fOk; }
template<class F, unsigned N>
template<class V >
bool ROOT::Math::CholeskyDecomp< F, N >::Solve ( V &  rhs) const [inline]

solves a linear system for the given right hand side

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 136 of file CholeskyDecomp.h.

References ROOT::Math::CholeskyDecomp< F, N >::fL, and ROOT::Math::CholeskyDecomp< F, N >::fOk.

   {
      using CholeskyDecompHelpers::_solver;
      if (fOk) _solver<F,N,V>()(rhs, fL); return fOk;
   }

Member Data Documentation

template<class F, unsigned N>
F ROOT::Math::CholeskyDecomp< F, N >::fL[N *(N+1)/2] [private]

lower triangular matrix L

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

Definition at line 69 of file CholeskyDecomp.h.

Referenced by ROOT::Math::CholeskyDecomp< F, N >::CholeskyDecomp(), ROOT::Math::CholeskyDecomp< F, N >::Invert(), and ROOT::Math::CholeskyDecomp< F, N >::Solve().

template<class F, unsigned N>
bool ROOT::Math::CholeskyDecomp< F, N >::fOk [private]