ROOT
Version v6.32
master
v6.34
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
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
from
__future__
import
print_function
14
import
ROOT
15
16
# Set up model
17
# ---------------------
18
19
# Construct observables x
20
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -10, 10)
21
22
# Construct gaussx(x,mx,1)
23
mx =
ROOT.RooRealVar
(
"mx"
,
"mx"
, 0, -10, 10)
24
gx =
ROOT.RooGaussian
(
"gx"
,
"gx"
, x, mx, 1.0)
25
26
# px = 1 (flat in x)
27
px =
ROOT.RooPolynomial
(
"px"
,
"px"
, x)
28
29
# model = f*gx + (1-f)px
30
f =
ROOT.RooRealVar
(
"f"
,
"f"
, 0.0, 1.0)
31
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [gx, px], [f])
32
33
# Generated 10000 events in (x,y) from pdf model
34
modelData =
model.generate
({x}, 10000)
35
36
# Fit full range
37
# ---------------------------
38
39
# Fit pdf to all data
40
r_full =
model.fitTo
(modelData, Save=
True
, PrintLevel=-1)
41
42
# Fit partial range
43
# ----------------------------------
44
45
# Define "signal" range in x as [-3,3]
46
x.setRange
(
"signal"
, -3, 3)
47
48
# Fit pdf only to data in "signal" range
49
r_sig =
model.fitTo
(modelData, Save=
True
, Range=
"signal"
, PrintLevel=-1)
50
51
# Plot/print results
52
# ---------------------------------------
53
54
# Make plot frame in x and add data and fitted model
55
frame =
x.frame
(Title=
"Fitting a sub range"
)
56
modelData.plotOn
(frame)
57
model.plotOn
(frame, Range=
"Full"
, LineColor=
"r"
, LineStyle=
"--"
)
# Add shape in full ranged dashed
58
model.plotOn
(frame)
# By default only fitted range is shown
59
60
# Print fit results
61
print(
"result of fit on all data "
)
62
r_full.Print
()
63
print(
"result of fit in in signal region (note increased error on signal fraction)"
)
64
r_sig.Print
()
65
66
# Draw frame on canvas
67
c =
ROOT.TCanvas
(
"rf203_ranges"
,
"rf203_ranges"
, 600, 600)
68
ROOT.gPad.SetLeftMargin
(0.15)
69
frame.GetYaxis
().SetTitleOffset(1.4)
70
frame.Draw
()
71
72
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
rf203_ranges.py
ROOT v6-32 - Reference Guide Generated on Wed Mar 26 2025 04:43:22 (GVA Time) using Doxygen 1.10.0