Logo ROOT   6.18/05
Reference Guide
minuit2FitBench2D.C File Reference

Detailed Description

View in nbviewer Open in SWAN

FCN=7196.63 FROM MIGRAD STATUS=CONVERGED 169 CALLS 170 TOTAL
EDM=3.06509e-09 STRATEGY= 1 ERROR MATRIX UNCERTAINTY 2.0 per cent
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 p0 5.28670e+01 2.67227e-01 1.71567e-03 1.41590e-04
2 p1 2.00562e+00 9.99521e-03 -6.72962e-05 -5.71949e-03
3 p2 -1.02279e+00 1.35440e-02 3.22465e-05 2.77027e-03
4 p3 2.89928e+00 7.94099e-03 -9.38703e-05 -4.55759e-03
5 p4 3.86656e+00 1.13330e-02 1.06568e-05 -2.11126e-03
Minuit, npass=0 : RT= 0.127 s, Cpu= 0.130 s
****************************************
Minimizer is Fumili
Chi2 = 7196.63
NDf = 7366
NCalls = 30
p0 = 52.8672 +/- 0.262932
p1 = 2.00562 +/- 0.00992527
p2 = -1.02279 +/- 0.0135821
p3 = 2.89927 +/- 0.00781808
p4 = 3.86655 +/- 0.011235
Fumili, npass=0 : RT= 0.189 s, Cpu= 0.190 s
****************************************
Minimizer is Minuit2
Chi2 = 7196.63
NDf = 7366
Edm = 4.35807e-08
NCalls = 176
p0 = 52.867 +/- 0.264728
p1 = 2.00562 +/- 0.00995152
p2 = -1.0228 +/- 0.0137398
p3 = 2.89928 +/- 0.00798634
p4 = 3.86656 +/- 0.0112306
Minuit2, npass=0 : RT= 0.108 s, Cpu= 0.110 s
****************************************
Minimizer is Minuit2 / Fumili
Chi2 = 7196.63
NDf = 7366
Edm = 2.49347e-07
NCalls = 99
p0 = 52.8669 +/- 0.265504
p1 = 2.00562 +/- 0.00992463
p2 = -1.0228 +/- 0.0135834
p3 = 2.89928 +/- 0.00792724
p4 = 3.86656 +/- 0.0113104
Fumili2, npass=0 : RT= 0.101 s, Cpu= 0.100 s
#include "TH1.h"
#include "TF1.h"
#include "TH2D.h"
#include "TF2.h"
#include "TCanvas.h"
#include "TStopwatch.h"
#include "TSystem.h"
#include "TRandom3.h"
#include "TVirtualFitter.h"
#include "TPaveLabel.h"
#include "TStyle.h"
TF2 *fitFcn;
TH2D *histo;
// Quadratic background function
Double_t gaus2D(Double_t *x, Double_t *par) {
double t1 = x[0] - par[1];
double t2 = x[1] - par[2];
return par[0]* exp( - 0.5 * ( t1*t1/( par[3]*par[3]) + t2*t2 /( par[4]*par[4] ) ) ) ;
}
// Sum of background and peak function
Double_t fitFunction(Double_t *x, Double_t *par) {
return gaus2D(x,par);
}
void fillHisto(int n =10000) {
gRandom = new TRandom3();
for (int i = 0; i < n; ++i) {
double x = gRandom->Gaus(2,3);
double y = gRandom->Gaus(-1,4);
histo->Fill(x,y,1.);
}
}
void DoFit(const char* fitter, TVirtualPad *pad, Int_t npass) {
TStopwatch timer;
pad->SetGrid();
fitFcn->SetParameters(100,0,0,2,7);
fitFcn->Update();
timer.Start();
histo->Fit("fitFcn","0");
timer.Stop();
histo->Draw();
Double_t cputime = timer.CpuTime();
printf("%s, npass=%d : RT=%7.3f s, Cpu=%7.3f s\n",fitter,npass,timer.RealTime(),cputime);
TPaveLabel *p = new TPaveLabel(0.5,0.7,0.85,0.8,Form("%s CPU= %g s",fitter,cputime),"brNDC");
p->Draw();
pad->Update();
}
void minuit2FitBench2D(int n = 100000) {
TCanvas *c1 = new TCanvas("c1","Fitting Demo",10,10,900,900);
c1->Divide(2,2);
// create a TF1 with the range from 0 to 3 and 6 parameters
fitFcn = new TF2("fitFcn",fitFunction,-10,10,-10,10,5);
//fitFcn->SetNpx(200);
histo = new TH2D("h2","2D Gauss",100,-10,10,100,-10,10);
fillHisto(n);
int npass=0;
//with Minuit
c1->cd(1);
DoFit("Minuit",gPad,npass);
//with Fumili
c1->cd(2);
DoFit("Fumili",gPad,npass);
//with Minuit2
c1->cd(3);
DoFit("Minuit2",gPad,npass);
//with Fumili2
c1->cd(4);
DoFit("Fumili2",gPad,npass);
}
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
double exp(double)
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
char * Form(const char *fmt,...)
R__EXTERN TStyle * gStyle
Definition: TStyle.h:406
#define gPad
Definition: TVirtualPad.h:286
The Canvas class.
Definition: TCanvas.h:31
virtual void Update()
Called by functions such as SetRange, SetNpx, SetParameters to force the deletion of the associated h...
Definition: TF1.cxx:3615
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
A 2-Dim function with parameters.
Definition: TF2.h:29
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition: TH1.cxx:1225
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:3791
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:289
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
A Pave (see TPave) with a text centered in the Pave.
Definition: TPaveLabel.h:20
virtual void Draw(Option_t *option="")
Draw this pavelabel with its current attributes.
Definition: TPaveLabel.cxx:77
Random number generator class based on M.
Definition: TRandom3.h:27
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:263
Stopwatch class.
Definition: TStopwatch.h:28
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
Definition: TStopwatch.cxx:110
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
Definition: TStopwatch.cxx:125
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
void SetStatY(Float_t y=0)
Definition: TStyle.h:375
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:1396
static void SetDefaultFitter(const char *name="")
static: set name of default fitter
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:50
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)=0
virtual void Update()=0
return c1
Definition: legend1.C:41
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
auto * t1
Definition: textangle.C:20
Author
Lorenzo Moneta

Definition in file minuit2FitBench2D.C.