library: libMatrix
#include "TMatrixFBase.h"

TMatrixFBase


class description - source file - inheritance tree (.pdf)

class TMatrixFBase : public TObject

Inheritance Chart:
TObject
<-
TMatrixFBase
<-
TMatrixF
<-
TMatrix
TMatrixFSym
 
    This is an abstract class, constructors will not be documented.
    Look at the header to check for available constructors.

    private:
Float_t* GetElements() protected:
virtual void Allocate(Int_t nrows, Int_t ncols, Int_t row_lwb = 0, Int_t col_lwb = 0, Int_t init = 0, Int_t nr_nonzero = -1) public:
virtual ~TMatrixFBase() virtual TMatrixFBase& Abs() virtual TMatrixFBase& Apply(const TElementActionF& action) virtual TMatrixFBase& Apply(const TElementPosActionF& action) static TClass* Class() virtual void Clear(Option_t* option = "") virtual Float_t ColNorm() const virtual Double_t Determinant() const virtual void Determinant(Double_t& d1, Double_t& d2) const virtual void Draw(Option_t* option = "") virtual Float_t E2Norm() const virtual void ExtractRow(Int_t row, Int_t col, Float_t* v, Int_t n = -1) const virtual const Int_t* GetColIndexArray() const virtual Int_t* GetColIndexArray() Int_t GetColLwb() const Int_t GetColUpb() const virtual void GetMatrix2Array(Float_t* data, Option_t* option = "") const virtual const Float_t* GetMatrixArray() const virtual Float_t* GetMatrixArray() Int_t GetNcols() const Int_t GetNoElements() const Int_t GetNrows() const virtual const Int_t* GetRowIndexArray() const virtual Int_t* GetRowIndexArray() Int_t GetRowLwb() const Int_t GetRowUpb() const virtual TMatrixFBase& GetSub(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, TMatrixFBase& target, Option_t* option = "S") const Float_t GetTol() const virtual TMatrixFBase& InsertRow(Int_t row, Int_t col, const Float_t* v, Int_t n = -1) void Invalidate() virtual TClass* IsA() const Bool_t IsOwner() const Bool_t IsSymmetric() const Bool_t IsValid() const void MakeValid() virtual Float_t Max() const virtual Float_t Min() const virtual Int_t NonZeros() const Float_t Norm1() const virtual TMatrixFBase& NormByDiag(const TVectorF& v, Option_t* option = "D") Float_t NormInf() const Bool_t operator!=(Float_t val) const virtual Float_t operator()(Int_t rown, Int_t coln) const virtual Float_t& operator()(Int_t rown, Int_t coln) Bool_t operator<(Float_t val) const Bool_t operator<=(Float_t val) const TMatrixFBase& operator=(const TMatrixFBase&) Bool_t operator==(Float_t val) const Bool_t operator>(Float_t val) const Bool_t operator>=(Float_t val) const virtual void Print(Option_t* name = "") const virtual TMatrixFBase& Randomize(Float_t alpha, Float_t beta, Double_t& seed) virtual TMatrixFBase& ResizeTo(Int_t nrows, Int_t ncols, Int_t nr_nonzeros = -1) virtual TMatrixFBase& ResizeTo(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Int_t nr_nonzeros = -1) virtual Float_t RowNorm() const virtual TMatrixFBase& SetColIndexArray(Int_t* data) virtual TMatrixFBase& SetMatrixArray(const Float_t* data, Option_t* option = "") virtual TMatrixFBase& SetRowIndexArray(Int_t* data) virtual TMatrixFBase& SetSub(Int_t row_lwb, Int_t col_lwb, const TMatrixFBase& source) Float_t SetTol(Float_t tol) virtual TMatrixFBase& Shift(Int_t row_shift, Int_t col_shift) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual TMatrixFBase& Sqr() virtual TMatrixFBase& Sqrt() virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b) virtual Float_t Sum() const virtual TMatrixFBase& UnitMatrix() virtual TMatrixFBase& Zero()

Data Members


    protected:
