ROOT
Version master
master
v6.34
v6.32
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
v6.04
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
rf109_chi2residpull.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit_main
3
## \notebook
4
## 'BASIC FUNCTIONALITY' RooFit tutorial macro #109
5
## Calculating chi^2 from histograms and curves in ROOT.RooPlots,
6
## making histogram of residual and pull distributions
7
##
8
## \macro_image
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \authors Clemens Lange, Wouter Verkerke (C version)
14
15
import
ROOT
16
17
# Set up model
18
# ---------------------
19
20
# Create observables
21
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -10, 10)
22
23
# Create Gaussian
24
sigma =
ROOT.RooRealVar
(
"sigma"
,
"sigma"
, 3, 0.1, 10)
25
mean =
ROOT.RooRealVar
(
"mean"
,
"mean"
, 0, -10, 10)
26
gauss =
ROOT.RooGaussian
(
"gauss"
,
"gauss"
, x, mean, sigma)
27
28
# Generate a sample of 1000 events with sigma=3
29
data =
gauss.generate
({x}, 10000)
30
31
# Change sigma to 3.15
32
sigma.setVal
(3.15)
33
34
# Plot data and slightly distorted model
35
# ---------------------------------------------------------------------------
36
37
# Overlay projection of gauss with sigma=3.15 on data with sigma=3.0
38
frame1 =
x.frame
(Title=
"Data with distorted Gaussian pdf"
, Bins=40)
39
data.plotOn
(frame1, DataError=
"SumW2"
)
40
gauss.plotOn
(frame1)
41
42
# Calculate chi^2
43
# ------------------------------
44
45
# Show the chi^2 of the curve w.r.t. the histogram
46
# If multiple curves or datasets live in the frame you can specify
47
# the name of the relevant curve and/or dataset in chiSquare()
48
print(
"chi^2 = "
,
frame1.chiSquare
())
49
50
# Show residual and pull dists
51
# -------------------------------------------------------
52
53
# Construct a histogram with the residuals of the data w.r.t. the curve
54
hresid =
frame1.residHist
()
55
56
# Construct a histogram with the pulls of the data w.r.t the curve
57
hpull =
frame1.pullHist
()
58
59
# Create a frame to draw the residual distribution and add the
60
# distribution to the frame
61
frame2 =
x.frame
(Title=
"Residual Distribution"
)
62
frame2.addPlotable
(hresid,
"P"
)
63
64
# Create a frame to draw the pull distribution and add the distribution to
65
# the frame
66
frame3 =
x.frame
(Title=
"Pull Distribution"
)
67
frame3.addPlotable
(hpull,
"P"
)
68
69
c =
ROOT.TCanvas
(
"rf109_chi2residpull"
,
"rf109_chi2residpull"
, 900, 300)
70
c.Divide
(3)
71
c.cd
(1)
72
ROOT.gPad.SetLeftMargin
(0.15)
73
frame1.GetYaxis
().SetTitleOffset(1.6)
74
frame1.Draw
()
75
c.cd
(2)
76
ROOT.gPad.SetLeftMargin
(0.15)
77
frame2.GetYaxis
().SetTitleOffset(1.6)
78
frame2.Draw
()
79
c.cd
(3)
80
ROOT.gPad.SetLeftMargin
(0.15)
81
frame3.GetYaxis
().SetTitleOffset(1.6)
82
frame3.Draw
()
83
84
c.SaveAs
(
"rf109_chi2residpull.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
roofit
rf109_chi2residpull.py
ROOT master - Reference Guide Generated on Fri Mar 7 2025 21:57:58 (GVA Time) using Doxygen 1.10.0