ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
minuit2FitBench2D.C
Go to the documentation of this file.
1 // @(#)root/minuit2:$Id$
2 // Author: L. Moneta 10/2005
3 
4 /**********************************************************************
5  * *
6  * Copyright (c) 2005 ROOT Foundation, CERN/PH-SFT *
7  * *
8  **********************************************************************/
9 
10 #include "TH1.h"
11 #include "TF1.h"
12 #include "TH2D.h"
13 #include "TF2.h"
14 #include "TCanvas.h"
15 #include "TStopwatch.h"
16 #include "TSystem.h"
17 #include "TRandom3.h"
18 #include "TVirtualFitter.h"
19 #include "TPaveLabel.h"
20 #include "TStyle.h"
21 
22 
25 
26 // Quadratic background function
28  double t1 = x[0] - par[1];
29  double t2 = x[1] - par[2];
30  return par[0]* exp( - 0.5 * ( t1*t1/( par[3]*par[3]) + t2*t2 /( par[4]*par[4] ) ) ) ;
31 }
32 
33 // Sum of background and peak function
35  return gaus2D(x,par);
36 }
37 
38 void fillHisto(int n =10000) {
39 
40  gRandom = new TRandom3();
41  for (int i = 0; i < n; ++i) {
42  double x = gRandom->Gaus(2,3);
43  double y = gRandom->Gaus(-1,4);
44  histo->Fill(x,y,1.);
45  }
46 }
47 
48 void DoFit(const char* fitter, TVirtualPad *pad, Int_t npass) {
51  pad->SetGrid();
52  fitFcn->SetParameters(100,0,0,2,7);
53  fitFcn->Update();
54 
55  timer.Start();
56  histo->Fit("fitFcn","0");
57  timer.Stop();
58 
59  histo->Draw();
60  Double_t cputime = timer.CpuTime();
61  printf("%s, npass=%d : RT=%7.3f s, Cpu=%7.3f s\n",fitter,npass,timer.RealTime(),cputime);
62  TPaveLabel *p = new TPaveLabel(0.5,0.7,0.85,0.8,Form("%s CPU= %g s",fitter,cputime),"brNDC");
63  p->Draw();
64  pad->Update();
65 }
66 
67 void minuit2FitBench2D(int n = 100000) {
69  TCanvas *c1 = new TCanvas("c1","Fitting Demo",10,10,900,900);
70  c1->Divide(2,2);
71  // create a TF1 with the range from 0 to 3 and 6 parameters
72  fitFcn = new TF2("fitFcn",fitFunction,-10,10,-10,10,5);
73  //fitFcn->SetNpx(200);
74  gStyle->SetOptFit();
75  gStyle->SetStatY(0.6);
76 
77 
78  histo = new TH2D("h2","2D Gauss",100,-10,10,100,-10,10);
79  fillHisto(n);
80 
81  int npass=0;
82 
83  //with Minuit
84  c1->cd(1);
85  DoFit("Minuit",gPad,npass);
86 
87  //with Fumili
88  c1->cd(2);
89  DoFit("Fumili",gPad,npass);
90 
91  //with Minuit2
92  c1->cd(3);
93  DoFit("Minuit2",gPad,npass);
94 
95  //with Fumili2
96  c1->cd(4);
97  DoFit("Fumili2",gPad,npass);
98 }
double par[1]
Definition: unuranDistr.cxx:38
const int npass
Definition: testPermute.cxx:21
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
Random number generator class based on M.
Definition: TRandom3.h:29
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
Definition: TStopwatch.cxx:108
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:56
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
TCanvas * c1
Definition: legend1.C:2
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
virtual void Update()=0
static void SetDefaultFitter(const char *name="")
static: set name of default fitter
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
Definition: TStopwatch.cxx:123
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: Rtypes.h:92
void SetStatY(Float_t y=0)
Definition: TStyle.h:392
void minuit2FitBench2D(int n=100000)
TLatex * t1
Definition: textangle.C:20
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition: TH1.cxx:1231
TStopwatch timer
Definition: pirndm.C:37
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:75
Double_t x[n]
Definition: legend1.C:17
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:59
A Pave (see TPave) with a text centered in the Pave.
Definition: TPaveLabel.h:32
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
TF2 * fitFcn
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
char * Form(const char *fmt,...)
virtual void Update()
Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated h...
Definition: TF1.cxx:3337
tuple pad
Definition: first.py:38
A 2-Dim function with parameters.
Definition: TF2.h:33
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
TH2D * histo
TText * t2
Definition: rootenv.C:28
The Canvas class.
Definition: TCanvas.h:48
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Double_t y[n]
Definition: legend1.C:17
Double_t gaus2D(Double_t *x, Double_t *par)
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1073
void fillHisto(int n=10000)
Double_t fitFunction(Double_t *x, Double_t *par)
#define gPad
Definition: TVirtualPad.h:288
double exp(double)
void DoFit(const char *fitter, TVirtualPad *pad, Int_t npass)
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:287
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3607
const Int_t n
Definition: legend1.C:16
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)=0
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:297
Stopwatch class.
Definition: TStopwatch.h:30