Logo ROOT   6.13/01
Reference Guide
List of all members | Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
TDecompSVD Class Reference

Single Value Decomposition class.

For an (m x n) matrix A with m >= n, the singular value decomposition is an (m x m) orthogonal matrix fU, an (m x n) diagonal matrix fS, and an (n x n) orthogonal matrix fV so that A = U*S*V'.

If the row/column index of A starts at (rowLwb,colLwb) then the decomposed matrices/vectors start at :

fU : (rowLwb,colLwb)
fV : (colLwb,colLwb)
fSig : (colLwb)

The diagonal matrix fS is stored in the singular values vector fSig . The singular values, fSig[k] = S[k][k], are ordered so that fSig[0] >= fSig[1] >= ... >= fSig[n-1].

The singular value decomposition always exists, so the decomposition will (as long as m >=n) never fail. If m < n, the user should add sufficient zero rows to A , so that m == n

Here fTol is used to set the threshold on the minimum allowed value of the singular values: min_singular = fTol*max(fSig[i])

Definition at line 23 of file TDecompSVD.h.

Public Types

enum  { kWorkMax = 100 }
 

Public Member Functions

 TDecompSVD ()
 
 TDecompSVD (Int_t nrows, Int_t ncols)
 Constructor for (nrows x ncols) matrix. More...
 
 TDecompSVD (Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb)
 Constructor for ([row_lwb..row_upb] x [col_lwb..col_upb]) matrix. More...
 
 TDecompSVD (const TMatrixD &m, Double_t tol=0.0)
 Constructor for general matrix A . More...
 
 TDecompSVD (const TDecompSVD &another)
 Copy constructor. More...
 
virtual ~TDecompSVD ()
 
virtual Double_t Condition ()
 Matrix condition number. More...
 
virtual Bool_t Decompose ()
 SVD decomposition of matrix If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular. More...
 
virtual void Det (Double_t &d1, Double_t &d2)
 Matrix determinant det = d1*TMath::Power(2.,d2) More...
 
const TMatrixD GetMatrix ()
 Reconstruct the original matrix using the decomposition parts. More...
 
virtual Int_t GetNcols () const
 
virtual Int_t GetNrows () const
 
const TVectorDGetSig ()
 
const TMatrixDGetU ()
 
const TMatrixDGetV ()
 
Bool_t Invert (TMatrixD &inv)
 For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit The user should always supply a matrix of size (m x m) ! If m > n , only the (n x m) part of the returned (pseudo inverse) matrix should be used . More...
 
TMatrixD Invert (Bool_t &status)
 For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit (n x m) Ainv is returned . More...
 
TMatrixD Invert ()
 
TDecompSVDoperator= (const TDecompSVD &source)
 Assignment operator. More...
 
void Print (Option_t *opt="") const
 Print class members. More...
 
virtual void SetMatrix (const TMatrixD &a)
 Set matrix to be decomposed. More...
 
virtual Bool_t Solve (TVectorD &b)
 Solve Ax=b assuming the SVD form of A is stored . More...
 
virtual TVectorD Solve (const TVectorD &b, Bool_t &ok)
 
virtual Bool_t Solve (TMatrixDColumn &b)
 Solve Ax=b assuming the SVD form of A is stored . More...
 
virtual Bool_t TransSolve (TVectorD &b)
 Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b. More...
 
virtual TVectorD TransSolve (const TVectorD &b, Bool_t &ok)
 
virtual Bool_t TransSolve (TMatrixDColumn &b)
 Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b. More...
 
- Public Member Functions inherited from TDecompBase
 TDecompBase ()
 Default constructor. More...
 
 TDecompBase (const TDecompBase &another)
 Copy constructor. More...
 
virtual ~TDecompBase ()
 
Int_t GetColLwb () const
 
Double_t GetCondition () const
 
Double_t GetDet1 () const
 
Double_t GetDet2 () const
 
Int_t GetRowLwb () const
 
Double_t GetTol () const
 
virtual Bool_t MultiSolve (TMatrixD &B)
 Solve set of equations with RHS in columns of B. More...
 
TDecompBaseoperator= (const TDecompBase &source)
 Assignment operator. More...
 
void Print (Option_t *opt="") const
 Print class members. More...
 
Double_t SetTol (Double_t tol)
 

Protected Member Functions

virtual const TMatrixDBaseGetDecompMatrix () const
 
- Protected Member Functions inherited from TDecompBase
Int_t Hager (Double_t &est, Int_t iter=5)
 
void ResetStatus ()
 

Static Protected Member Functions

