Logo ROOT   6.18/05
Reference Guide
TDecompLU.h
Go to the documentation of this file.
1// @(#)root/matrix:$Id$
2// Authors: Fons Rademakers, Eddy Offermann Dec 2003
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#ifndef ROOT_TDecompLU
13#define ROOT_TDecompLU
14
15///////////////////////////////////////////////////////////////////////////
16// //
17// LU Decomposition class //
18// //
19///////////////////////////////////////////////////////////////////////////
20
21#include "TDecompBase.h"
22
23class TDecompLU : public TDecompBase
24{
25protected :
26
27 Int_t fImplicitPivot; // control to determine implicit row scale before
28 // deciding on the pivot (Crout method)
29 Int_t fNIndex; // size of row permutation index
30 Int_t *fIndex; //[fNIndex] row permutation index
31 Double_t fSign; // = +/- 1 reflecting even/odd row permutations, resp.
32 TMatrixD fLU; // decomposed matrix so that a = l u where
33 // l is stored lower left and u upper right side
34
35 static Bool_t DecomposeLUCrout(TMatrixD &lu,Int_t *index,Double_t &sign,Double_t tol,Int_t &nrZeros);
36 static Bool_t DecomposeLUGauss(TMatrixD &lu,Int_t *index,Double_t &sign,Double_t tol,Int_t &nrZeros);
37
38 virtual const TMatrixDBase &GetDecompMatrix() const { return fLU; }
39
40public :
41
42 TDecompLU();
43 explicit TDecompLU(Int_t nrows);
44 TDecompLU(Int_t row_lwb,Int_t row_upb);
45 TDecompLU(const TMatrixD &m,Double_t tol = 0.0,Int_t implicit = 1);
46 TDecompLU(const TDecompLU &another);
47 virtual ~TDecompLU() {if (fIndex) delete [] fIndex; fIndex = 0; }
48
49 const TMatrixD GetMatrix ();
50 virtual Int_t GetNrows () const { return fLU.GetNrows(); }
51 virtual Int_t GetNcols () const { return fLU.GetNcols(); }
52 const TMatrixD &GetLU () { if ( !TestBit(kDecomposed) ) Decompose();
53 return fLU; }
54
55 virtual void SetMatrix (const TMatrixD &a);
56
57 virtual Bool_t Decompose ();
58 virtual Bool_t Solve ( TVectorD &b);
59 virtual TVectorD Solve (const TVectorD& b,Bool_t &ok) { TVectorD x = b; ok = Solve(x); return x; }
60 virtual Bool_t Solve ( TMatrixDColumn &b);
61 virtual Bool_t TransSolve ( TVectorD &b);
62 virtual TVectorD TransSolve (const TVectorD& b,Bool_t &ok) { TVectorD x = b; ok = TransSolve(x); return x; }
64 virtual void Det (Double_t &d1,Double_t &d2);
65
66 static Bool_t InvertLU (TMatrixD &a,Double_t tol,Double_t *det=0);
68 TMatrixD Invert (Bool_t &status);
69 TMatrixD Invert () { Bool_t status; return Invert(status); }
70
71 void Print(Option_t *opt ="") const; // *MENU*
72
73 TDecompLU &operator= (const TDecompLU &source);
74
75 ClassDef(TDecompLU,1) // Matrix Decompositition LU
76};
77
78#endif
#define b(i)
Definition: RSha256.hxx:100
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const char Option_t
Definition: RtypesCore.h:62
#define ClassDef(name, id)
Definition: Rtypes.h:326
Decomposition Base class.
Definition: TDecompBase.h:34
LU Decomposition class.
Definition: TDecompLU.h:24
virtual void Det(Double_t &d1, Double_t &d2)
Calculate determinant det = d1*TMath::Power(2.,d2)
Definition: TDecompLU.cxx:509
virtual 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.
Definition: TDecompLU.cxx:233
const TMatrixD & GetLU()
Definition: TDecompLU.h:52
virtual void SetMatrix(const TMatrixD &a)
Set matrix to be decomposed.
Definition: TDecompLU.cxx:202
static Bool_t DecomposeLUGauss(TMatrixD &lu, Int_t *index, Double_t &sign, Double_t tol, Int_t &nrZeros)
LU decomposition using Gaussian Elimination with partial pivoting (See Golub & Van Loan,...
Definition: TDecompLU.cxx:708
virtual Int_t GetNcols() const
Definition: TDecompLU.h:51
TDecompLU()
Default constructor.
Definition: TDecompLU.cxx:45
static Bool_t InvertLU(TMatrixD &a, Double_t tol, Double_t *det=0)
Calculate matrix inversion through in place forward/backward substitution.
Definition: TDecompLU.cxx:775
virtual 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.
Definition: TDecompLU.cxx:371
virtual Int_t GetNrows() const
Definition: TDecompLU.h:50
Double_t fSign
Definition: TDecompLU.h:31
virtual ~TDecompLU()
Definition: TDecompLU.h:47
const TMatrixD GetMatrix()
Reconstruct the original matrix using the decomposition parts.
Definition: TDecompLU.cxx:151
Int_t fNIndex
Definition: TDecompLU.h:29
virtual TVectorD TransSolve(const TVectorD &b, Bool_t &ok)
Definition: TDecompLU.h:62
static 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.
Definition: TDecompLU.cxx:599
TMatrixD Invert()
Definition: TDecompLU.h:69
virtual Bool_t Decompose()
Matrix A is decomposed in components U and L so that P * A = U * L If the decomposition succeeds,...
Definition: TDecompLU.cxx:126
TMatrixD fLU
Definition: TDecompLU.h:32
void Print(Option_t *opt="") const
Print internals of this object.
Definition: TDecompLU.cxx:559
virtual TVectorD Solve(const TVectorD &b, Bool_t &ok)
Definition: TDecompLU.h:59
Int_t * fIndex
Definition: TDecompLU.h:30
Int_t fImplicitPivot
Definition: TDecompLU.h:27
virtual const TMatrixDBase & GetDecompMatrix() const
Definition: TDecompLU.h:38
TDecompLU & operator=(const TDecompLU &source)
assignment operator
Definition: TDecompLU.cxx:573
Int_t GetNrows() const
Definition: TMatrixTBase.h:124
Int_t GetNcols() const
Definition: TMatrixTBase.h:127
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
Double_t x[n]
Definition: legend1.C:17
void inv(rsa_NUMBER *, rsa_NUMBER *, rsa_NUMBER *)
Definition: rsaaux.cxx:949
auto * m
Definition: textangle.C:8
auto * a
Definition: textangle.C:12