ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
solveLinear.C
Go to the documentation of this file.
1 //Author: Eddy Offermann
2 // This macro shows several ways to perform a linear least-squares
3 // analysis . To keep things simple we fit a straight line to 4
4 // data points
5 // The first 4 methods use the linear algebra package to find
6 // x such that min (A x - b)^T (A x - b) where A and b
7 // are calculated with the data points and the functional expression :
8 //
9 // 1. Normal equations:
10 // Expanding the expression (A x - b)^T (A x - b) and taking the
11 // derivative wrt x leads to the "Normal Equations":
12 // A^T A x = A^T b where A^T A is a positive definite matrix. Therefore,
13 // a Cholesky decomposition scheme can be used to calculate its inverse .
14 // This leads to the solution x = (A^T A)^-1 A^T b . All this is done in
15 // routine NormalEqn . We made it a bit more complicated by giving the
16 // data weights .
17 // Numerically this is not the best way to proceed because effctively the
18 // condition number of A^T A is twice as large as that of A, making inversion
19 // more difficult
20 //
21 // 2. SVD
22 // One can show that solving A x = b for x with A of size (m x n) and m > n
23 // through a Singular Value Decomposition is equivalent to miminizing
24 // (A x - b)^T (A x - b)
25 // Numerically , this is the most stable method of all 5
26 //
27 // 3. Pseudo Inverse
28 // Here we calulate the generalized matrix inverse ("pseudo inverse") by
29 // solving A X = Unit for matrix X through an SVD . The formal expression for
30 // is X = (A^T A)^-1 A^T . Then we multiply it by b .
31 // Numerically, not as good as 2 and not as fast . In general it is not a
32 // good idea to solve a set of linear equations with a matrix inversion .
33 //
34 // 4. Pseudo Inverse , brute force
35 // The pseudo inverse is calculated brute force through a series of matrix
36 // manipulations . It shows nicely some operations in the matrix package,
37 // but is otherwise a big "no no" .
38 //
39 // 5. Least-squares analysis with Minuit
40 // An objective function L is minimized by Minuit, where
41 // L = sum_i { (y - c_0 -c_1 * x / e)^2 }
42 // Minuit will calculate numerically the derivative of L wrt c_0 and c_1 .
43 // It has not been told that these derivatives are linear in the parameters
44 // c_0 and c_1 .
45 // For ill-conditioned linear problems it is better to use the fact it is
46 // a linear fit as in 2 .
47 //
48 // Another interesting thing is the way we assign data to the vectors and
49 // matrices through adoption .
50 // This allows data assignment without physically moving bytes around .
51 //
52 // USAGE
53 // -----
54 // This macro can be execued via CINT or via ACLIC
55 // - via CINT, do
56 // root > .x solveLinear.C
57 // - via ACLIC
58 // root > gSystem->Load("libMatrix");
59 // root > gSystem->Load("libGpad");
60 // root > .x solveLinear.C+
61 //
62 #ifndef __CINT__
63 #include "Riostream.h"
64 #include "TMatrixD.h"
65 #include "TVectorD.h"
66 #include "TGraphErrors.h"
67 #include "TDecompChol.h"
68 #include "TDecompSVD.h"
69 #include "TF1.h"
70 #endif
71 
72 
73 void solveLinear(Double_t eps = 1.e-12)
74 {
75 #ifdef __CINT__
76  gSystem->Load("libMatrix");
77 #endif
78  cout << "Perform the fit y = c0 + c1 * x in four different ways" << endl;
79 
80  const Int_t nrVar = 2;
81  const Int_t nrPnts = 4;
82 
83  Double_t ax[] = {0.0,1.0,2.0,3.0};
84  Double_t ay[] = {1.4,1.5,3.7,4.1};
85  Double_t ae[] = {0.5,0.2,1.0,0.5};
86 
87  // Make the vectors 'Use" the data : they are not copied, the vector data
88  // pointer is just set appropriately
89 
90  TVectorD x; x.Use(nrPnts,ax);
91  TVectorD y; y.Use(nrPnts,ay);
92  TVectorD e; e.Use(nrPnts,ae);
93 
94  TMatrixD A(nrPnts,nrVar);
95  TMatrixDColumn(A,0) = 1.0;
96  TMatrixDColumn(A,1) = x;
97 
98  cout << " - 1. solve through Normal Equations" << endl;
99 
100  const TVectorD c_norm = NormalEqn(A,y,e);
101 
102  cout << " - 2. solve through SVD" << endl;
103  // numerically preferred method
104 
105  // first bring the weights in place
106  TMatrixD Aw = A;
107  TVectorD yw = y;
108  for (Int_t irow = 0; irow < A.GetNrows(); irow++) {
109  TMatrixDRow(Aw,irow) *= 1/e(irow);
110  yw(irow) /= e(irow);
111  }
112 
113  TDecompSVD svd(Aw);
114  Bool_t ok;
115  const TVectorD c_svd = svd.Solve(yw,ok);
116 
117  cout << " - 3. solve with pseudo inverse" << endl;
118 
119  const TMatrixD pseudo1 = svd.Invert();
120  TVectorD c_pseudo1 = yw;
121  c_pseudo1 *= pseudo1;
122 
123  cout << " - 4. solve with pseudo inverse, calculated brute force" << endl;
124 
126  const TMatrixD pseudo2 = AtA.Invert() * Aw.T();
127  TVectorD c_pseudo2 = yw;
128  c_pseudo2 *= pseudo2;
129 
130  cout << " - 5. Minuit through TGraph" << endl;
131 
132  TGraphErrors *gr = new TGraphErrors(nrPnts,ax,ay,0,ae);
133  TF1 *f1 = new TF1("f1","pol1",0,5);
134  gr->Fit("f1","Q");
135  TVectorD c_graph(nrVar);
136  c_graph(0) = f1->GetParameter(0);
137  c_graph(1) = f1->GetParameter(1);
138 
139  // Check that all 4 answers are identical within a certain
140  // tolerance . The 1e-12 is somewhat arbitrary . It turns out that
141  // the TGraph fit is different by a few times 1e-13.
142 
143  Bool_t same = kTRUE;
144  same &= VerifyVectorIdentity(c_norm,c_svd,0,eps);
145  same &= VerifyVectorIdentity(c_norm,c_pseudo1,0,eps);
146  same &= VerifyVectorIdentity(c_norm,c_pseudo2,0,eps);
147  same &= VerifyVectorIdentity(c_norm,c_graph,0,eps);
148  if (same)
149  cout << " All solutions are the same within tolerance of " << eps << endl;
150  else
151  cout << " Some solutions differ more than the allowed tolerance of " << eps << endl;
152 }
TVectorD NormalEqn(const TMatrixD &A, const TVectorD &b)
Solve min {(A .
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0)
Fit this graph with function with name fname.
Definition: TGraph.cxx:1024
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1766
TMatrixT< Element > & T()
Definition: TMatrixT.h:150
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
static double A[]
TVectorT< Element > & Use(Int_t lwb, Int_t upb, Element *data)
Use the array data to fill the vector lwb..upb].
Definition: TVectorT.cxx:346
Double_t x[n]
Definition: legend1.C:17
TMatrixTRow< Double_t > TMatrixDRow
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
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
void solveLinear(Double_t eps=1.e-12)
Definition: solveLinear.C:73
TGraphErrors * gr
Definition: legend1.C:25
virtual Bool_t Solve(TVectorD &b)
Solve Ax=b assuming the SVD form of A is stored .
Definition: TDecompSVD.cxx:609
double Double_t
Definition: RtypesCore.h:55
Double_t y[n]
Definition: legend1.C:17
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:352
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
TMatrixTColumn< Double_t > TMatrixDColumn
1-Dim function class
Definition: TF1.h:149
TF1 * f1
Definition: legend1.C:11
A TGraphErrors is a TGraph with error bars.
Definition: TGraphErrors.h:28
Bool_t VerifyVectorIdentity(const TVectorT< Element > &m1, const TVectorT< Element > &m2, Int_t verbose, Element maxDevAllow)
Verify that elements of the two vectors are equal within maxDevAllow .
Definition: TVectorT.cxx:2294
const Bool_t kTRUE
Definition: Rtypes.h:91