static Bool_t Bidiagonalize (TMatrixD &v, TMatrixD &u, TVectorD &sDiag, TVectorD &oDiag)
 Bidiagonalize the (m x n) - matrix a (stored in v) through a series of Householder transformations applied to the left (Q^T) and to the right (H) of a , so that A = Q . More...
 
static void Diag_1 (TMatrixD &v, TVectorD &sDiag, TVectorD &oDiag, Int_t k)
 Step 1 in the matrix diagonalization. More...
 
static void Diag_2 (TVectorD &sDiag, TVectorD &oDiag, Int_t k, Int_t l)
 Step 2 in the matrix diagonalization. More...
 
static void Diag_3 (TMatrixD &v, TMatrixD &u, TVectorD &sDiag, TVectorD &oDiag, Int_t k, Int_t l)
 Step 3 in the matrix diagonalization. More...
 
static Bool_t Diagonalize (TMatrixD &v, TMatrixD &u, TVectorD &sDiag, TVectorD &oDiag)
 Diagonalizes in an iterative fashion the bidiagonal matrix C as described through sDiag and oDiag, so that S' = U'^T . More...
 
static void SortSingular (TMatrixD &v, TMatrixD &u, TVectorD &sDiag)
 Perform a permutation transformation on the diagonal matrix S', so that matrix S'' = U''^T . More...
 
- Static Protected Member Functions inherited from TDecompBase
static void DiagProd (const TVectorD &diag, Double_t tol, Double_t &d1, Double_t &d2)
 

Protected Attributes

TVectorD fSig
 
TMatrixD fU
 
TMatrixD fV
 
- Protected Attributes inherited from TDecompBase
Int_t fColLwb
 
Double_t fCondition
 
Double_t fDet1
 
Double_t fDet2
 
Int_t fRowLwb
 
Double_t fTol
 

Additional Inherited Members

- Protected Types inherited from TDecompBase
enum  { kWorkMax = 100 }
 
enum  EMatrixDecompStat {
  kInit = BIT(14), kPatternSet = BIT(15), kValuesSet = BIT(16), kMatrixSet = BIT(17),
  kDecomposed = BIT(18), kDetermined = BIT(19), kCondition = BIT(20), kSingular = BIT(21)
}
 

#include <TDecompSVD.h>

Inheritance diagram for TDecompSVD:
[legend]

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
kWorkMax 

Definition at line 43 of file TDecompSVD.h.

Constructor & Destructor Documentation

◆ TDecompSVD() [1/5]

TDecompSVD::TDecompSVD ( )
inline

Definition at line 45 of file TDecompSVD.h.

◆ TDecompSVD() [2/5]

TDecompSVD::TDecompSVD ( Int_t  nrows,
Int_t  ncols 
)

Constructor for (nrows x ncols) matrix.

Definition at line 52 of file TDecompSVD.cxx.

◆ TDecompSVD() [3/5]

TDecompSVD::TDecompSVD ( Int_t  row_lwb,
Int_t  row_upb,
Int_t  col_lwb,
Int_t  col_upb 
)

Constructor for ([row_lwb..row_upb] x [col_lwb..col_upb]) matrix.

Definition at line 66 of file TDecompSVD.cxx.

◆ TDecompSVD() [4/5]

TDecompSVD::TDecompSVD ( const TMatrixD m,
Double_t  tol = 0.0 
)

Constructor for general matrix A .

Definition at line 85 of file TDecompSVD.cxx.

◆ TDecompSVD() [5/5]

TDecompSVD::TDecompSVD ( const TDecompSVD another)

Copy constructor.

Definition at line 114 of file TDecompSVD.cxx.

◆ ~TDecompSVD()

virtual TDecompSVD::~TDecompSVD ( )
inlinevirtual

Definition at line 50 of file TDecompSVD.h.

Member Function Documentation

◆ Bidiagonalize()

Bool_t TDecompSVD::Bidiagonalize ( TMatrixD v,
TMatrixD u,
TVectorD sDiag,
TVectorD oDiag 
)
staticprotected

Bidiagonalize the (m x n) - matrix a (stored in v) through a series of Householder transformations applied to the left (Q^T) and to the right (H) of a , so that A = Q .

C . H^T with matrix C bidiagonal. Q and H are orthogonal matrices .

Output: v - (n x n) - matrix H in the (n x n) part of v u - (m x m) - matrix Q^T sDiag - diagonal of the (m x n) C oDiag - off-diagonal elements of matrix C

Test code for the output:

