ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
invertMatrix.C
Go to the documentation of this file.
1 //Author: Eddy Offermann
2 
3 // This macro shows several ways to invert a matrix . Each method
4 // is a trade-off between accuracy of the inversion and speed.
5 // Which method to chose depends on "how well-behaved" the matrix is.
6 // This is best checked through a call to Condition(), available in each
7 // decomposition class. A second possibilty (less preferred) would be to
8 // check the determinant
9 //
10 // USAGE
11 // -----
12 // This macro can be execued via CINT or via ACLIC
13 // - via CINT, do
14 // root > .x invertMatrix.C
15 // - via ACLIC
16 // root > gSystem->Load("libMatrix");
17 // root > .x invertMatrix.C+
18 
19 #ifndef __CINT__
20 #include "Riostream.h"
21 #include "TMath.h"
22 #include "TMatrixD.h"
23 #include "TMatrixDLazy.h"
24 #include "TVectorD.h"
25 #include "TDecompLU.h"
26 #include "TDecompSVD.h"
27 #endif
28 
29 void invertMatrix(Int_t msize=6)
30 {
31 #ifdef __CINT__
32  gSystem->Load("libMatrix");
33 #endif
34  if (msize < 2 || msize > 10) {
35  cout << "2 <= msize <= 10" <<endl;
36  return;
37  }
38  cout << "--------------------------------------------------------" <<endl;
39  cout << "Inversion results for a ("<<msize<<","<<msize<<") matrix" <<endl;
40  cout << "For each inversion procedure we check the maxmimum size " <<endl;
41  cout << "of the off-diagonal elements of Inv(A) * A " <<endl;
42  cout << "--------------------------------------------------------" <<endl;
43 
44  TMatrixD H_square = THilbertMatrixD(msize,msize);
45 
46 // 1. InvertFast(Double_t *det=0)
47 // It is identical to Invert() for sizes > 6 x 6 but for smaller sizes, the
48 // inversion is performed according to Cramer's rule by explicitly calculating
49 // all Jacobi's sub-determinants . For instance for a 6 x 6 matrix this means:
50 // # of 5 x 5 determinant : 36
51 // # of 4 x 4 determinant : 75
52 // # of 3 x 3 determinant : 80
53 // # of 2 x 2 determinant : 45 (see TMatrixD/FCramerInv.cxx)
54 //
55 // The only "quality" control in this process is to check whether the 6 x 6
56 // determinant is unequal 0 . But speed gains are significant compared to Invert() ,
57 // upto an order of magnitude for sizes <= 4 x 4
58 //
59 // The inversion is done "in place", so the original matrix will be overwritten
60 // If a pointer to a Double_t is supplied the determinant is calculated
61 //
62 
63  cout << "1. Use .InvertFast(&det)" <<endl;
64  if (msize > 6)
65  cout << " for ("<<msize<<","<<msize<<") this is identical to .Invert(&det)" <<endl;
66 
67  Double_t det1;
68  TMatrixD H1 = H_square;
69  H1.InvertFast(&det1);
70 
71  // Get the maximum off-diagonal matrix value . One way to do this is to set the
72  // diagonal to zero .
73 
74  TMatrixD U1(H1,TMatrixD::kMult,H_square);
75  TMatrixDDiag diag1(U1); diag1 = 0.0;
76  const Double_t U1_max_offdiag = (U1.Abs()).Max();
77  cout << " Maximum off-diagonal = " << U1_max_offdiag << endl;
78  cout << " Determinant = " << det1 <<endl;
79 
80 // 2. Invert(Double_t *det=0)
81 // Again the inversion is performed in place .
82 // It consists out of a sequence of calls to the decomposition classes . For instance
83 // for the general dense matrix TMatrixD the LU decomposition is invoked:
84 // - The matrix is decomposed using a scheme according to Crout which involves
85 // "implicit partial pivoting", see for instance Num. Recip. (we have also available
86 // a decomposition scheme that does not the scaling and is therefore even slightly
87 // faster but less stable)
88 // With each decomposition, a tolerance has to be specified . If this tolerance
89 // requirement is not met, the matrix is regarded as being singular. The value
90 // passed to this decomposition, is the data member fTol of the matrix . Its
91 // default value is DBL_EPSILON, which is defined as the smallest nuber so that
92 // 1+DBL_EPSILON > 1
93 // - The last step is a standard forward/backward substitution .
94 //
95 // It is important to realize that both InvertFast() and Invert() are "one-shot" deals , speed
96 // comes at a price . If something goes wrong because the matrix is (near) singular, you have
97 // overwritten your original matrix and no factorization is available anymore to get more
98 // information like condition number or change the tolerance number .
99 //
100 // All other calls in the matrix classes involving inversion like the ones with the "smart"
101 // constructors (kInverted,kInvMult...) use this inversion method .
102 //
103 
104  cout << "2. Use .Invert(&det)" <<endl;
105 
106  Double_t det2;
107  TMatrixD H2 = H_square;
108  H2.Invert(&det2);
109 
110  TMatrixD U2(H2,TMatrixD::kMult,H_square);
111  TMatrixDDiag diag2(U2); diag2 = 0.0;
112  const Double_t U2_max_offdiag = (U2.Abs()).Max();
113  cout << " Maximum off-diagonal = " << U2_max_offdiag << endl;
114  cout << " Determinant = " << det2 <<endl;
115 
116 // 3. Inversion through LU decomposition
117 // The (default) algorithms used are similar to 2. (Not identical because in 2, the whole
118 // calculation is done "in-place". Here the orginal matrix is copied (so more memory
119 // management => slower) and several operations can be performed without having to repeat
120 // the decomposition step .
121 // Inverting a matrix is nothing else than solving a set of equations where the rhs is given
122 // by the unit matrix, so the steps to take are identical to those solving a linear equation :
123 //
124 
125  cout << "3. Use TDecompLU" <<endl;
126 
127  TMatrixD H3 = H_square;
128  TDecompLU lu(H_square);
129 
130  // Any operation that requires a decomposition will trigger it . The class keeps
131  // an internal state so that following operations will not perform the decomposition again
132  // unless the matrix is changed through SetMatrix(..)
133  // One might want to proceed more cautiously by invoking first Decompose() and check its
134  // return value before proceeding....
135 
136  lu.Invert(H3);
137  Double_t d1_lu; Double_t d2_lu;
138  lu.Det(d1_lu,d2_lu);
139  Double_t det3 = d1_lu*TMath::Power(2.,d2_lu);
140 
141  TMatrixD U3(H3,TMatrixD::kMult,H_square);
142  TMatrixDDiag diag3(U3); diag3 = 0.0;
143  const Double_t U3_max_offdiag = (U3.Abs()).Max();
144  cout << " Maximum off-diagonal = " << U3_max_offdiag << endl;
145  cout << " Determinant = " << det3 <<endl;
146 
147 // 4. Inversion through SVD decomposition
148 // For SVD and QRH, the (n x m) matrix does only have to fulfill n >=m . In case n > m
149 // a pseudo-inverse is calculated
150  cout << "4. Use TDecompSVD on non-square matrix" <<endl;
151 
152  TMatrixD H_nsquare = THilbertMatrixD(msize,msize-1);
153 
154  TDecompSVD svd(H_nsquare);
155 
156  TMatrixD H4 = svd.Invert();
157  Double_t d1_svd; Double_t d2_svd;
158  svd.Det(d1_svd,d2_svd);
159  Double_t det4 = d1_svd*TMath::Power(2.,d2_svd);
160 
161  TMatrixD U4(H4,TMatrixD::kMult,H_nsquare);
162  TMatrixDDiag diag4(U4); diag4 = 0.0;
163  const Double_t U4_max_offdiag = (U4.Abs()).Max();
164  cout << " Maximum off-diagonal = " << U4_max_offdiag << endl;
165  cout << " Determinant = " << det4 <<endl;
166 }
virtual void Det(Double_t &d1, Double_t &d2)
Matrix determinant det = d1*TMath::Power(2.,d2)
Definition: TDecompSVD.cxx:840
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1766
int Int_t
Definition: RtypesCore.h:41
virtual void Det(Double_t &d1, Double_t &d2)
Calculate determinant det = d1*TMath::Power(2.,d2)
Definition: TDecompLU.cxx:509
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
Definition: TMatrixT.cxx:1388
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
Bool_t Invert(TMatrixD &inv)
For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit The user should alw...
Definition: TDecompSVD.cxx:875
THilbertMatrixT< Double_t > THilbertMatrixD
Definition: TMatrixDLazy.h:41
double Double_t
Definition: RtypesCore.h:55
Bool_t Invert(TMatrixD &inv)
For a matrix A(m,m), its inverse A_inv is defined as A * A_inv = A_inv * A = unit (m x m) Ainv is ret...
Definition: TDecompLU.cxx:526
void invertMatrix(Int_t msize=6)
Definition: invertMatrix.C:29
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
TMatrixT< Element > & InvertFast(Double_t *det=0)
Invert the matrix and calculate its determinant, however upto (6x6) a fast Cramer inversion is used ...
Definition: TMatrixT.cxx:1402
virtual TMatrixTBase< Element > & Abs()
Take an absolute value of a matrix, i.e. apply Abs() to each element.