Logo ROOT   6.13/01
Reference Guide
TMatrixDEigen.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_TMatrixDEigen
13 #define ROOT_TMatrixDEigen
14 
15 //////////////////////////////////////////////////////////////////////////
16 // //
17 // TMatrixDEigen //
18 // //
19 // Eigenvalues and eigenvectors of a real matrix. //
20 // //
21 //////////////////////////////////////////////////////////////////////////
22 
23 #include "TMatrixD.h"
24 #include "TVectorD.h"
25 
27 {
28 protected :
29 
30  static void MakeHessenBerg (TMatrixD &v,TVectorD &ortho,TMatrixD &H);
31  static void MakeSchurr (TMatrixD &v,TVectorD &d, TVectorD &e,TMatrixD &H);
32  static void Sort (TMatrixD &v,TVectorD &d, TVectorD &e);
33 
34  TMatrixD fEigenVectors; // Eigen-vectors of matrix
35  TVectorD fEigenValuesRe; // Eigen-values
36  TVectorD fEigenValuesIm; // Eigen-values
37 
38 public :
39 
40  enum {kWorkMax = 100}; // size of work array
41 
43  : fEigenVectors(), fEigenValuesRe(), fEigenValuesIm() {};
44  TMatrixDEigen(const TMatrixD &a);
45  TMatrixDEigen(const TMatrixDEigen &another);
46  virtual ~TMatrixDEigen() {}
47 
48 // If matrix A has shape (rowLwb,rowUpb,rowLwb,rowUpb), then each eigen-vector
49 // must have an index running between (rowLwb,rowUpb) .
50 // For convenience, the column index of the eigen-vector matrix
51 // also runs from rowLwb to rowUpb so that the returned matrix
52 // has also index/shape (rowLwb,rowUpb,rowLwb,rowUpb) .
53 // The same is true for the eigen-value vectors an matrix .
54 
55  const TMatrixD &GetEigenVectors () const { return fEigenVectors; }
56  const TVectorD &GetEigenValuesRe() const { return fEigenValuesRe; }
57  const TVectorD &GetEigenValuesIm() const { return fEigenValuesIm; }
58  const TMatrixD GetEigenValues () const;
59 
60  TMatrixDEigen &operator= (const TMatrixDEigen &source);
61 
62  ClassDef(TMatrixDEigen,1) // Eigen-Vectors/Values of a Matrix
63 };
64 #endif
constexpr Double_t H()
Planck's constant in .
Definition: TMath.h:189
TMatrixD fEigenVectors
Definition: TMatrixDEigen.h:34
TMatrixDEigen.
Definition: TMatrixDEigen.h:26
const TVectorD & GetEigenValuesIm() const
Definition: TMatrixDEigen.h:57
static void MakeHessenBerg(TMatrixD &v, TVectorD &ortho, TMatrixD &H)
Nonsymmetric reduction to Hessenberg form.
TVectorD fEigenValuesRe
Definition: TMatrixDEigen.h:35
TVectorD fEigenValuesIm
Definition: TMatrixDEigen.h:36
const TMatrixD GetEigenValues() const
Computes the block diagonal eigenvalue matrix.
SVector< double, 2 > v
Definition: Dict.h:5
static void MakeSchurr(TMatrixD &v, TVectorD &d, TVectorD &e, TMatrixD &H)
Nonsymmetric reduction from Hessenberg to real Schur form.
static void Sort(TMatrixD &v, TVectorD &d, TVectorD &e)
Sort eigenvalues and corresponding vectors in descending order of Re^2+Im^2 of the complex eigenvalue...
const TMatrixD & GetEigenVectors() const
Definition: TMatrixDEigen.h:55
virtual ~TMatrixDEigen()
Definition: TMatrixDEigen.h:46
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
const TVectorD & GetEigenValuesRe() const
Definition: TMatrixDEigen.h:56
TMatrixDEigen & operator=(const TMatrixDEigen &source)
Assignment operator.