ROOT  6.06/09
Reference Guide
TDecompBK.h
Go to the documentation of this file.
1 // @(#)root/matrix:$Id$
2 // Authors: Fons Rademakers, Eddy Offermann Sep 2004
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_TDecompBK
13 #define ROOT_TDecompBK
14 
15 ///////////////////////////////////////////////////////////////////////////
16 // //
17 // Bunch-Kaufman Decomposition class //
18 // //
19 ///////////////////////////////////////////////////////////////////////////
20 
21 #ifndef ROOT_TDecompBase
22 #include "TDecompBase.h"
23 #endif
24 #ifndef ROOT_TMatrixDSym
25 #include "TMatrixDSym.h"
26 #endif
27 #ifndef ROOT_TVectorD
28 #include "TVectorD.h"
29 #endif
30 
31 class TDecompBK : public TDecompBase
32 {
33 protected :
34 
35  Int_t fNIpiv; // size of row permutation index
36  Int_t *fIpiv; //[fNIpiv] row permutation index
37  TMatrixD fU; // decomposed matrix so that a = u d u^T
38 
39  virtual const TMatrixDBase &GetDecompMatrix() const { return fU; }
40 
41 public :
42 
43  TDecompBK();
44  explicit TDecompBK(Int_t nrows);
45  TDecompBK(Int_t row_lwb,Int_t row_upb);
46  TDecompBK(const TMatrixDSym &m,Double_t tol = 0.0);
47  TDecompBK(const TDecompBK &another);
48  virtual ~TDecompBK() {if (fIpiv) delete [] fIpiv; fIpiv = 0; }
49 
50  virtual Int_t GetNrows () const { return fU.GetNrows(); }
51  virtual Int_t GetNcols () const { return fU.GetNcols(); }
52  const TMatrixD &GetU () { if ( !TestBit(kDecomposed) ) Decompose();
53  return fU; }
54 
55  virtual void SetMatrix (const TMatrixDSym &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) { return Solve(b); }
62  virtual TVectorD TransSolve (const TVectorD& b,Bool_t &ok) { TVectorD x = b; ok = Solve(x); return x; }
63  virtual Bool_t TransSolve ( TMatrixDColumn &b) { return Solve(b); }
64  virtual void Det (Double_t &/*d1*/,Double_t &/*d2*/)
65  { MayNotUse("Det(Double_t&,Double_t&)"); }
66 
69  TMatrixDSym Invert () { Bool_t status; return Invert(status); }
70 
71  void Print(Option_t *opt ="") const; // *MENU*
72 
73  TDecompBK &operator= (const TDecompBK &source);
74 
75  ClassDef(TDecompBK,1) // Matrix Decomposition Bunch-Kaufman
76 };
77 
78 #endif
virtual Bool_t TransSolve(TMatrixDColumn &b)
Definition: TDecompBK.h:63
const char Option_t
Definition: RtypesCore.h:62
Int_t fNIpiv
Definition: TDecompBK.h:35
virtual Bool_t Solve(TVectorD &b)
Solve Ax=b assuming the BK form of A is stored in fU . Solution returned in b.
Definition: TDecompBK.cxx:335
void MayNotUse(const char *method) const
Use this method to signal that a method (defined in a base class) may not be called in a derived clas...
Definition: TObject.cxx:971
virtual const TMatrixDBase & GetDecompMatrix() const
Definition: TDecompBK.h:39
double inv(double x)
For comparisons.
Definition: inv.h:58
const TMatrixD & GetU()
Definition: TDecompBK.h:52
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
Double_t x[n]
Definition: legend1.C:17
#define ClassDef(name, id)
Definition: Rtypes.h:254
void Print(Option_t *opt="") const
Print the class members.
Definition: TDecompBK.cxx:644
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
const double tol
virtual TVectorD TransSolve(const TVectorD &b, Bool_t &ok)
Definition: TDecompBK.h:62
TDecompBK & operator=(const TDecompBK &source)
Assigment operator.
Definition: TDecompBK.cxx:656
TMatrixDSym Invert()
Definition: TDecompBK.h:69
virtual Int_t GetNrows() const
Definition: TDecompBK.h:50
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
TMarker * m
Definition: textangle.C:8
virtual void Det(Double_t &, Double_t &)
Matrix determinant det = d1*TMath::Power(2.,d2)
Definition: TDecompBK.h:64
TMatrixD fU
Definition: TDecompBK.h:37
double Double_t
Definition: RtypesCore.h:55
virtual void SetMatrix(const TMatrixDSym &a)
Set the matrix to be decomposed, decomposition status is reset.
Definition: TDecompBK.cxx:310
Int_t * fIpiv
Definition: TDecompBK.h:36
Int_t GetNcols() const
Definition: TMatrixTBase.h:137
virtual ~TDecompBK()
Definition: TDecompBK.h:48
virtual Int_t GetNcols() const
Definition: TDecompBK.h:51
virtual Bool_t TransSolve(TVectorD &b)
Definition: TDecompBK.h:61
virtual TVectorD Solve(const TVectorD &b, Bool_t &ok)
Definition: TDecompBK.h:59
virtual Bool_t Decompose()
Matrix A is decomposed in components U and D so that A = U*D*U^T If the decomposition succeeds...
Definition: TDecompBK.cxx:127