ROOT logo
//Author: Eddy Offermann

// This macro shows several ways to invert a matrix . Each  method
// is a trade-off between accuracy of the inversion and speed.
// Which method to chose depends on "how well-behaved" the matrix is.
// This is best checked through a call to Condition(), available in each
// decomposition class. A second possibilty (less preferred) would be to
// check the determinant
//
//   USAGE
//   -----
// This macro can be execued via CINT or via ACLIC
// - via CINT, do
//    root > .x invertMatrix.C
// - via ACLIC
//    root > gSystem->Load("libMatrix");
//    root > .x invertMatrix.C+

#ifndef __CINT__
#include "Riostream.h"
#include "TMath.h"
#include "TMatrixD.h"
#include "TMatrixDLazy.h"
#include "TVectorD.h"
#include "TDecompLU.h"
#include "TDecompSVD.h"
#endif

void invertMatrix(Int_t msize=6)
{
#ifdef __CINT__
  gSystem->Load("libMatrix");
#endif
  if (msize < 2 || msize > 10) {
    cout << "2 <= msize <= 10" <<endl;
    return;
  }
  cout << "--------------------------------------------------------" <<endl;
  cout << "Inversion results for a ("<<msize<<","<<msize<<") matrix" <<endl;
  cout << "For each inversion procedure we check the maxmimum size " <<endl;
  cout << "of the off-diagonal elements of Inv(A) * A              " <<endl;
  cout << "--------------------------------------------------------" <<endl;

  TMatrixD H_square = THilbertMatrixD(msize,msize);

//  1. InvertFast(Double_t *det=0)
//   It is identical to Invert() for sizes > 6 x 6 but for smaller sizes, the
//   inversion is performed according to Cramer's rule by explicitly calculating
//   all Jacobi's sub-determinants . For instance for a 6 x 6 matrix this means:
//    # of 5 x 5 determinant : 36 
//    # of 4 x 4 determinant : 75 
//    # of 3 x 3 determinant : 80 
//    # of 2 x 2 determinant : 45    (see TMatrixD/FCramerInv.cxx)
//
//    The only "quality" control in this process is to check whether the 6 x 6
//    determinant is unequal 0 . But speed gains are significant compared to Invert() ,
//    upto an order of magnitude for sizes <= 4 x 4
//
//    The inversion is done "in place", so the original matrix will be overwritten
//    If a pointer to a Double_t is supplied the determinant is calculated
//

  cout << "1. Use .InvertFast(&det)" <<endl;
  if (msize > 6)
    cout << " for ("<<msize<<","<<msize<<") this is identical to .Invert(&det)" <<endl;

  Double_t det1;
  TMatrixD H1 = H_square;
  H1.InvertFast(&det1);

  // Get the maximum off-diagonal matrix value . One way to do this is to set the
  // diagonal to zero .

  TMatrixD U1(H1,TMatrixD::kMult,H_square);
  TMatrixDDiag diag1(U1); diag1 = 0.0;
  const Double_t U1_max_offdiag = (U1.Abs()).Max();
  cout << "  Maximum off-diagonal = " << U1_max_offdiag << endl;
  cout << "  Determinant          = " << det1 <<endl;

// 2. Invert(Double_t *det=0)
//   Again the inversion is performed in place .
//   It consists out of a sequence of calls to the decomposition classes . For instance
//   for the general dense matrix TMatrixD the LU decomposition is invoked:
//    - The matrix is decomposed using a scheme according to Crout which involves
//      "implicit partial pivoting", see for instance Num. Recip. (we have also available
//      a decomposition scheme that does not the scaling and is therefore even slightly
//      faster but less stable)
//      With each decomposition, a tolerance has to be specified . If this tolerance
//      requirement is not met, the matrix is regarded as being singular. The value
//      passed to this decomposition, is the data member fTol of the matrix . Its
//      default value is DBL_EPSILON, which is defined as the smallest nuber so that
//      1+DBL_EPSILON > 1
//    - The last step is a standard forward/backward substitution .
//
//   It is important to realize that both InvertFast() and Invert() are "one-shot" deals , speed
//   comes at a price . If something goes wrong because the matrix is (near) singular, you have
//   overwritten your original matrix and  no factorization is available anymore to get more
//   information like condition number or change the tolerance number .
//
//   All other calls in the matrix classes involving inversion like the ones with the "smart"
//   constructors (kInverted,kInvMult...) use this inversion method .
//

  cout << "2. Use .Invert(&det)" <<endl;

  Double_t det2;
  TMatrixD H2 = H_square;
  H2.Invert(&det2);

  TMatrixD U2(H2,TMatrixD::kMult,H_square);
  TMatrixDDiag diag2(U2); diag2 = 0.0;
  const Double_t U2_max_offdiag = (U2.Abs()).Max();
  cout << "  Maximum off-diagonal = " << U2_max_offdiag << endl;
  cout << "  Determinant          = " << det2 <<endl;

// 3. Inversion through LU decomposition
//   The (default) algorithms used are similar to 2. (Not identical because in 2, the whole
//   calculation is done "in-place". Here the orginal matrix is copied (so more memory
//   management => slower) and several operations can be performed without having to repeat
//   the decomposition step .
//   Inverting a matrix is nothing else than solving a set of equations where the rhs is given
//   by the unit matrix, so the steps to take are identical to those solving a linear equation :
//

  cout << "3. Use TDecompLU" <<endl;

  TMatrixD H3 = H_square;
  TDecompLU lu(H_square);

  // Any operation that requires a decomposition will trigger it . The class keeps
  // an internal state so that following operations will not perform the decomposition again
  // unless the matrix is changed through SetMatrix(..)
  // One might want to proceed more cautiously by invoking first Decompose() and check its
  // return value before proceeding....

  lu.Invert(H3);
  Double_t d1_lu; Double_t d2_lu;
  lu.Det(d1_lu,d2_lu);
  Double_t det3 = d1_lu*TMath::Power(2.,d2_lu);

  TMatrixD U3(H3,TMatrixD::kMult,H_square);
  TMatrixDDiag diag3(U3); diag3 = 0.0;
  const Double_t U3_max_offdiag = (U3.Abs()).Max();
  cout << "  Maximum off-diagonal = " << U3_max_offdiag << endl;
  cout << "  Determinant          = " << det3 <<endl;

// 4. Inversion through SVD decomposition
//   For SVD and QRH, the (n x m) matrix does only have to fulfill n >=m . In case n > m
//   a pseudo-inverse is calculated
  cout << "4. Use TDecompSVD on non-square matrix" <<endl;

  TMatrixD H_nsquare = THilbertMatrixD(msize,msize-1);

  TDecompSVD svd(H_nsquare);

  TMatrixD H4 = svd.Invert();
  Double_t d1_svd; Double_t d2_svd;
  svd.Det(d1_svd,d2_svd);
  Double_t det4 = d1_svd*TMath::Power(2.,d2_svd);

  TMatrixD U4(H4,TMatrixD::kMult,H_nsquare);
  TMatrixDDiag diag4(U4); diag4 = 0.0;
  const Double_t U4_max_offdiag = (U4.Abs()).Max();
  cout << "  Maximum off-diagonal = " << U4_max_offdiag << endl;
  cout << "  Determinant          = " << det4 <<endl;
}
 invertMatrix.C:1
 invertMatrix.C:2
 invertMatrix.C:3
 invertMatrix.C:4
 invertMatrix.C:5
 invertMatrix.C:6
 invertMatrix.C:7
 invertMatrix.C:8
 invertMatrix.C:9
 invertMatrix.C:10
 invertMatrix.C:11
 invertMatrix.C:12
 invertMatrix.C:13
 invertMatrix.C:14
 invertMatrix.C:15
 invertMatrix.C:16
 invertMatrix.C:17
 invertMatrix.C:18
 invertMatrix.C:19
 invertMatrix.C:20
 invertMatrix.C:21
 invertMatrix.C:22
 invertMatrix.C:23
 invertMatrix.C:24
 invertMatrix.C:25
 invertMatrix.C:26
 invertMatrix.C:27
 invertMatrix.C:28
 invertMatrix.C:29
 invertMatrix.C:30
 invertMatrix.C:31
 invertMatrix.C:32
 invertMatrix.C:33
 invertMatrix.C:34
 invertMatrix.C:35
 invertMatrix.C:36
 invertMatrix.C:37
 invertMatrix.C:38
 invertMatrix.C:39
 invertMatrix.C:40
 invertMatrix.C:41
 invertMatrix.C:42
 invertMatrix.C:43
 invertMatrix.C:44
 invertMatrix.C:45
 invertMatrix.C:46
 invertMatrix.C:47
 invertMatrix.C:48
 invertMatrix.C:49
 invertMatrix.C:50
 invertMatrix.C:51
 invertMatrix.C:52
 invertMatrix.C:53
 invertMatrix.C:54
 invertMatrix.C:55
 invertMatrix.C:56
 invertMatrix.C:57
 invertMatrix.C:58
 invertMatrix.C:59
 invertMatrix.C:60
 invertMatrix.C:61
 invertMatrix.C:62
 invertMatrix.C:63
 invertMatrix.C:64
 invertMatrix.C:65
 invertMatrix.C:66
 invertMatrix.C:67
 invertMatrix.C:68
 invertMatrix.C:69
 invertMatrix.C:70
 invertMatrix.C:71
 invertMatrix.C:72
 invertMatrix.C:73
 invertMatrix.C:74
 invertMatrix.C:75
 invertMatrix.C:76
 invertMatrix.C:77
 invertMatrix.C:78
 invertMatrix.C:79
 invertMatrix.C:80
 invertMatrix.C:81
 invertMatrix.C:82
 invertMatrix.C:83
 invertMatrix.C:84
 invertMatrix.C:85
 invertMatrix.C:86
 invertMatrix.C:87
 invertMatrix.C:88
 invertMatrix.C:89
 invertMatrix.C:90
 invertMatrix.C:91
 invertMatrix.C:92
 invertMatrix.C:93
 invertMatrix.C:94
 invertMatrix.C:95
 invertMatrix.C:96
 invertMatrix.C:97
 invertMatrix.C:98
 invertMatrix.C:99
 invertMatrix.C:100
 invertMatrix.C:101
 invertMatrix.C:102
 invertMatrix.C:103
 invertMatrix.C:104
 invertMatrix.C:105
 invertMatrix.C:106
 invertMatrix.C:107
 invertMatrix.C:108
 invertMatrix.C:109
 invertMatrix.C:110
 invertMatrix.C:111
 invertMatrix.C:112
 invertMatrix.C:113
 invertMatrix.C:114
 invertMatrix.C:115
 invertMatrix.C:116
 invertMatrix.C:117
 invertMatrix.C:118
 invertMatrix.C:119
 invertMatrix.C:120
 invertMatrix.C:121
 invertMatrix.C:122
 invertMatrix.C:123
 invertMatrix.C:124
 invertMatrix.C:125
 invertMatrix.C:126
 invertMatrix.C:127
 invertMatrix.C:128
 invertMatrix.C:129
 invertMatrix.C:130
 invertMatrix.C:131
 invertMatrix.C:132
 invertMatrix.C:133
 invertMatrix.C:134
 invertMatrix.C:135
 invertMatrix.C:136
 invertMatrix.C:137
 invertMatrix.C:138
 invertMatrix.C:139
 invertMatrix.C:140
 invertMatrix.C:141
 invertMatrix.C:142
 invertMatrix.C:143
 invertMatrix.C:144
 invertMatrix.C:145
 invertMatrix.C:146
 invertMatrix.C:147
 invertMatrix.C:148
 invertMatrix.C:149
 invertMatrix.C:150
 invertMatrix.C:151
 invertMatrix.C:152
 invertMatrix.C:153
 invertMatrix.C:154
 invertMatrix.C:155
 invertMatrix.C:156
 invertMatrix.C:157
 invertMatrix.C:158
 invertMatrix.C:159
 invertMatrix.C:160
 invertMatrix.C:161
 invertMatrix.C:162
 invertMatrix.C:163
 invertMatrix.C:164
 invertMatrix.C:165
 invertMatrix.C:166
 invertMatrix.C:167