Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
graph2dfit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_fit
3 /// \notebook
4 /// Fitting a TGraph2D
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \author Olivier Couet
11 
12 #include <TMath.h>
13 #include <TGraph2D.h>
14 #include <TRandom.h>
15 #include <TStyle.h>
16 #include <TCanvas.h>
17 #include <TF2.h>
18 #include <TH1.h>
19 
20 TCanvas* graph2dfit()
21 {
22  gStyle->SetOptStat(0);
23  gStyle->SetOptFit();
24 
25  TCanvas *c = new TCanvas("c","Graph2D example",0,0,600,800);
26  c->Divide(2,3);
27 
28  Double_t rnd, x, y, z;
29  Double_t e = 0.3;
30  Int_t nd = 400;
31  Int_t np = 10000;
32 
33  TRandom r;
34  Double_t fl = 6;
35  TF2 *f2 = new TF2("f2","1000*(([0]*sin(x)/x)*([1]*sin(y)/y))+200",
36  -fl,fl,-fl,fl);
37  f2->SetParameters(1,1);
38  TGraph2D *dt = new TGraph2D();
39 
40  // Fill the 2D graph
41  Double_t zmax = 0;
42  for (Int_t N=0; N<nd; N++) {
43  f2->GetRandom2(x,y);
44  // Generate a random number in [-e,e]
45  rnd = 2*r.Rndm()*e-e;
46  z = f2->Eval(x,y)*(1+rnd);
47  if (z>zmax) zmax = z;
48  dt->SetPoint(N,x,y,z);
49  }
50 
51  Double_t hr = 350;
52  TH1D *h1 = new TH1D("h1",
53  "#splitline{Difference between Original}{#splitline{function and Function}{with noise}}",
54  100, -hr, hr);
55  TH1D *h2 = new TH1D("h2",
56  "#splitline{Difference between Original}{#splitline{function and Delaunay triangles}{interpolation}}",
57  100, -hr, hr);
58  TH1D *h3 = new TH1D("h3",
59  "#splitline{Difference between Original}{function and Minuit fit}",
60  500, -hr, hr);
61 
62  f2->SetParameters(0.5,1.5);
63  dt->Fit(f2);
64  TF2 *fit2 = (TF2*)dt->FindObject("f2");
65 
66  f2->SetParameters(1,1);
67 
68  for (Int_t N=0; N<np; N++) {
69  f2->GetRandom2(x,y);
70  // Generate a random number in [-e,e]
71  rnd = 2*r.Rndm()*e-e;
72  z = f2->Eval(x,y)*(1+rnd);
73  h1->Fill(f2->Eval(x,y)-z);
74  z = dt->Interpolate(x,y);
75  h2->Fill(f2->Eval(x,y)-z);
76  z = fit2->Eval(x,y);
77  h3->Fill(f2->Eval(x,y)-z);
78  }
79 
80  c->cd(1);
81  f2->SetTitle("Original function with Graph2D points on top");
82  f2->SetMaximum(zmax);
83  gStyle->SetHistTopMargin(0);
84  f2->Draw("surf1");
85  dt->Draw("same p0");
86 
87  c->cd(3);
88  dt->SetMargin(0.1);
89  dt->SetFillColor(36);
90  dt->SetTitle("Histogram produced with Delaunay interpolation");
91  dt->Draw("surf4");
92 
93  c->cd(5);
94  fit2->SetTitle("Minuit fit result on the Graph2D points");
95  fit2->Draw("surf1");
96 
97  h1->SetFillColor(47);
98  h2->SetFillColor(38);
99  h3->SetFillColor(29);
100 
101  c->cd(2); h1->Fit("gaus","Q") ; h1->Draw();
102  c->cd(4); h2->Fit("gaus","Q") ; h2->Draw();
103  c->cd(6); h3->Fit("gaus","Q") ; h3->Draw();
104  c->cd();
105  return c;
106 }