solveLinear.C

Go to the documentation of this file.
00001 //Author: Eddy Offermann
00002 // This macro shows several ways to perform a linear least-squares
00003 // analysis . To keep things simple we fit a straight line to 4
00004 // data points
00005 // The first 4 methods use the linear algebra package to find
00006 //  x  such that min (A x - b)^T (A x - b) where A and b
00007 //  are calculated with the data points  and the functional expression :
00008 //
00009 //  1. Normal equations:
00010 //   Expanding the expression (A x - b)^T (A x - b) and taking the
00011 //   derivative wrt x leads to the "Normal Equations":
00012 //   A^T A x = A^T b  where A^T A is a positive definite matrix. Therefore,
00013 //   a Cholesky decomposition scheme can be used to calculate its inverse .
00014 //   This leads to the solution x = (A^T A)^-1 A^T b . All this is done in
00015 //   routine NormalEqn . We made it a bit more complicated by giving the
00016 //   data weights .
00017 //   Numerically this is not the best way to proceed because effctively the
00018 //   condition number of A^T A is twice as large as that of A, making inversion
00019 //   more difficult
00020 //
00021 //  2. SVD
00022 //   One can show that solving A x = b for x with A of size (m x n) and m > n
00023 //   through a Singular Value Decomposition is equivalent to miminizing 
00024 //   (A x - b)^T (A x - b)
00025 //   Numerically , this is the most stable method of all 5
00026 //
00027 //  3. Pseudo Inverse
00028 //   Here we calulate the generalized matrix inverse ("pseudo inverse") by
00029 //   solving A X = Unit for matrix X through an SVD . The formal expression for
00030 //   is X = (A^T A)^-1 A^T . Then we multiply it by b .
00031 //   Numerically, not as good as 2 and not as fast . In general it is not a
00032 //   good idea to solve a set of linear equations with a matrix inversion .
00033 //
00034 //  4. Pseudo Inverse , brute force
00035 //   The pseudo inverse is calculated brute force through a series of matrix
00036 //   manipulations . It shows nicely some operations in the matrix package,
00037 //   but is otherwise a big "no no" .
00038 //
00039 //  5. Least-squares analysis with Minuit
00040 //   An objective function L is minimized by Minuit, where
00041 //    L = sum_i { (y - c_0 -c_1 * x / e)^2 }
00042 //   Minuit will calculate numerically the derivative of L wrt c_0 and c_1 .
00043 //   It has not been told that these derivatives are linear in the parameters
00044 //   c_0 and c_1 .
00045 //   For ill-conditioned linear problems it is better to use the fact it is
00046 //   a linear fit as in 2 .
00047 //
00048 // Another interesting thing is the way we assign data to the vectors and 
00049 // matrices through adoption .
00050 // This allows data assignment without physically moving bytes around .
00051 //
00052 //   USAGE
00053 //   -----
00054 // This macro can be execued via CINT or via ACLIC
00055 // - via CINT, do
00056 //    root > .x solveLinear.C
00057 // - via ACLIC
00058 //    root > gSystem->Load("libMatrix");
00059 //    root > gSystem->Load("libGpad");
00060 //    root > .x solveLinear.C+
00061 //
00062 #ifndef __CINT__
00063 #include "Riostream.h"
00064 #include "TMatrixD.h"
00065 #include "TVectorD.h"
00066 #include "TGraphErrors.h"
00067 #include "TDecompChol.h"
00068 #include "TDecompSVD.h"
00069 #include "TF1.h"
00070 #endif
00071 
00072 
00073 void solveLinear(Double_t eps = 1.e-12)
00074 {
00075 #ifdef __CINT__
00076   gSystem->Load("libMatrix");
00077 #endif
00078   cout << "Perform the fit  y = c0 + c1 * x in four different ways" << endl;
00079 
00080   const Int_t nrVar  = 2;
00081   const Int_t nrPnts = 4;
00082 
00083   Double_t ax[] = {0.0,1.0,2.0,3.0};
00084   Double_t ay[] = {1.4,1.5,3.7,4.1};
00085   Double_t ae[] = {0.5,0.2,1.0,0.5};
00086 
00087   // Make the vectors 'Use" the data : they are not copied, the vector data
00088   // pointer is just set appropriately
00089 
00090   TVectorD x; x.Use(nrPnts,ax);
00091   TVectorD y; y.Use(nrPnts,ay);
00092   TVectorD e; e.Use(nrPnts,ae);
00093 
00094   TMatrixD A(nrPnts,nrVar);
00095   TMatrixDColumn(A,0) = 1.0;
00096   TMatrixDColumn(A,1) = x;
00097 
00098   cout << " - 1. solve through Normal Equations" << endl;
00099 
00100   const TVectorD c_norm = NormalEqn(A,y,e);
00101 
00102   cout << " - 2. solve through SVD" << endl;
00103   // numerically  preferred method
00104 
00105   // first bring the weights in place
00106   TMatrixD Aw = A;
00107   TVectorD yw = y;
00108   for (Int_t irow = 0; irow < A.GetNrows(); irow++) {
00109     TMatrixDRow(Aw,irow) *= 1/e(irow);
00110     yw(irow) /= e(irow);
00111   }
00112 
00113   TDecompSVD svd(Aw);
00114   Bool_t ok;
00115   const TVectorD c_svd = svd.Solve(yw,ok);
00116 
00117   cout << " - 3. solve with pseudo inverse" << endl;
00118 
00119   const TMatrixD pseudo1  = svd.Invert();
00120   TVectorD c_pseudo1 = yw;
00121   c_pseudo1 *= pseudo1;
00122 
00123   cout << " - 4. solve with pseudo inverse, calculated brute force" << endl;
00124 
00125   TMatrixDSym AtA(TMatrixDSym::kAtA,Aw);
00126   const TMatrixD pseudo2 = AtA.Invert() * Aw.T();
00127   TVectorD c_pseudo2 = yw;
00128   c_pseudo2 *= pseudo2;
00129 
00130   cout << " - 5. Minuit through TGraph" << endl;
00131 
00132   TGraphErrors *gr = new TGraphErrors(nrPnts,ax,ay,0,ae);
00133   TF1 *f1 = new TF1("f1","pol1",0,5);
00134   gr->Fit("f1","Q");
00135   TVectorD c_graph(nrVar);
00136   c_graph(0) = f1->GetParameter(0);
00137   c_graph(1) = f1->GetParameter(1);
00138 
00139   // Check that all 4 answers are identical within a certain
00140   // tolerance . The 1e-12 is somewhat arbitrary . It turns out that
00141   // the TGraph fit is different by a few times 1e-13.
00142 
00143   Bool_t same = kTRUE;
00144   same &= VerifyVectorIdentity(c_norm,c_svd,0,eps);
00145   same &= VerifyVectorIdentity(c_norm,c_pseudo1,0,eps);
00146   same &= VerifyVectorIdentity(c_norm,c_pseudo2,0,eps);
00147   same &= VerifyVectorIdentity(c_norm,c_graph,0,eps);
00148   if (same)
00149     cout << " All solutions are the same within tolerance of " << eps << endl;
00150   else
00151     cout << " Some solutions differ more than the allowed tolerance of " << eps << endl;
00152 }

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