Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
vectorizedFit.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Tutorial for creating a Vectorized TF1 function using a formula expression and use it for fitting an histogram

To create a vectorized function (if ROOT has been compiled with support for vectorization) is very easy. One needs to create the TF1 object with the option "VEC" or call the method TF1::SetVectorized

Doing Serial Gaussian Fit
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 40217.9
NDf = 39409
Edm = 3.38976e-08
NCalls = 75
Constant = 120.018 +/- 0.105817
Mean = 0.00114402 +/- 0.000709328
Sigma = 0.979817 +/- 0.000519995 (limited)
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 20355.3
Chi2 = 40710.5
NDf = 39997
Edm = 8.97826e-10
NCalls = 65
Constant = 120.024 +/- 0.105551
Mean = 0.000138332 +/- 0.000716607
Sigma = 0.99985 +/- 0.000537073 (limited)
Real time 0:00:00, CP time 0.450
Doing Vectorized Gaussian Fit
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 40217.9
NDf = 39409
Edm = 3.38976e-08
NCalls = 75
Constant = 120.018 +/- 0.105817
Mean = 0.00114402 +/- 0.000709328
Sigma = 0.979817 +/- 0.000519995 (limited)
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 20355.3
Chi2 = 40710.5
NDf = 39997
Edm = 8.97826e-10
NCalls = 65
Constant = 120.024 +/- 0.105551
Mean = 0.000138332 +/- 0.000716607
Sigma = 0.99985 +/- 0.000537073 (limited)
Real time 0:00:00, CP time 0.360
Doing Serial Polynomial Fit
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 37690
NDf = 38075
Edm = 3.47827e-15
NCalls = 72
A = 0.202001 +/- 0.00176461
B = 0.268032 +/- 0.0153893
C = 1.05504 +/- 0.0248331
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 20527.5
Chi2 = 41055
NDf = 39997
Edm = 8.16351e-08
NCalls = 90
A = 0.149763 +/- 0.00165111
B = 0.880262 +/- 0.0135143
C = 0.6066 +/- 0.0181308
Real time 0:00:00, CP time 0.330
Doing Vectorized Polynomial Fit
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 37690
NDf = 38075
Edm = 3.15424e-16
NCalls = 70
A = 0.202001 +/- 0.00176461
B = 0.268032 +/- 0.0153893
C = 1.05504 +/- 0.0248331
****************************************
Minimizer is Minuit2 / Migrad
MinFCN = 20527.5
Chi2 = 41055
NDf = 39997
Edm = 8.16253e-08
NCalls = 90
A = 0.149763 +/- 0.0016511
B = 0.880262 +/- 0.0135143
C = 0.6066 +/- 0.0181308
Real time 0:00:00, CP time 0.300
#include <TCanvas.h>
#include <TF1.h>
#include <TH1D.h>
#include <TStopwatch.h>
#include <TStyle.h>
#include <iostream>
void vectorizedFit() {
gStyle->SetOptFit(111111);
int nbins = 40000;
auto h1 = new TH1D("h1","h1",nbins,-3,3);
h1->FillRandom("gaus",nbins*50);
auto c1 = new TCanvas("Fit","Fit",800,1000);
c1->Divide(1,2);
c1->cd(1);
std::cout << "Doing Serial Gaussian Fit " << std::endl;
auto f1 = new TF1("f1","gaus");
f1->SetNpx(nbins*10);
w.Start();
h1->Fit(f1);
h1->Fit(f1,"L+");
w.Print();
std::cout << "Doing Vectorized Gaussian Fit " << std::endl;
auto f2 = new TF1("f2","gaus",-3,3,"VEC");
// alternatively you can also use the TF1::SetVectorized function
//f2->SetVectorized(true);
w.Start();
h1->Fit(f2);
h1->Fit(f2,"L+");
w.Print();
// rebin histograms and scale it back to the function
h1->Rebin(nbins/100);
h1->Scale(100./nbins);
((TF1 *)h1->GetListOfFunctions()->At(0))->SetTitle("Chi2 Fit");
((TF1 *)h1->GetListOfFunctions()->At(1))->SetTitle("Likelihood Fit");
//c1->cd(1)->BuildLegend();
/// Do a polynomial fit now
c1->cd(2);
auto f3 = new TF1("f3","[A]*x^2+[B]*x+[C]",0,10);
f3->SetParameters(0.5,3,2);
f3->SetNpx(nbins*10);
// generate the events
auto h2 = new TH1D("h2","h2",nbins,0,10);
h2->FillRandom("f3",10*nbins);
std::cout << "Doing Serial Polynomial Fit " << std::endl;
f3->SetParameters(2,2,2);
w.Start();
h2->Fit(f3);
h2->Fit(f3,"L+");
w.Print();
std::cout << "Doing Vectorized Polynomial Fit " << std::endl;
auto f4 = new TF1("f4","[A]*x*x+[B]*x+[C]",0,10);
f4->SetVectorized(true);
f4->SetParameters(2,2,2);
w.Start();
h2->Fit(f4);
h2->Fit(f4,"L+");
w.Print();
// rebin histograms and scale it back to the function
h2->Rebin(nbins/100);
h2->Scale(100./nbins);
((TF1 *)h2->GetListOfFunctions()->At(0))->SetTitle("Chi2 Fit");
((TF1 *)h2->GetListOfFunctions()->At(1))->SetTitle("Likelihood Fit");
((TF1 *)h2->GetListOfFunctions()->At(1))->SetLineColor(kBlue);
//c1->cd(2)->BuildLegend();
}
@ kBlue
Definition Rtypes.h:66
Option_t Option_t SetLineColor
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
static void SetDefaultMinimizer(const char *type, const char *algo=nullptr)
Set the default Minimizer type and corresponding algorithms.
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:233
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition TF1.cxx:3433
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:670
virtual void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr)
Fill histogram following distribution in function fname.
Definition TH1.cxx:3519
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:3898
virtual TH1 * Rebin(Int_t ngroup=2, const char *newname="", const Double_t *xbins=nullptr)
Rebin this histogram.
Definition TH1.cxx:6275
TList * GetListOfFunctions() const
Definition TH1.h:244
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition TH1.cxx:6604
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:355
Stopwatch class.
Definition TStopwatch.h:28
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:1593
return c1
Definition legend1.C:41
TH1F * h1
Definition legend1.C:5
TF1 * f1
Definition legend1.C:11
Author
Lorenzo Moneta

Definition in file vectorizedFit.C.