ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
line3Dfit.C
Go to the documentation of this file.
1 //Fitting of a TGraph2D with a 3D straight line
2 //
3 // run this macro by doing:
4 //
5 // root>.x line3Dfit.C+
6 //
7 
8 //Author: L. Moneta
9 //
10 
11 
12 #include <TMath.h>
13 #include <TGraph2D.h>
14 #include <TRandom2.h>
15 #include <TStyle.h>
16 #include <TCanvas.h>
17 #include <TF2.h>
18 #include <TH1.h>
19 #include <Math/Functor.h>
20 #include <TPolyLine3D.h>
21 #include <Math/Vector3D.h>
22 #include <Fit/Fitter.h>
23 
24 using namespace ROOT::Math;
25 
26 
27 // define the parameteric line equation
28 void line(double t, const double *p, double &x, double &y, double &z) {
29  // a parameteric line is define from 6 parameters but 4 are independent
30  // x0,y0,z0,z1,y1,z1 which are the coordinates of two points on the line
31  // can choose z0 = 0 if line not parallel to x-y plane and z1 = 1;
32  x = p[0] + p[1]*t;
33  y = p[2] + p[3]*t;
34  z = t;
35 }
36 
37 // calculate distance line-point
38 double distance2(double x,double y,double z, const double *p) {
39  // distance line point is D= | (xp-x0) cross ux |
40  // where ux is direction of line and x0 is a point in the line (like t = 0)
41  XYZVector xp(x,y,z);
42  XYZVector x0(p[0], p[2], 0. );
43  XYZVector x1(p[0] + p[1], p[2] + p[3], 1. );
44  XYZVector u = (x1-x0).Unit();
45  double d2 = ((xp-x0).Cross(u)) .Mag2();
46  return d2;
47 }
48 bool first = true;
49 
50 
51 // function Object to be minimized
52 struct SumDistance2 {
53  // the TGraph is a data member of the object
54  TGraph2D * fGraph;
55 
56  SumDistance2(TGraph2D * g) : fGraph(g) {}
57 
58  // implementation of the function to be minimized
59  double operator() (const double * par) {
60 
61  assert(fGraph != 0);
62  double * x = fGraph->GetX();
63  double * y = fGraph->GetY();
64  double * z = fGraph->GetZ();
65  int npoints = fGraph->GetN();
66  double sum = 0;
67  for (int i = 0; i < npoints; ++i) {
68  double d = distance2(x[i],y[i],z[i],par);
69  sum += d;
70 #ifdef DEBUG
71  if (first) std::cout << "point " << i << "\t"
72  << x[i] << "\t"
73  << y[i] << "\t"
74  << z[i] << "\t"
75  << std::sqrt(d) << std::endl;
76 #endif
77  }
78  if (first)
79  std::cout << "Total Initial distance square = " << sum << std::endl;
80  first = false;
81  return sum;
82  }
83 
84 };
85 
87 {
88  gStyle->SetOptStat(0);
89  gStyle->SetOptFit();
90 
91 
92  //double e = 0.1;
93  Int_t nd = 10000;
94 
95 
96 // double xmin = 0; double ymin = 0;
97 // double xmax = 10; double ymax = 10;
98 
99  TGraph2D * gr = new TGraph2D();
100 
101  // Fill the 2D graph
102  double p0[4] = {10,20,1,2};
103 
104  // generate graph with the 3d points
105  for (Int_t N=0; N<nd; N++) {
106  double x,y,z = 0;
107  // Generate a random number
108  double t = gRandom->Uniform(0,10);
109  line(t,p0,x,y,z);
110  double err = 1;
111  // do a gaussian smearing around the points in all coordinates
112  x += gRandom->Gaus(0,err);
113  y += gRandom->Gaus(0,err);
114  z += gRandom->Gaus(0,err);
115  gr->SetPoint(N,x,y,z);
116  //dt->SetPointError(N,0,0,err);
117  }
118  // fit the graph now
119 
120  ROOT::Fit::Fitter fitter;
121 
122 
123  // make the functor objet
124  SumDistance2 sdist(gr);
125 #ifdef __CINT__
126  ROOT::Math::Functor fcn(&sdist,4,"SumDistance2");
127 #else
128  ROOT::Math::Functor fcn(sdist,4);
129 #endif
130  // set the function and the initial parameter values
131  double pStart[4] = {1,1,1,1};
132  fitter.SetFCN(fcn,pStart);
133  // set step sizes different than default ones (0.3 times parameter values)
134  for (int i = 0; i < 4; ++i) fitter.Config().ParSettings(i).SetStepSize(0.01);
135 
136  bool ok = fitter.FitFCN();
137  if (!ok) {
138  Error("line3Dfit","Line3D Fit failed");
139  return 1;
140  }
141 
142  const ROOT::Fit::FitResult & result = fitter.Result();
143 
144  std::cout << "Total final distance square " << result.MinFcnValue() << std::endl;
145  result.Print(std::cout);
146 
147 
148  gr->Draw("p0");
149 
150  // get fit parameters
151  const double * parFit = result.GetParams();
152 
153 
154  // draw the fitted line
155  int n = 1000;
156  double t0 = 0;
157  double dt = 10;
158  TPolyLine3D *l = new TPolyLine3D(n);
159  for (int i = 0; i <n;++i) {
160  double t = t0+ dt*i/n;
161  double x,y,z;
162  line(t,parFit,x,y,z);
163  l->SetPoint(i,x,y,z);
164  }
165  l->SetLineColor(kRed);
166  l->Draw("same");
167 
168  // draw original line
169  TPolyLine3D *l0 = new TPolyLine3D(n);
170  for (int i = 0; i <n;++i) {
171  double t = t0+ dt*i/n;
172  double x,y,z;
173  line(t,p0,x,y,z);
174  l0->SetPoint(i,x,y,z);
175  }
176  l0->SetLineColor(kBlue);
177  l0->Draw("same");
178  return 0;
179 }
180 
181 int main() {
182  return line3Dfit();
183 }
const FitConfig & Config() const
access to the fit configuration (const method)
Definition: Fitter.h:382
void line(double t, const double *p, double &x, double &y, double &z)
Definition: line3Dfit.C:28
double par[1]
Definition: unuranDistr.cxx:38
void fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
Definition: Ifit.C:26
virtual void Draw(Option_t *option="")
Specific drawing options can be used to paint a TGraph2D:
Definition: TGraph2D.cxx:719
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
Definition: Rtypes.h:61
Documentation for class Functor class.
Definition: Functor.h:394
#define assert(cond)
Definition: unittest.h:542
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
bool FitFCN(unsigned int npar, Function &fcn, const double *params=0, unsigned int dataSize=0, bool chi2fit=false)
Fit using the a generic FCN function as a C++ callable object implementing double () (const double *)...
Definition: Fitter.h:538
#define N
int Int_t
Definition: RtypesCore.h:41
A 3-dimensional polyline.
Definition: TPolyLine3D.h:41
double sqrt(double)
Double_t x[n]
Definition: legend1.C:17
int d
Definition: tornado.py:11
const FitResult & Result() const
get fit result
Definition: Fitter.h:354
double distance2(double x, double y, double z, const double *p)
Definition: line3Dfit.C:38
Float_t z[5]
Definition: Ifit.C:16
Int_t line3Dfit()
Definition: line3Dfit.C:86
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
TThread * t[5]
Definition: threadsh1.C:13
void SetStepSize(double err)
set the step size
const double * GetParams() const
parameter values (return const pointer)
Definition: FitResult.h:181
Fitter class, entry point for performing all type of fits.
Definition: Fitter.h:94
const ParameterSettings & ParSettings(unsigned int i) const
get the parameter settings for the i-th parameter (const method)
Definition: FitConfig.h:80
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1204
int main()
Definition: line3Dfit.C:181
bool first
Definition: line3Dfit.C:48
TLine * l
Definition: textangle.C:4
T Mag2(const SVector< T, D > &rhs)
Vector magnitude square Template to compute .
Definition: Functions.h:231
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
TGraphErrors * gr
Definition: legend1.C:25
virtual void Draw(Option_t *option="")
Draw this 3-D polyline with its current attributes.
SVector< T, 3 > Cross(const SVector< T, 3 > &lhs, const SVector< T, 3 > &rhs)
Vector Cross Product (only for 3-dim vectors) .
Definition: Functions.h:324
class containg the result of the fit and all the related information (fitted parameter values...
Definition: FitResult.h:52
TRObject operator()(const T1 &t1) const
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
static const double x1[5]
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Sets point number n.
Definition: TGraph2D.cxx:1707
Double_t y[n]
Definition: legend1.C:17
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
Definition: FitResult.cxx:446
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:606
double MinFcnValue() const
Return value of the objective function (chi2 or likelihood) used in the fit.
Definition: FitResult.h:125
bool SetFCN(unsigned int npar, Function &fcn, const double *params=0, unsigned int dataSize=0, bool chi2fit=false)
Set a generic FCN function as a C++ callable object implementing double () (const double *) Note that...
Definition: Fitter.h:543
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1252
double result[121]
Definition: Rtypes.h:61
SVector< T, D > Unit(const SVector< T, D > &rhs)
Unit.
Definition: Functions.h:383
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:50
const Int_t n
Definition: legend1.C:16
void Error(ErrorHandler_t func, int code, const char *va_(fmt),...)
Write error message and call a handler, if required.