const Int_t nRow = v.GetNrows();
const Int_t nCol = v.GetNcols();
TMatrixD H(v); H.ResizeTo(nCol,nCol);
TMatrixD E1(nCol,nCol); E1.UnitMatrix();
TMatrixD Ht(TMatrixDBase::kTransposed,H);
Bool_t ok = kTRUE;
ok &= VerifyMatrixIdentity(Ht * H,E1,kTRUE,1.0e-13);
ok &= VerifyMatrixIdentity(H * Ht,E1,kTRUE,1.0e-13);
TMatrixD E2(nRow,nRow); E2.UnitMatrix();
TMatrixD Qt(u);
TMatrixD Q(TMatrixDBase::kTransposed,Qt);
ok &= VerifyMatrixIdentity(Q * Qt,E2,kTRUE,1.0e-13);
TMatrixD C(nRow,nCol);
TMatrixDDiag(C) = sDiag;
for (Int_t i = 0; i < nCol-1; i++)
C(i,i+1) = oDiag(i+1);
TMatrixD A = Q*C*Ht;
ok &= VerifyMatrixIdentity(A,a,kTRUE,1.0e-13);

Definition at line 192 of file TDecompSVD.cxx.

◆ Condition()

Double_t TDecompSVD::Condition ( )
virtual

Matrix condition number.

Reimplemented from TDecompBase.

Definition at line 823 of file TDecompSVD.cxx.

◆ Decompose()

Bool_t TDecompSVD::Decompose ( )
virtual

SVD decomposition of matrix If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular.

Implements TDecompBase.

Definition at line 123 of file TDecompSVD.cxx.

◆ Det()

void TDecompSVD::Det ( Double_t &  d1,
Double_t &  d2 
)
virtual

Matrix determinant det = d1*TMath::Power(2.,d2)

Reimplemented from TDecompBase.

Definition at line 846 of file TDecompSVD.cxx.

◆ Diag_1()

void TDecompSVD::Diag_1 ( TMatrixD v,
TVectorD sDiag,
TVectorD oDiag,
Int_t  k 
)
staticprotected

Step 1 in the matrix diagonalization.

Definition at line 374 of file TDecompSVD.cxx.

◆ Diag_2()

void TDecompSVD::Diag_2 ( TVectorD sDiag,
TVectorD oDiag,
Int_t  k,
Int_t  l 
)
staticprotected

Step 2 in the matrix diagonalization.

Definition at line 398 of file TDecompSVD.cxx.

◆ Diag_3()

void TDecompSVD::Diag_3 ( TMatrixD v,
TMatrixD u,
TVectorD sDiag,
TVectorD oDiag,
Int_t  k,
Int_t  l 
)
staticprotected

Step 3 in the matrix diagonalization.

Definition at line 416 of file TDecompSVD.cxx.

◆ Diagonalize()

Bool_t TDecompSVD::Diagonalize ( TMatrixD v,
TMatrixD u,
TVectorD sDiag,
TVectorD oDiag 
)
staticprotected

Diagonalizes in an iterative fashion the bidiagonal matrix C as described through sDiag and oDiag, so that S' = U'^T .

C . V' is diagonal. U' and V' are orthogonal matrices .

Output:

  • v - (n x n) - matrix H . V' in the (n x n) part of v
  • u - (m x m) - matrix U'^T . Q^T
  • sDiag - diagonal of the (m x n) S'

    return convergence flag: 0 -> no convergence 1 -> convergence

Test code for the output:

const Int_t nRow = v.GetNrows();
const Int_t nCol = v.GetNcols();
TMatrixD tmp = v; tmp.ResizeTo(nCol,nCol);
TMatrixD Vprime = Ht*tmp;
TMatrixD Vprimet(TMatrixDBase::kTransposed,Vprime);
TMatrixD Uprimet = u*Q;
TMatrixD Uprime(TMatrixDBase::kTransposed,Uprimet);
TMatrixD Sprime(nRow,nCol);
TMatrixDDiag(Sprime) = sDiag;
ok &= VerifyMatrixIdentity(Uprimet * C * Vprime,Sprime,kTRUE,1.0e-13);
ok &= VerifyMatrixIdentity(Q*Uprime * Sprime * Vprimet * Ht,a,kTRUE,1.0e-13);

Definition at line 307 of file TDecompSVD.cxx.

◆ GetDecompMatrix()

virtual const TMatrixDBase& TDecompSVD::GetDecompMatrix ( ) const
inlineprotectedvirtual

Implements TDecompBase.

Definition at line 39 of file TDecompSVD.h.

◆ GetMatrix()

const TMatrixD TDecompSVD::GetMatrix ( )

Reconstruct the original matrix using the decomposition parts.

Definition at line 557 of file TDecompSVD.cxx.

◆ GetNcols()

Int_t TDecompSVD::GetNcols ( ) const
virtual

Implements TDecompBase.

Definition at line 870 of file TDecompSVD.cxx.

