ROOT
Version 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
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
Loading...
Searching...
No Matches
rf203_ranges.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit_main
3
## \notebook
4
## Addition and convolution: fitting and plotting in sub ranges
5
##
6
## \macro_image
7
## \macro_code
8
## \macro_output
9
##
10
## \date February 2018
11
## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13
import
ROOT
14
15
# Set up model
16
# ---------------------
17
18
# Construct observables x
19
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -10, 10)
20
21
# Construct gaussx(x,mx,1)
22
mx =
ROOT.RooRealVar
(
"mx"
,
"mx"
, 0, -10, 10)
23
gx =
ROOT.RooGaussian
(
"gx"
,
"gx"
, x, mx, 1.0)
24
25
# px = 1 (flat in x)
26
px =
ROOT.RooPolynomial
(
"px"
,
"px"
, x)
27
28
# model = f*gx + (1-f)px
29
f =
ROOT.RooRealVar
(
"f"
,
"f"
, 0.0, 1.0)
30
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [gx, px], [f])
31
32
# Generated 10000 events in (x,y) from pdf model
33
modelData =
model.generate
({x}, 10000)
34
35
# Fit full range
36
# ---------------------------
37
38
# Fit pdf to all data
39
r_full =
model.fitTo
(modelData, Save=
True
, PrintLevel=-1)
40
41
# Fit partial range
42
# ----------------------------------
43
44
# Define "signal" range in x as [-3,3]
45
x.setRange
(
"signal"
, -3, 3)
46
47
# Fit pdf only to data in "signal" range
48
r_sig =
model.fitTo
(modelData, Save=
True
, Range=
"signal"
, PrintLevel=-1)
49
50
# Plot/print results
51
# ---------------------------------------
52
53
# Make plot frame in x and add data and fitted model
54
frame =
x.frame
(Title=
"Fitting a sub range"
)
55
modelData.plotOn
(frame)
56
model.plotOn
(frame, Range=
"Full"
, LineColor=
"r"
, LineStyle=
"--"
)
# Add shape in full ranged dashed
57
model.plotOn
(frame)
# By default only fitted range is shown
58
59
# Print fit results
60
print(
"result of fit on all data "
)
61
r_full.Print
()
62
print(
"result of fit in in signal region (note increased error on signal fraction)"
)
63
r_sig.Print
()
64
65
# Draw frame on canvas
66
c =
ROOT.TCanvas
(
"rf203_ranges"
,
"rf203_ranges"
, 600, 600)
67
ROOT.gPad.SetLeftMargin
(0.15)
68
frame.GetYaxis
().SetTitleOffset(1.4)
69
frame.Draw
()
70
71
c.SaveAs
(
"rf203_ranges.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
rf203_ranges.py
ROOT master - Reference Guide Generated on Wed Apr 2 2025 07:35:44 (GVA Time) using Doxygen 1.10.0