Int_t fNrows number of rows Int_t fNcols number of columns Int_t fRowLwb lower bound of the row index Int_t fColLwb lower bound of the col index Int_t fNelems number of elements in matrix Int_t fNrowIndex length of row index array (= fNrows+1) wich is only used for sparse matrices Float_t fTol sqrt(epsilon); epsilon is smallest number number so that 1+epsilon > 1 Bool_t fIsOwner !default kTRUE, when Use array kFALSE public:
static const enum TMatrixFBase:: kSizeMax static const enum TMatrixFBase:: kWorkMax static const TMatrixFBase::EMatrixStatusBits kStatus static const TMatrixFBase::EMatrixCreatorsOp1 kZero static const TMatrixFBase::EMatrixCreatorsOp1 kUnit static const TMatrixFBase::EMatrixCreatorsOp1 kTransposed static const TMatrixFBase::EMatrixCreatorsOp1 kInverted static const TMatrixFBase::EMatrixCreatorsOp1 kAtA static const TMatrixFBase::EMatrixCreatorsOp2 kMult static const TMatrixFBase::EMatrixCreatorsOp2 kTransposeMult static const TMatrixFBase::EMatrixCreatorsOp2 kInvMult static const TMatrixFBase::EMatrixCreatorsOp2 kMultTranspose static const TMatrixFBase::EMatrixCreatorsOp2 kPlus static const TMatrixFBase::EMatrixCreatorsOp2 kMinus

