ROOT logo
ROOT » MATH » MATRIX » TDecompLU

class TDecompLU: public TDecompBase


LU Decomposition class

Decompose  a general n x n matrix A into P A = L U

where P is a permutation matrix, L is unit lower triangular and U
is upper triangular.
L is stored in the strict lower triangular part of the matrix fLU.
The diagonal elements of L are unity and are not stored.
U is stored in the diagonal and upper triangular part of the matrix
fU.
P is stored in the index array fIndex : j = fIndex[i] indicates that
row j and row i should be swapped .

fSign gives the sign of the permutation, (-1)^n, where n is the
number of interchanges in the permutation.

fLU has the same indexing range as matrix A .

The decomposition fails if a diagonal element of abs(fLU) is == 0,
The matrix fUL is made invalid .


Function Members (Methods)

public:
TDecompLU()
TDecompLU(Int_t nrows)
TDecompLU(const TDecompLU& another)
TDecompLU(Int_t row_lwb, Int_t row_upb)
TDecompLU(const TMatrixD& m, Double_t tol = 0.0, Int_t implicit = 1)
virtual~TDecompLU()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTObject::Clear(Option_t* = "")
virtual TObject*TObject::Clone(const char* newname = "") const
virtual Int_tTObject::Compare(const TObject* obj) const
virtual Double_tTDecompBase::Condition()
virtual voidTObject::Copy(TObject& object) const
virtual Bool_tDecompose()
virtual voidTObject::Delete(Option_t* option = "")MENU
virtual voidDet(Double_t& d1, Double_t& d2)
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
Int_tTDecompBase::GetColLwb() const
Double_tTDecompBase::GetCondition() const
Double_tTDecompBase::GetDet1() const
Double_tTDecompBase::GetDet2() const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
const TMatrixD&GetLU()
const TMatrixDGetMatrix()
virtual const char*TObject::GetName() const
virtual Int_tGetNcols() const
virtual Int_tGetNrows() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
Int_tTDecompBase::GetRowLwb() const
virtual const char*TObject::GetTitle() const
Double_tTDecompBase::GetTol() const
virtual UInt_tTObject::GetUniqueID() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTObject::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() constMENU
TMatrixDInvert()
Bool_tInvert(TMatrixD& inv)
TMatrixDInvert(Bool_t& status)
voidTObject::InvertBit(UInt_t f)
static Bool_tInvertLU(TMatrixD& a, Double_t tol, Double_t* det = 0)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTObject::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTObject::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
virtual Bool_tTDecompBase::MultiSolve(TMatrixD& B)
virtual Bool_tTObject::Notify()
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
TDecompLU&operator=(const TDecompLU& source)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidPrint(Option_t* opt = "") constMENU
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(basic_ostream<char,char_traits<char> >& out, Option_t* option = "")
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
virtual voidSetMatrix(const TMatrixD& a)
static voidTObject::SetObjectStat(Bool_t stat)
Double_tTDecompBase::SetTol(Double_t newTol)
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual Bool_tSolve(TVectorD& b)
virtual Bool_tSolve(TMatrixDColumn& b)
virtual TVectorDSolve(const TVectorD& b, Bool_t& ok)
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual Bool_tTransSolve(TVectorD& b)
virtual Bool_tTransSolve(TMatrixDColumn& b)
virtual TVectorDTransSolve(const TVectorD& b, Bool_t& ok)
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
static Bool_tDecomposeLUCrout(TMatrixD& lu, Int_t* index, Double_t& sign, Double_t tol, Int_t& nrZeros)
static Bool_tDecomposeLUGauss(TMatrixD& lu, Int_t* index, Double_t& sign, Double_t tol, Int_t& nrZeros)
static voidTDecompBase::DiagProd(const TVectorD& diag, Double_t tol, Double_t& d1, Double_t& d2)
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
virtual const TMatrixDBase&GetDecompMatrix() const
Int_tTDecompBase::Hager(Double_t& est, Int_t iter = 5)
voidTObject::MakeZombie()
voidTDecompBase::ResetStatus()

Data Members

