Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
minuit2FitBench2D.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_fit
3 /// \notebook -js
4 ///
5 /// \macro_image
6 /// \macro_output
7 /// \macro_code
8 ///
9 /// \author Lorenzo Moneta
10 
11 #include "TH1.h"
12 #include "TF1.h"
13 #include "TH2D.h"
14 #include "TF2.h"
15 #include "TCanvas.h"
16 #include "TStopwatch.h"
17 #include "TSystem.h"
18 #include "TRandom3.h"
19 #include "TVirtualFitter.h"
20 #include "TPaveLabel.h"
21 #include "TStyle.h"
22 
23 
24 TF2 *fitFcn;
25 TH2D *histo;
26 
27 // Quadratic background function
28 Double_t gaus2D(Double_t *x, Double_t *par) {
29  double t1 = x[0] - par[1];
30  double t2 = x[1] - par[2];
31  return par[0]* exp( - 0.5 * ( t1*t1/( par[3]*par[3]) + t2*t2 /( par[4]*par[4] ) ) ) ;
32 }
33 
34 // Sum of background and peak function
35 Double_t fitFunction(Double_t *x, Double_t *par) {
36  return gaus2D(x,par);
37 }
38 
39 void fillHisto(int n =10000) {
40 
41  gRandom = new TRandom3();
42  for (int i = 0; i < n; ++i) {
43  double x = gRandom->Gaus(2,3);
44  double y = gRandom->Gaus(-1,4);
45  histo->Fill(x,y,1.);
46  }
47 }
48 
49 void DoFit(const char* fitter, TVirtualPad *pad, Int_t npass) {
50  TStopwatch timer;
51  TVirtualFitter::SetDefaultFitter(fitter);
52  pad->SetGrid();
53  fitFcn->SetParameters(100,0,0,2,7);
54  fitFcn->Update();
55 
56  timer.Start();
57  histo->Fit("fitFcn","0");
58  timer.Stop();
59 
60  histo->Draw();
61  Double_t cputime = timer.CpuTime();
62  printf("%s, npass=%d : RT=%7.3f s, Cpu=%7.3f s\n",fitter,npass,timer.RealTime(),cputime);
63  TPaveLabel *p = new TPaveLabel(0.5,0.7,0.85,0.8,Form("%s CPU= %g s",fitter,cputime),"brNDC");
64  p->Draw();
65  pad->Update();
66 }
67 
68 void minuit2FitBench2D(int n = 100000) {
69  TH1::AddDirectory(kFALSE);
70  TCanvas *c1 = new TCanvas("c1","Fitting Demo",10,10,900,900);
71  c1->Divide(2,2);
72  // create a TF1 with the range from 0 to 3 and 6 parameters
73  fitFcn = new TF2("fitFcn",fitFunction,-10,10,-10,10,5);
74  //fitFcn->SetNpx(200);
75  gStyle->SetOptFit();
76  gStyle->SetStatY(0.6);
77 
78 
79  histo = new TH2D("h2","2D Gauss",100,-10,10,100,-10,10);
80  fillHisto(n);
81 
82  int npass=0;
83 
84  //with Minuit
85  c1->cd(1);
86  DoFit("Minuit",gPad,npass);
87 
88  //with Fumili
89  c1->cd(2);
90  DoFit("Fumili",gPad,npass);
91 
92  //with Minuit2
93  c1->cd(3);
94  DoFit("Minuit2",gPad,npass);
95 
96  //with Fumili2
97  c1->cd(4);
98  DoFit("Fumili2",gPad,npass);
99 }