Class Description

                                                                      
 Linear Algebra Package                                               
 ----------------------                                               
                                                                      
 The present package implements all the basic algorithms dealing      
 with vectors, matrices, matrix columns, rows, diagonals, etc.        
 In addition eigen-Vector analysis and several matrix decomposition   
 have been added (LU,QRH,Cholesky and SVD) . The decompositions are   
 used in matrix inversion, equation solving.                          
                                                                      
 Matrix elements are arranged in memory in a ROW-wise fashion         
 fashion . For (n x m) matrices where n*m < kSizeMax (=25 currently)  
 storage space is avaialble on the stack, thus avoiding expensive     
 allocation/deallocation of heap space . However, this introduces of  
 course kSizeMax overhead for each matrix object . If this is an      
 issue recompile with a new appropriate value (>=0) for kSizeMax      
                                                                      
 Another way to assign and store matrix data is through Use           
 see for instance stressLinear.cxx file .                             
                                                                      
 Unless otherwise specified, matrix and vector indices always start   
 with 0, spanning up to the specified limit-1. However, there are     
 constructors to which one can specify aribtrary lower and upper      
 bounds, e.g. TMatrixF m(1,10,1,5) defines a matrix that ranges,      
 nad that can be addresses, from 1..10, 1..5 (a(1,1)..a(10,5)).       
                                                                      
 The present package provides all facilities to completely AVOID      
 returning matrices. Use "TMatrixF A(TMatrixF::kTransposed,B);"       
 and other fancy constructors as much as possible. If one really needs
 to return a matrix, return a TMatLazy object instead. The            
 conversion is completely transparent to the end user, e.g.           
 "TMatrixF m = THaarMatrixF(5);" and _is_ efficient.                  
                                                                      
 Since TMatrixF et al. are fully integrated in ROOT they of course    
 can be stored in a ROOT database.                                    
                                                                      
 For usage examples see $ROOTSYS/test/stressLinear.cxx                
                                                                      
 Acknowledgements                                                     
 ----------------                                                     
 1. Oleg E. Kiselyov                                                  
  First implementations were based on the his code . We have diverged 
  quite a bit since then but the ideas/code for lazy matrix and       
  "nested function" are 100% his .                                    
  You can see him and his code in action at http://okmij.org/ftp      
 2. Chris R. Birchenhall,                                             
  We adapted his idea of the implementation for the decomposition     
  classes instead of our messy installation of matrix inversion       
  His installation of matrix condition number, using an iterative     
  scheme using the Hage algorithm is worth looking at !               
  Chris has a nice writeup (matdoc.ps) on his matrix classes at       
   ftp://ftp.mcc.ac.uk/pub/matclass/                                  
 3. Mark Fischler and Steven Haywood of CLHEP                         
  They did the slave labor of spelling out all sub-determinants       
   for Cramer inversion  of (4x4),(5x5) and (6x6) matrices            
  The stack storage for small matrices was also taken from them       
 4. Roldan Pozo of TNT (http://math.nist.gov/tnt/)                    
  He converted the EISPACK routines for the eigen-vector analysis to  
  C++ . We started with his implementation                            
 5. Siegmund Brandt (http://siux00.physik.uni-siegen.de/~brandt/datan 
  We adapted his (very-well) documented SVD routines                  
                                                                      
 How to efficiently use this package                                  
 -----------------------------------                                  
                                                                      
 1. Never return complex objects (matrices or vectors)                
    Danger: For example, when the following snippet:                  
        TMatrixF foo(int n)                                           
        {                                                             
           TMatrixF foom(n,n); fill_in(foom); return foom;            
        }                                                             
        TMatrixF m = foo(5);                                          
    runs, it constructs matrix foo:foom, copies it onto stack as a    
    return value and destroys foo:foom. Return value (a matrix)       
    from foo() is then copied over to m (via a copy constructor),     
    and the return value is destroyed. So, the matrix constructor is  
    called 3 times and the destructor 2 times. For big matrices,      
    the cost of multiple constructing/copying/destroying of objects   
    may be very large. *Some* optimized compilers can cut down on 1   
    copying/destroying, but still it leaves at least two calls to     
    the constructor. Note, TMatLazy (see below) can construct         
    TMatrixF m "inplace", with only a _single_ call to the            
    constructor.                                                      
                                                                      
 2. Use "two-address instructions"                                    
        "void TMatrixF::operator += (const TMatrixF &B);"             
    as much as possible.                                              
    That is, to add two matrices, it's much more efficient to write   
        A += B;                                                       
    than                                                              
        TMatrixF C = A + B;                                           
    (if both operand should be preserved,                             
        TMatrixF C = A; C += B;                                       
    is still better).                                                 
                                                                      
 3. Use glorified constructors when returning of an object seems      
    inevitable:                                                       
        "TMatrixF A(TMatrixF::kTransposed,B);"                        
        "TMatrixF C(A,TMatrixF::kTransposeMult,B);"                   
                                                                      
    like in the following snippet (from $ROOTSYS/test/vmatrix.cxx)    
    that verifies that for an orthogonal matrix T, T'T = TT' = E.     
                                                                      
    TMatrixF haar = THaarMatrixF(5);                                  
    TMatrixF unit(TMatrixF::kUnit,haar);                              
    TMatrixF haar_t(TMatrixF::kTransposed,haar);                      
    TMatrixF hth(haar,TMatrixF::kTransposeMult,haar);                 
    TMatrixF hht(haar,TMatrixF::kMult,haar_t);                        
    TMatrixF hht1 = haar; hht1 *= haar_t;                             
    VerifyMatrixIdentity(unit,hth);                                   
    VerifyMatrixIdentity(unit,hht);                                   
    VerifyMatrixIdentity(unit,hht1);                                  
                                                                      
 4. Accessing row/col/diagonal of a matrix without much fuss          
    (and without moving a lot of stuff around):                       
                                                                      
        TMatrixF m(n,n); TVectorF v(n); TMatrixFDiag(m) += 4;         
        v = TMatrixFRow(m,0);                                         
        TMatrixColumn m1(m,1); m1(2) = 3; // the same as m(2,1)=3;    
    Note, constructing of, say, TMatrixFDiag does *not* involve any   
    copying of any elements of the source matrix.                     
                                                                      
 5. It's possible (and encouraged) to use "nested" functions          
    For example, creating of a Hilbert matrix can be done as follows: 
                                                                      
    void foo(const TMatrixF &m)                                       
    {                                                                 
      TMatrixF m1(TMatrixF::kZero,m);                                 
      struct MakeHilbert : public TElementPosActionF {                
        void Operation(Float_t &element) { element = 1./(fI+fJ-1); }  
      };                                                              
      m1.Apply(MakeHilbert());                                        
    }                                                                 
                                                                      
    of course, using a special method THilbertMatrixF() is            
    still more optimal, but not by a whole lot. And that's right,     
    class MakeHilbert is declared *within* a function and local to    
    that function. It means one can define another MakeHilbert class  
    (within another function or outside of any function, that is, in  
    the global scope), and it still will be OK. Note, this currently  
    is not yet supported by the interpreter CINT.                     
                                                                      
    Another example is applying of a simple function to each matrix   
    element:                                                          
                                                                      
    void foo(TMatrixF &m,TMatrixF &m1)                                
    {                                                                 
      typedef  double (*dfunc_t)(double);                             
      class ApplyFunction : public TElementActionF {                  
        dfunc_t fFunc;                                                
        void Operation(Float_t &element) { element=fFunc(element); }  
      public:                                                         
        ApplyFunction(dfunc_t func):fFunc(func) {}                    
      };                                                              
      ApplyFunction x(TMath::Sin);                                    
      m.Apply(x);                                                     
    }                                                                 
                                                                      
    Validation code $ROOTSYS/test/vmatrix.cxx and vvector.cxx contain 
    a few more examples of that kind.                                 
                                                                      
 6. Lazy matrices: instead of returning an object return a "recipe"   
    how to make it. The full matrix would be rolled out only when     
    and where it's needed:                                            
       TMatrixF haar = THaarMatrixF(5);                               
    THaarMatrixF() is a *class*, not a simple function. However       
    similar this looks to a returning of an object (see note #1       
    above), it's dramatically different. THaarMatrixF() constructs a  
    TMatLazy, an object of just a few bytes long. A
    "TMatrixF(const TMatLazy &recipe)" constructor follows the        
    recipe and makes the matrix haar() right in place. No matrix      
    element is moved whatsoever!                                      
                                                                      


TMatrixFBase& SetMatrixArray(const Float_t *data,Option_t *option)
 Copy array data to matrix . It is assumed that array is of size >= fNelems
 (=)))) fNrows*fNcols
 option indicates how the data is stored in the array:
 option =
          'F'   : column major (Fortran) m[i][j] = array[i+j*fNrows]
          else  : row major    (C)       m[i][j] = array[i*fNcols+j] (default)