public:
enum TDecompBase::EMatrixDecompStat { kInit
kPatternSet
kValuesSet
kMatrixSet
kDecomposed
kDetermined
kCondition
kSingular
};
enum TDecompBase::[unnamed] { kWorkMax
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
protected:
Int_tTDecompBase::fColLwbColumn lower bound of decomposed matrix
Double_tTDecompBase::fConditionmatrix condition number
Double_tTDecompBase::fDet1determinant mantissa
Double_tTDecompBase::fDet2determinant exponent for powers of 2
Int_tfImplicitPivotcontrol to determine implicit row scale before
Int_t*fIndex[fNIndex] row permutation index
TMatrixDfLUdecomposed matrix so that a = l u where
Int_tfNIndexsize of row permutation index
Int_tTDecompBase::fRowLwbRow lower bound of decomposed matrix
Double_tfSign= +/- 1 reflecting even/odd row permutations, resp.
Double_tTDecompBase::fTolsqrt(epsilon); epsilon is smallest number number so that 1+epsilon > 1

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

TDecompLU()
 Default constructor
TDecompLU(Int_t nrows)
 Constructor for (nrows x nrows) matrix
TDecompLU(Int_t row_lwb, Int_t row_upb)
 Constructor for ([row_lwb..row_upb] x [row_lwb..row_upb]) matrix
TDecompLU(const TMatrixD& m, Double_t tol = 0.0, Int_t implicit = 1)
 Constructor for matrix a
TDecompLU(const TDecompLU &another)
 Copy constructor
Bool_t Decompose()
 Matrix A is decomposed in components U and L so that P * A = U * L
 If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular
const TMatrixD GetMatrix()
 Reconstruct the original matrix using the decomposition parts
void SetMatrix(const TMatrixD& a)
 Set matrix to be decomposed
Bool_t Solve(TVectorD &b)
 Solve Ax=b assuming the LU form of A is stored in fLU, but assume b has *not*
 been transformed.  Solution returned in b.
Bool_t Solve(TMatrixDColumn &cb)
 Solve Ax=b assuming the LU form of A is stored in fLU, but assume b has *not*
 been transformed.  Solution returned in b.
Bool_t TransSolve(TVectorD &b)
 Solve A^T x=b assuming the LU form of A^T is stored in fLU, but assume b has *not*
 been transformed.  Solution returned in b.
Bool_t TransSolve(TMatrixDColumn &cb)
 Solve A^T x=b assuming the LU form of A^T is stored in fLU, but assume b has *not*
 been transformed.  Solution returned in b.
void Det(Double_t& d1, Double_t& d2)
 Calculate determinant det = d1*TMath::Power(2.,d2)
Bool_t Invert(TMatrixD &inv)
 For a matrix A(m,m), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
 (m x m) Ainv is returned .
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 .
void Print(Option_t* opt = "") const
Print internals of this object
TDecompLU & operator=(const TDecompLU& source)
assignement operator
Bool_t DecomposeLUCrout(TMatrixD& lu, Int_t* index, Double_t& sign, Double_t tol, Int_t& nrZeros)
 Crout/Doolittle algorithm of LU decomposing a square matrix, with implicit partial
 pivoting.  The decomposition is stored in fLU: U is explicit in the upper triag
 and L is in multiplier form in the subdiagionals .
 Row permutations are mapped out in fIndex. fSign, used for calculating the
 determinant, is +/- 1 for even/odd row permutations. .
Bool_t DecomposeLUGauss(TMatrixD& lu, Int_t* index, Double_t& sign, Double_t tol, Int_t& nrZeros)
 LU decomposition using Gaussain Elimination with partial pivoting (See Golub &
 Van Loan, Matrix Computations, Algorithm 3.4.1) of a square matrix .
 The decomposition is stored in fLU: U is explicit in the upper triag and L is in
 multiplier form in the subdiagionals . Row permutations are mapped out in fIndex.
 fSign, used for calculating the determinant, is +/- 1 for even/odd row permutations.
 Since this algorithm uses partial pivoting without scaling like in Crout/Doolitle.
 it is somewhat faster but less precise .
Bool_t InvertLU(TMatrixD& a, Double_t tol, Double_t* det = 0)
 Calculate matrix inversion through in place forward/backward substitution
const TMatrixDBase & GetDecompMatrix() const
{ return fLU; }
virtual ~TDecompLU()
{if (fIndex) delete [] fIndex; fIndex = 0; }
Int_t GetNrows() const
{ return fLU.GetNrows(); }
Int_t GetNcols() const
{ return fLU.GetNcols(); }
const TMatrixD & GetLU()
Bool_t Solve( TVectorD &b)
Bool_t TransSolve( TVectorD &b)
Bool_t Invert(TMatrixD &inv)