ROOT logo
// @(#)root/matrix:$Id: TDecompLU.h 20882 2007-11-19 11:31:26Z rdm $
// Authors: Fons Rademakers, Eddy Offermann   Dec 2003

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

#ifndef ROOT_TDecompLU
#define ROOT_TDecompLU

///////////////////////////////////////////////////////////////////////////
//                                                                       //
// LU Decomposition class                                                //
//                                                                       //
///////////////////////////////////////////////////////////////////////////

#ifndef ROOT_TDecompBase
#include "TDecompBase.h"
#endif

class TDecompLU : public TDecompBase
{
protected :

   Int_t     fImplicitPivot; // control to determine implicit row scale before
                             //  deciding on the pivot (Crout method)
   Int_t     fNIndex;        // size of row permutation index
   Int_t    *fIndex;         //[fNIndex] row permutation index
   Double_t  fSign;          // = +/- 1 reflecting even/odd row permutations, resp.
   TMatrixD  fLU;            // decomposed matrix so that a = l u where
                             // l is stored lower left and u upper right side

   static Bool_t DecomposeLUCrout(TMatrixD &lu,Int_t *index,Double_t &sign,Double_t tol,Int_t &nrZeros);
   static Bool_t DecomposeLUGauss(TMatrixD &lu,Int_t *index,Double_t &sign,Double_t tol,Int_t &nrZeros);

   virtual const TMatrixDBase &GetDecompMatrix() const { return fLU; }

public :

   TDecompLU();
   explicit TDecompLU(Int_t nrows);
   TDecompLU(Int_t row_lwb,Int_t row_upb);
   TDecompLU(const TMatrixD &m,Double_t tol = 0.0,Int_t implicit = 1);
   TDecompLU(const TDecompLU &another);
   virtual ~TDecompLU() {if (fIndex) delete [] fIndex; fIndex = 0; }

           const TMatrixD  GetMatrix ();
   virtual       Int_t     GetNrows  () const { return fLU.GetNrows(); }
   virtual       Int_t     GetNcols  () const { return fLU.GetNcols(); }
           const TMatrixD &GetLU     ()       { if ( !TestBit(kDecomposed) ) Decompose();
                                                return fLU; }

   virtual       void      SetMatrix (const TMatrixD &a);

   virtual Bool_t   Decompose  ();
   virtual Bool_t   Solve      (      TVectorD &b);
   virtual TVectorD Solve      (const TVectorD& b,Bool_t &ok) { TVectorD x = b; ok = Solve(x); return x; }
   virtual Bool_t   Solve      (      TMatrixDColumn &b);
   virtual Bool_t   TransSolve (      TVectorD &b);
   virtual TVectorD TransSolve (const TVectorD& b,Bool_t &ok) { TVectorD x = b; ok = TransSolve(x); return x; }
   virtual Bool_t   TransSolve (      TMatrixDColumn &b);
   virtual void     Det        (Double_t &d1,Double_t &d2);

   static  Bool_t   InvertLU  (TMatrixD &a,Double_t tol,Double_t *det=0);
   Bool_t           Invert    (TMatrixD &inv);
   TMatrixD         Invert    (Bool_t &status);
   TMatrixD         Invert    () { Bool_t status; return Invert(status); }

   void Print(Option_t *opt ="") const; // *MENU*

   TDecompLU &operator= (const TDecompLU &source);

   ClassDef(TDecompLU,1) // Matrix Decompositition LU
};

#endif
 TDecompLU.h:1
 TDecompLU.h:2
 TDecompLU.h:3
 TDecompLU.h:4
 TDecompLU.h:5
 TDecompLU.h:6
 TDecompLU.h:7
 TDecompLU.h:8
 TDecompLU.h:9
 TDecompLU.h:10
 TDecompLU.h:11
 TDecompLU.h:12
 TDecompLU.h:13
 TDecompLU.h:14
 TDecompLU.h:15
 TDecompLU.h:16
 TDecompLU.h:17
 TDecompLU.h:18
 TDecompLU.h:19
 TDecompLU.h:20
 TDecompLU.h:21
 TDecompLU.h:22
 TDecompLU.h:23
 TDecompLU.h:24
 TDecompLU.h:25
 TDecompLU.h:26
 TDecompLU.h:27
 TDecompLU.h:28
 TDecompLU.h:29
 TDecompLU.h:30
 TDecompLU.h:31
 TDecompLU.h:32
 TDecompLU.h:33
 TDecompLU.h:34
 TDecompLU.h:35
 TDecompLU.h:36
 TDecompLU.h:37
 TDecompLU.h:38
 TDecompLU.h:39
 TDecompLU.h:40
 TDecompLU.h:41
 TDecompLU.h:42
 TDecompLU.h:43
 TDecompLU.h:44
 TDecompLU.h:45
 TDecompLU.h:46
 TDecompLU.h:47
 TDecompLU.h:48
 TDecompLU.h:49
 TDecompLU.h:50
 TDecompLU.h:51
 TDecompLU.h:52
 TDecompLU.h:53
 TDecompLU.h:54
 TDecompLU.h:55
 TDecompLU.h:56
 TDecompLU.h:57
 TDecompLU.h:58
 TDecompLU.h:59
 TDecompLU.h:60
 TDecompLU.h:61
 TDecompLU.h:62
 TDecompLU.h:63
 TDecompLU.h:64
 TDecompLU.h:65
 TDecompLU.h:66
 TDecompLU.h:67
 TDecompLU.h:68
 TDecompLU.h:69
 TDecompLU.h:70
 TDecompLU.h:71
 TDecompLU.h:72
 TDecompLU.h:73
 TDecompLU.h:74
 TDecompLU.h:75
 TDecompLU.h:76
 TDecompLU.h:77
 TDecompLU.h:78
 TDecompLU.h:79
 TDecompLU.h:80