Bool_t IsSymmetric() const

void GetMatrix2Array(Float_t *data,Option_t *option) const
 Copy matrix data to array . It is assumed that array is of size >= fNelems
 (=)))) fNrows*fNcols
 option indicates how the data is stored in the array:
 option =
          'F'   : column major (Fortran) array[i+j*fNrows] = m[i][j]
          else  : row major    (C)       array[i*fNcols+j] = m[i][j] (default)

TMatrixFBase& InsertRow(Int_t rown,Int_t coln,const Float_t *v,Int_t n)

void ExtractRow(Int_t rown,Int_t coln,Float_t *v,Int_t n) const
 Store in array v, n matrix elements of row rown starting at column coln

TMatrixFBase& Shift(Int_t row_shift,Int_t col_shift)
 Shift the row index by adding row_shift and the column index by adding
 col_shift, respectively. So [rowLwb..rowUpb][colLwb..colUpb] becomes
 [rowLwb+row_shift..rowUpb+row_shift][colLwb+col_shift..colUpb+col_shift]

TMatrixFBase& Zero()

TMatrixFBase& Abs()
 Take an absolute value of a matrix, i.e. apply Abs() to each element.

TMatrixFBase& Sqr()
 Square each element of the matrix.

TMatrixFBase& Sqrt()
 Take square root of all elements.

TMatrixFBase& UnitMatrix()
 Make a unit matrix (matrix need not be a square one).