◆ GetNrows()

Int_t TDecompSVD::GetNrows ( ) const
virtual

Implements TDecompBase.

Definition at line 865 of file TDecompSVD.cxx.

◆ GetSig()

const TVectorD& TDecompSVD::GetSig ( )
inline

Definition at line 59 of file TDecompSVD.h.

◆ GetU()

const TMatrixD& TDecompSVD::GetU ( )
inline

Definition at line 55 of file TDecompSVD.h.

◆ GetV()

const TMatrixD& TDecompSVD::GetV ( )
inline

Definition at line 57 of file TDecompSVD.h.

◆ Invert() [1/3]

Bool_t TDecompSVD::Invert ( TMatrixD inv)

For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit The user should always supply a matrix of size (m x m) ! If m > n , only the (n x m) part of the returned (pseudo inverse) matrix should be used .

Definition at line 881 of file TDecompSVD.cxx.

◆ Invert() [2/3]

TMatrixD TDecompSVD::Invert ( Bool_t &  status)

For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit (n x m) Ainv is returned .

Definition at line 903 of file TDecompSVD.cxx.

◆ Invert() [3/3]

TMatrixD TDecompSVD::Invert ( )
inline

Definition at line 82 of file TDecompSVD.h.

◆ operator=()

TDecompSVD & TDecompSVD::operator= ( const TDecompSVD source)

Assignment operator.

Definition at line 930 of file TDecompSVD.cxx.

◆ Print()

void TDecompSVD::Print ( Option_t *  opt = "") const

Print class members.

Definition at line 919 of file TDecompSVD.cxx.

◆ SetMatrix()

void TDecompSVD::SetMatrix ( const TMatrixD a)
virtual

Set matrix to be decomposed.

Definition at line 582 of file TDecompSVD.cxx.

◆ Solve() [1/3]

Bool_t TDecompSVD::Solve ( TVectorD b)
virtual

Solve Ax=b assuming the SVD form of A is stored .

Solution returned in b. If A is of size (m x n), input vector b should be of size (m), however, the solution, returned in b, will be in the first (n) elements .

For m > n , x is the least-squares solution of min(A . x - b)

Implements TDecompBase.

Definition at line 615 of file TDecompSVD.cxx.

◆ Solve() [2/3]

virtual TVectorD TDecompSVD::Solve ( const TVectorD b,
Bool_t &  ok 
)
inlinevirtual

Implements TDecompBase.

Definition at line 66 of file TDecompSVD.h.

◆ Solve() [3/3]

Bool_t TDecompSVD::Solve ( TMatrixDColumn cb)
virtual

Solve Ax=b assuming the SVD form of A is stored .

Solution returned in the matrix column cb b. If A is of size (m x n), input vector b should be of size (m), however, the solution, returned in b, will be in the first (n) elements .

For m > n , x is the least-squares solution of min(A . x - b)

Implements TDecompBase.

Definition at line 670 of file TDecompSVD.cxx.

◆ SortSingular()

void TDecompSVD::SortSingular ( TMatrixD v,
TMatrixD u,
TVectorD sDiag 
)
staticprotected

Perform a permutation transformation on the diagonal matrix S', so that matrix S'' = U''^T .

S' . V'' has diagonal elements ordered such that they do not increase.

Output:

  • v - (n x n) - matrix H . V' . V'' in the (n x n) part of v
  • u - (m x m) - matrix U''^T . U'^T . Q^T
  • sDiag - diagonal of the (m x n) S''

Definition at line 497 of file TDecompSVD.cxx.

◆ TransSolve() [1/3]

Bool_t TDecompSVD::TransSolve ( TVectorD b)
virtual

Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b.

Implements TDecompBase.

Definition at line 723 of file TDecompSVD.cxx.

◆ TransSolve() [2/3]

virtual TVectorD TDecompSVD::TransSolve ( const TVectorD b,
Bool_t &  ok 
)
inlinevirtual

Implements TDecompBase.

Definition at line 72 of file TDecompSVD.h.

◆ TransSolve() [3/3]

Bool_t TDecompSVD::TransSolve ( TMatrixDColumn b)
virtual

Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b.

Implements TDecompBase.

Definition at line 772 of file TDecompSVD.cxx.

Member Data Documentation

◆ fSig

TVectorD TDecompSVD::fSig
protected

Definition at line 30 of file TDecompSVD.h.

◆ fU

TMatrixD TDecompSVD::fU
protected

Definition at line 28 of file TDecompSVD.h.

◆ fV

TMatrixD TDecompSVD::fV
protected

Definition at line 29 of file TDecompSVD.h.

Libraries for TDecompSVD:
[legend]

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