|
ROOT 6.13/01 Reference Guide |
QR Decomposition class.
Decompose a general (m x n) matrix A into A = fQ fR H where
If row/column index of A starts at (rowLwb,colLwb) then the decomposed matrices start from :
Errors arise from formation of reflectors i.e. singularity . Note it attempts to handle the cases where the nRow <= nCol .
Definition at line 25 of file TDecompQRH.h.
Public Types | |
| enum | { kWorkMax = 100 } |
Public Member Functions | |
| TDecompQRH () | |
| TDecompQRH (Int_t nrows, Int_t ncols) | |
| Constructor for (nrows x ncols) matrix. More... | |
| TDecompQRH (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... | |
| TDecompQRH (const TMatrixD &m, Double_t tol=0.0) | |
| Constructor for general matrix A . More... | |
| TDecompQRH (const TDecompQRH &another) | |
| Copy constructor. More... | |
| virtual | ~TDecompQRH () |
| virtual Bool_t | Decompose () |
| QR decomposition of matrix a by Householder transformations, see Golub & Loan first edition p41 & Sec 6.2. More... | |
| virtual void | Det (Double_t &d1, Double_t &d2) |
| This routine calculates the absolute (!) value of the determinant |det| = d1*TMath::Power(2.,d2) More... | |
| virtual Int_t | GetNcols () const |
| virtual Int_t | GetNrows () const |
| virtual const TMatrixD & | GetQ () |
| virtual const TMatrixD & | GetR () |
| virtual const TVectorD & | GetUp () |
| virtual const TVectorD & | GetW () |
| 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 () |
| TDecompQRH & | operator= (const TDecompQRH &source) |
| Assignment operator. More... | |
| void | Print (Option_t *opt="") const |
| Print the 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 QR form of A is stored in fR,fQ and fW, but assume b has not been transformed. More... | |
| virtual TVectorD | Solve (const TVectorD &b, Bool_t &ok) |
| virtual Bool_t | Solve (TMatrixDColumn &b) |
| Solve Ax=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed. More... | |
| virtual Bool_t | TransSolve (TVectorD &b) |
| Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed. More... | |
| virtual TVectorD | TransSolve (const TVectorD &b, Bool_t &ok) |
| virtual Bool_t | TransSolve (TMatrixDColumn &b) |
| Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed. More... | |
Public Member Functions inherited from TDecompBase | |
| TDecompBase () | |
| Default constructor. More... | |
| TDecompBase (const TDecompBase &another) | |
| Copy constructor. More... | |
| virtual | ~TDecompBase () |
| virtual Double_t | Condition () |
| Matrix condition number. More... | |
| 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... | |
| TDecompBase & | operator= (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 TMatrixDBase & | GetDecompMatrix () 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 | QRH (TMatrixD &q, TVectorD &diagR, TVectorD &up, TVectorD &w, Double_t tol) |
| Decomposition function . 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 | |
| TMatrixD | fQ |
| TMatrixD | fR |
| TVectorD | fUp |
| TVectorD | fW |
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 <TDecompQRH.h>
| anonymous enum |
| Enumerator | |
|---|---|
| kWorkMax | |
Definition at line 41 of file TDecompQRH.h.
|
inline |
Definition at line 43 of file TDecompQRH.h.
| TDecompQRH::TDecompQRH | ( | Int_t | nrows, |
| Int_t | ncols | ||
| ) |
Constructor for (nrows x ncols) matrix.
Definition at line 48 of file TDecompQRH.cxx.
| TDecompQRH::TDecompQRH | ( | 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 69 of file TDecompQRH.cxx.
| TDecompQRH::TDecompQRH | ( | const TMatrixD & | m, |
| Double_t | tol = 0.0 |
||
| ) |
Constructor for general matrix A .
Definition at line 96 of file TDecompQRH.cxx.
| TDecompQRH::TDecompQRH | ( | const TDecompQRH & | another | ) |
Copy constructor.
Definition at line 130 of file TDecompQRH.cxx.
|
inlinevirtual |
Definition at line 48 of file TDecompQRH.h.
|
virtual |
QR decomposition of matrix a by Householder transformations, see Golub & Loan first edition p41 & Sec 6.2.
First fR is returned in upper triang of fQ and diagR. fQ returned in 'u-form' in lower triang of fQ and fW, the latter containing the "Householder betas". If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular
Implements TDecompBase.
Definition at line 143 of file TDecompQRH.cxx.
|
virtual |
This routine calculates the absolute (!) value of the determinant |det| = d1*TMath::Power(2.,d2)
Reimplemented from TDecompBase.
Definition at line 486 of file TDecompQRH.cxx.
|
inlineprotectedvirtual |
Implements TDecompBase.
Definition at line 37 of file TDecompQRH.h.
|
inlinevirtual |
Implements TDecompBase.
Definition at line 51 of file TDecompQRH.h.
|
inlinevirtual |
Implements TDecompBase.
Definition at line 50 of file TDecompQRH.h.
|
inlinevirtual |
Definition at line 52 of file TDecompQRH.h.
|
inlinevirtual |
Definition at line 54 of file TDecompQRH.h.
|
inlinevirtual |
Definition at line 56 of file TDecompQRH.h.
|
inlinevirtual |
Definition at line 58 of file TDecompQRH.h.
| Bool_t TDecompQRH::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 508 of file TDecompQRH.cxx.
| TMatrixD TDecompQRH::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 526 of file TDecompQRH.cxx.
|
inline |
Definition at line 74 of file TDecompQRH.h.
| TDecompQRH & TDecompQRH::operator= | ( | const TDecompQRH & | source | ) |
Assignment operator.
Definition at line 554 of file TDecompQRH.cxx.
| void TDecompQRH::Print | ( | Option_t * | opt = "" | ) | const |
Print the class members.
Definition at line 542 of file TDecompQRH.cxx.
|
staticprotected |
Decomposition function .
Definition at line 182 of file TDecompQRH.cxx.
|
virtual |
Set matrix to be decomposed.
Definition at line 215 of file TDecompQRH.cxx.
|
virtual |
Solve Ax=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed.
Solution returned in b.
Implements TDecompBase.
Definition at line 249 of file TDecompQRH.cxx.
Implements TDecompBase.
Definition at line 65 of file TDecompQRH.h.
|
virtual |
Solve Ax=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed.
Solution returned in b.
Implements TDecompBase.
Definition at line 304 of file TDecompQRH.cxx.
|
virtual |
Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed.
Solution returned in b.
Implements TDecompBase.
Definition at line 363 of file TDecompQRH.cxx.
Implements TDecompBase.
Definition at line 68 of file TDecompQRH.h.
|
virtual |
Solve A^T x=b assuming the QR form of A is stored in fR,fQ and fW, but assume b has not been transformed.
Solution returned in b.
Implements TDecompBase.
Definition at line 423 of file TDecompQRH.cxx.
|
protected |
Definition at line 30 of file TDecompQRH.h.
|
protected |
Definition at line 31 of file TDecompQRH.h.
|
protected |
Definition at line 32 of file TDecompQRH.h.
|
protected |
Definition at line 33 of file TDecompQRH.h.