TMatrixFBase& NormByDiag(const TVectorF &v,Option_t *option)
 option:
 "D"   :  b(i,j) = a(i,j)/sqrt(abs*(v(i)*v(j)))  (default)
 else  :  b(i,j) = a(i,j)*sqrt(abs*(v(i)*v(j)))  (default)

Float_t RowNorm() const
 Row matrix norm, MAX{ SUM{ |M(i,j)|, over j}, over i}.
 The norm is induced by the infinity vector norm.

Float_t ColNorm() const
 Column matrix norm, MAX{ SUM{ |M(i,j)|, over i}, over j}.
 The norm is induced by the 1 vector norm.

Float_t E2Norm() const
 Square of the Euclidian norm, SUM{ m(i,j)^2 }.

Int_t NonZeros() const
 Compute the number of elements != 0.0

Float_t Sum() const
 Compute sum of elements

Float_t Min() const
 return minimum matrix element value

Float_t Max() const
 return maximum vector element value

void Draw(Option_t *option)
 Draw this matrix using an intermediate histogram
 The histogram is named "TMatrixF" by default and no title

void Print(Option_t *name) const
 Print the matrix as a table of elements .

TMatrixFBase& Apply(const TElementActionF &action)

TMatrixFBase& Apply(const TElementPosActionF &action)
 Apply action to each element of the matrix. To action the location
 of the current element is passed.

TMatrixFBase& Randomize(Float_t alpha,Float_t beta,Double_t &seed)
 randomize matrix element values

void Streamer(TBuffer &R__b)
 Stream an object of class TMatrixFBase.



Inline Functions


                  void ~TMatrixFBase()
              Float_t* GetElements()
                  void Allocate(Int_t nrows, Int_t ncols, Int_t row_lwb = 0, Int_t col_lwb = 0, Int_t init = 0, Int_t nr_nonzero = -1)
                 Int_t GetRowLwb() const
                 Int_t GetRowUpb() const
                 Int_t GetNrows() const
                 Int_t GetColLwb() const
                 Int_t GetColUpb() const
                 Int_t GetNcols() const
                 Int_t GetNoElements() const
               Float_t GetTol() const
        const Float_t* GetMatrixArray() const
              Float_t* GetMatrixArray()
          const Int_t* GetRowIndexArray() const
                Int_t* GetRowIndexArray()
          const Int_t* GetColIndexArray() const
                Int_t* GetColIndexArray()
         TMatrixFBase& SetRowIndexArray(Int_t* data)
         TMatrixFBase& SetColIndexArray(Int_t* data)
               Float_t SetTol(Float_t tol)
                  void Clear(Option_t* option = "")
                  void Invalidate()
                  void MakeValid()
                Bool_t IsValid() const
                Bool_t IsOwner() const
         TMatrixFBase& GetSub(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, TMatrixFBase& target, Option_t* option = "S") const
         TMatrixFBase& SetSub(Int_t row_lwb, Int_t col_lwb, const TMatrixFBase& source)
         TMatrixFBase& ResizeTo(Int_t nrows, Int_t ncols, Int_t nr_nonzeros = -1)
         TMatrixFBase& ResizeTo(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb, Int_t nr_nonzeros = -1)
              Double_t Determinant() const
                  void Determinant(Double_t& d1, Double_t& d2) const
               Float_t NormInf() const
               Float_t Norm1() const
               Float_t operator()(Int_t rown, Int_t coln) const
              Float_t& operator()(Int_t rown, Int_t coln)
                Bool_t operator==(Float_t val) const
                Bool_t operator!=(Float_t val) const
                Bool_t operator<(Float_t val) const
                Bool_t operator<=(Float_t val) const
                Bool_t operator>(Float_t val) const
                Bool_t operator>=(Float_t val) const
               TClass* Class()
               TClass* IsA() const
                  void ShowMembers(TMemberInspector& insp, char* parent)
                  void StreamerNVirtual(TBuffer& b)
         TMatrixFBase& operator=(const TMatrixFBase&)


Last update: root/matrix:$Name: $:$Id: TMatrixFBase.cxx,v 1.17 2004/10/23 20:19:05 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.