Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
FittingDemo.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Example for fitting signal/background.

This example can be executed with:

root > .x FittingDemo.C (using the cling interpreter)
root > .x FittingDemo.C+ (using the native complier via ACLIC)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 58.9284
NDf = 54
Edm = 9.73335e-07
NCalls = 606
p0 = -0.864746 +/- 0.887832
p1 = 45.8428 +/- 2.65453
p2 = -13.3213 +/- 0.980305
p3 = 13.8087 +/- 2.24457
p4 = 0.172313 +/- 0.0374057
p5 = 0.987278 +/- 0.0112928
Minuit2Minimizer: Minimize with max-calls 1780 convergence for edm < 0.01 strategy 1
Number of iterations 10
----------> Iteration 0
FVAL = 60.856975016 Edm = 2.04157626492 Nfcn = 23
Error matrix change = 1
Parameters : p0 = -0.864746 p1 = 45.8428 p2 = -13.3213 p3 = 13.8087 p4 = 0.2 p5 = 1
----------> Iteration 1
FVAL = 59.0984738211 Edm = 0.192336813475 Nfcn = 40
Error matrix change = 0.589159
Parameters : p0 = -0.913462 p1 = 45.8154 p2 = -13.3295 p3 = 14.3068 p4 = 0.176677 p5 = 0.990619
----------> Iteration 2
FVAL = 58.9599783699 Edm = 0.018918033963 Nfcn = 54
Error matrix change = 0.342284
Parameters : p0 = -0.906747 p1 = 45.8363 p2 = -13.3167 p3 = 14.116 p4 = 0.177064 p5 = 0.986977
----------> Iteration 3
FVAL = 58.9376425149 Edm = 0.00958566000238 Nfcn = 68
Error matrix change = 0.364189
Parameters : p0 = -0.924256 p1 = 45.83 p2 = -13.3146 p3 = 13.9161 p4 = 0.174171 p5 = 0.987051
----------> Iteration 4
FVAL = 58.9318510093 Edm = 0.00157028886174 Nfcn = 82
Error matrix change = 0.246637
Parameters : p0 = -0.905334 p1 = 45.839 p2 = -13.3101 p3 = 13.8652 p4 = 0.172941 p5 = 0.987268
----------> Iteration 5
FVAL = 58.9287717045 Edm = 0.000439744580493 Nfcn = 96
Error matrix change = 0.274881
Parameters : p0 = -0.867943 p1 = 45.8302 p2 = -13.3151 p3 = 13.8162 p4 = 0.172792 p5 = 0.987388
----------> Iteration 6
FVAL = 58.9284179879 Edm = 9.10495513083e-07 Nfcn = 110
Error matrix change = 0.141986
Parameters : p0 = -0.862363 p1 = 45.8286 p2 = -13.316 p3 = 13.8144 p4 = 0.172388 p5 = 0.987291
----------> Iteration 7
FVAL = 58.9284179879 Edm = 3.19179962054e-05 Nfcn = 150
Error matrix change = 0
Parameters : p0 = -0.862363 p1 = 45.8286 p2 = -13.316 p3 = 13.8144 p4 = 0.172388 p5 = 0.987291
----------> Iteration 8
FVAL = 58.9283860619 Edm = 4.4843686043e-13 Nfcn = 163
Error matrix change = 3.09334e-05
Parameters : p0 = -0.864713 p1 = 45.8434 p2 = -13.3214 p3 = 13.8074 p4 = 0.172309 p5 = 0.987281
----------> Iteration 9
FVAL = 58.9283860619 Edm = 4.4843686043e-13 Nfcn = 163
Error matrix change = 3.09334e-05
Parameters : p0 = -0.864713 p1 = 45.8434 p2 = -13.3214 p3 = 13.8074 p4 = 0.172309 p5 = 0.987281
Minuit2Minimizer : Valid minimum - status = 0
FVAL = 58.9283860619491762
Edm = 4.48436860430199946e-13
Nfcn = 163
p0 = -0.864713 +/- 0.891795
p1 = 45.8434 +/- 2.64223
p2 = -13.3214 +/- 0.976969
p3 = 13.8074 +/- 2.1776
p4 = 0.172309 +/- 0.0358282
p5 = 0.987281 +/- 0.0112683
Covariance Matrix:
p0 p1 p2 p3 p4 p5
p0 0.7953 -1.2054 0.34842 -0.15946 -0.0037284 0.00042265
p1 -1.2054 6.9814 -2.5255 -3.0272 -0.037043 -0.0018129
p2 0.34842 -2.5255 0.95447 1.1587 0.014337 0.00058758
p3 -0.15946 -3.0272 1.1587 4.7419 0.05537 0.0017371
p4 -0.0037284 -0.037043 0.014337 0.05537 0.0012837 2.8395e-05
p5 0.00042265 -0.0018129 0.00058758 0.0017371 2.8395e-05 0.00012697
Correlation Matrix:
p0 p1 p2 p3 p4 p5
p0 1 -0.51155 0.39991 -0.08211 -0.11669 0.042058
p1 -0.51155 1 -0.97834 -0.52614 -0.3913 -0.060891
p2 0.39991 -0.97834 1 0.54464 0.40958 0.053373
p3 -0.08211 -0.52614 0.54464 1 0.7097 0.070794
p4 -0.11669 -0.3913 0.40958 0.7097 1 0.070334
p5 0.042058 -0.060891 0.053373 0.070794 0.070334 1
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 58.9284
NDf = 54
Edm = 4.48437e-13
NCalls = 163
p0 = -0.864713 +/- 0.891795
p1 = 45.8434 +/- 2.64223
p2 = -13.3214 +/- 0.976969
p3 = 13.8074 +/- 2.1776
p4 = 0.172309 +/- 0.0358282
p5 = 0.987281 +/- 0.0112683
#include "TH1.h"
#include "TMath.h"
#include "TF1.h"
#include "TLegend.h"
#include "TCanvas.h"
// Quadratic background function
double background(double *x, double *par) {
return par[0] + par[1]*x[0] + par[2]*x[0]*x[0];
}
// Lorenzian Peak function
double lorentzianPeak(double *x, double *par) {
return (0.5*par[0]*par[1]/TMath::Pi()) /
TMath::Max( 1.e-10,(x[0]-par[2])*(x[0]-par[2])
+ .25*par[1]*par[1]);
}
// Sum of background and peak function
double fitFunction(double *x, double *par) {
return background(x,par) + lorentzianPeak(x,&par[3]);
}
void FittingDemo() {
//Bevington Exercise by Peter Malzacher, modified by Rene Brun
const int nBins = 60;
double data[nBins] = { 6, 1,10,12, 6,13,23,22,15,21,
23,26,36,25,27,35,40,44,66,81,
75,57,48,45,46,41,35,36,53,32,
40,37,38,31,36,44,42,37,32,32,
43,44,35,33,33,39,29,41,32,44,
26,39,29,35,32,21,21,15,25,15};
TCanvas *c1 = new TCanvas("c1","Fitting Demo",10,10,700,500);
c1->SetFillColor(33);
c1->SetFrameFillColor(41);
c1->SetGrid();
TH1F *histo = new TH1F("histo",
"Lorentzian Peak on Quadratic Background",60,0,3);
histo->SetMarkerStyle(21);
histo->SetMarkerSize(0.8);
histo->SetStats(false);
for(int i=0; i < nBins; i++) histo->SetBinContent(i+1,data[i]);
// create a TF1 with the range from 0 to 3 and 6 parameters
TF1 *fitFcn = new TF1("fitFcn",fitFunction,0,3,6);
fitFcn->SetNpx(500);
fitFcn->SetLineWidth(4);
// first try without starting values for the parameters
// This defaults to 1 for each param.
// this results in an ok fit for the polynomial function
// however the non-linear part (lorenzian) does not
// respond well.
fitFcn->SetParameters(1,1,1,1,1,1);
histo->Fit("fitFcn","0");
// second try: set start values for some parameters
fitFcn->SetParameter(4,0.2); // width
fitFcn->SetParameter(5,1); // peak
histo->Fit("fitFcn","V+","ep");
// improve the picture:
TF1 *backFcn = new TF1("backFcn",background,0,3,3);
backFcn->SetLineColor(kRed);
TF1 *signalFcn = new TF1("signalFcn",lorentzianPeak,0,3,3);
signalFcn->SetLineColor(kBlue);
signalFcn->SetNpx(500);
double par[6];
// writes the fit results into the par array
fitFcn->GetParameters(par);
backFcn->SetParameters(par);
backFcn->Draw("same");
signalFcn->SetParameters(&par[3]);
signalFcn->Draw("same");
// draw the legend
TLegend *legend=new TLegend(0.6,0.65,0.88,0.85);
legend->SetTextFont(72);
legend->SetTextSize(0.04);
legend->AddEntry(histo,"Data","lpe");
legend->AddEntry(backFcn,"Background fit","l");
legend->AddEntry(signalFcn,"Signal fit","l");
legend->AddEntry(fitFcn,"Global Fit","l");
legend->Draw();
}
@ kRed
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
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
virtual Double_t * GetParameters() const
Definition TF1.h:546
virtual void SetParameters(const Double_t *params)
Definition TF1.h:670
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:660
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
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 void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9190
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8958
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
constexpr Double_t Pi()
Definition TMath.h:37
Author
Rene Brun

Definition in file FittingDemo.C.