invertMatrix.C

Go to the documentation of this file.
00001 //Author: Eddy Offermann
00002 
00003 // This macro shows several ways to invert a matrix . Each  method
00004 // is a trade-off between accuracy of the inversion and speed.
00005 // Which method to chose depends on "how well-behaved" the matrix is.
00006 // This is best checked through a call to Condition(), available in each
00007 // decomposition class. A second possibilty (less preferred) would be to
00008 // check the determinant
00009 //
00010 //   USAGE
00011 //   -----
00012 // This macro can be execued via CINT or via ACLIC
00013 // - via CINT, do
00014 //    root > .x invertMatrix.C
00015 // - via ACLIC
00016 //    root > gSystem->Load("libMatrix");
00017 //    root > .x invertMatrix.C+
00018 
00019 #ifndef __CINT__
00020 #include "Riostream.h"
00021 #include "TMath.h"
00022 #include "TMatrixD.h"
00023 #include "TMatrixDLazy.h"
00024 #include "TVectorD.h"
00025 #include "TDecompLU.h"
00026 #include "TDecompSVD.h"
00027 #endif
00028 
00029 void invertMatrix(Int_t msize=6)
00030 {
00031 #ifdef __CINT__
00032   gSystem->Load("libMatrix");
00033 #endif
00034   if (msize < 2 || msize > 10) {
00035     cout << "2 <= msize <= 10" <<endl;
00036     return;
00037   }
00038   cout << "--------------------------------------------------------" <<endl;
00039   cout << "Inversion results for a ("<<msize<<","<<msize<<") matrix" <<endl;
00040   cout << "For each inversion procedure we check the maxmimum size " <<endl;
00041   cout << "of the off-diagonal elements of Inv(A) * A              " <<endl;
00042   cout << "--------------------------------------------------------" <<endl;
00043 
00044   TMatrixD H_square = THilbertMatrixD(msize,msize);
00045 
00046 //  1. InvertFast(Double_t *det=0)
00047 //   It is identical to Invert() for sizes > 6 x 6 but for smaller sizes, the
00048 //   inversion is performed according to Cramer's rule by explicitly calculating
00049 //   all Jacobi's sub-determinants . For instance for a 6 x 6 matrix this means:
00050 //    # of 5 x 5 determinant : 36 
00051 //    # of 4 x 4 determinant : 75 
00052 //    # of 3 x 3 determinant : 80 
00053 //    # of 2 x 2 determinant : 45    (see TMatrixD/FCramerInv.cxx)
00054 //
00055 //    The only "quality" control in this process is to check whether the 6 x 6
00056 //    determinant is unequal 0 . But speed gains are significant compared to Invert() ,
00057 //    upto an order of magnitude for sizes <= 4 x 4
00058 //
00059 //    The inversion is done "in place", so the original matrix will be overwritten
00060 //    If a pointer to a Double_t is supplied the determinant is calculated
00061 //
00062 
00063   cout << "1. Use .InvertFast(&det)" <<endl;
00064   if (msize > 6)
00065     cout << " for ("<<msize<<","<<msize<<") this is identical to .Invert(&det)" <<endl;
00066 
00067   Double_t det1;
00068   TMatrixD H1 = H_square;
00069   H1.InvertFast(&det1);
00070 
00071   // Get the maximum off-diagonal matrix value . One way to do this is to set the
00072   // diagonal to zero .
00073 
00074   TMatrixD U1(H1,TMatrixD::kMult,H_square);
00075   TMatrixDDiag diag1(U1); diag1 = 0.0;
00076   const Double_t U1_max_offdiag = (U1.Abs()).Max();
00077   cout << "  Maximum off-diagonal = " << U1_max_offdiag << endl;
00078   cout << "  Determinant          = " << det1 <<endl;
00079 
00080 // 2. Invert(Double_t *det=0)
00081 //   Again the inversion is performed in place .
00082 //   It consists out of a sequence of calls to the decomposition classes . For instance
00083 //   for the general dense matrix TMatrixD the LU decomposition is invoked:
00084 //    - The matrix is decomposed using a scheme according to Crout which involves
00085 //      "implicit partial pivoting", see for instance Num. Recip. (we have also available
00086 //      a decomposition scheme that does not the scaling and is therefore even slightly
00087 //      faster but less stable)
00088 //      With each decomposition, a tolerance has to be specified . If this tolerance
00089 //      requirement is not met, the matrix is regarded as being singular. The value
00090 //      passed to this decomposition, is the data member fTol of the matrix . Its
00091 //      default value is DBL_EPSILON, which is defined as the smallest nuber so that
00092 //      1+DBL_EPSILON > 1
00093 //    - The last step is a standard forward/backward substitution .
00094 //
00095 //   It is important to realize that both InvertFast() and Invert() are "one-shot" deals , speed
00096 //   comes at a price . If something goes wrong because the matrix is (near) singular, you have
00097 //   overwritten your original matrix and  no factorization is available anymore to get more
00098 //   information like condition number or change the tolerance number .
00099 //
00100 //   All other calls in the matrix classes involving inversion like the ones with the "smart"
00101 //   constructors (kInverted,kInvMult...) use this inversion method .
00102 //
00103 
00104   cout << "2. Use .Invert(&det)" <<endl;
00105 
00106   Double_t det2;
00107   TMatrixD H2 = H_square;
00108   H2.Invert(&det2);
00109 
00110   TMatrixD U2(H2,TMatrixD::kMult,H_square);
00111   TMatrixDDiag diag2(U2); diag2 = 0.0;
00112   const Double_t U2_max_offdiag = (U2.Abs()).Max();
00113   cout << "  Maximum off-diagonal = " << U2_max_offdiag << endl;
00114   cout << "  Determinant          = " << det2 <<endl;
00115 
00116 // 3. Inversion through LU decomposition
00117 //   The (default) algorithms used are similar to 2. (Not identical because in 2, the whole
00118 //   calculation is done "in-place". Here the orginal matrix is copied (so more memory
00119 //   management => slower) and several operations can be performed without having to repeat
00120 //   the decomposition step .
00121 //   Inverting a matrix is nothing else than solving a set of equations where the rhs is given
00122 //   by the unit matrix, so the steps to take are identical to those solving a linear equation :
00123 //
00124 
00125   cout << "3. Use TDecompLU" <<endl;
00126 
00127   TMatrixD H3 = H_square;
00128   TDecompLU lu(H_square);
00129 
00130   // Any operation that requires a decomposition will trigger it . The class keeps
00131   // an internal state so that following operations will not perform the decomposition again
00132   // unless the matrix is changed through SetMatrix(..)
00133   // One might want to proceed more cautiously by invoking first Decompose() and check its
00134   // return value before proceeding....
00135 
00136   lu.Invert(H3);
00137   Double_t d1_lu; Double_t d2_lu;
00138   lu.Det(d1_lu,d2_lu);
00139   Double_t det3 = d1_lu*TMath::Power(2.,d2_lu);
00140 
00141   TMatrixD U3(H3,TMatrixD::kMult,H_square);
00142   TMatrixDDiag diag3(U3); diag3 = 0.0;
00143   const Double_t U3_max_offdiag = (U3.Abs()).Max();
00144   cout << "  Maximum off-diagonal = " << U3_max_offdiag << endl;
00145   cout << "  Determinant          = " << det3 <<endl;
00146 
00147 // 4. Inversion through SVD decomposition
00148 //   For SVD and QRH, the (n x m) matrix does only have to fulfill n >=m . In case n > m
00149 //   a pseudo-inverse is calculated
00150   cout << "4. Use TDecompSVD on non-square matrix" <<endl;
00151 
00152   TMatrixD H_nsquare = THilbertMatrixD(msize,msize-1);
00153 
00154   TDecompSVD svd(H_nsquare);
00155 
00156   TMatrixD H4 = svd.Invert();
00157   Double_t d1_svd; Double_t d2_svd;
00158   svd.Det(d1_svd,d2_svd);
00159   Double_t det4 = d1_svd*TMath::Power(2.,d2_svd);
00160 
00161   TMatrixD U4(H4,TMatrixD::kMult,H_nsquare);
00162   TMatrixDDiag diag4(U4); diag4 = 0.0;
00163   const Double_t U4_max_offdiag = (U4.Abs()).Max();
00164   cout << "  Maximum off-diagonal = " << U4_max_offdiag << endl;
00165   cout << "  Determinant          = " << det4 <<endl;
00166 }

Generated on Tue Jul 5 15:44:51 2011 for ROOT_528-00b_version by  doxygen 1.5.1