ROOT
v6-32
Reference Guide
Loading...
Searching...
No Matches
rf312_multirangefit.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook -nodraw
4
## Multidimensional models: performing fits in multiple (disjoint) ranges in one or more dimensions
5
##
6
## \macro_code
7
## \macro_output
8
##
9
## \date February 2018
10
## \authors Clemens Lange, Wouter Verkerke (C++ version)
11
12
import
ROOT
13
14
15
# Create 2D pdf and data
16
# -------------------------------------------
17
18
# Define observables x,y
19
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -10, 10)
20
y =
ROOT.RooRealVar
(
"y"
,
"y"
, -10, 10)
21
22
# Construct the signal pdf gauss(x)*gauss(y)
23
mx =
ROOT.RooRealVar
(
"mx"
,
"mx"
, 1, -10, 10)
24
my =
ROOT.RooRealVar
(
"my"
,
"my"
, 1, -10, 10)
25
26
gx =
ROOT.RooGaussian
(
"gx"
,
"gx"
, x, mx, 1.0)
27
gy =
ROOT.RooGaussian
(
"gy"
,
"gy"
, y, my, 1.0)
28
29
sig =
ROOT.RooProdPdf
(
"sig"
,
"sig"
, gx, gy)
30
31
# Construct the background pdf (flat in x,y)
32
px =
ROOT.RooPolynomial
(
"px"
,
"px"
, x)
33
py =
ROOT.RooPolynomial
(
"py"
,
"py"
, y)
34
bkg =
ROOT.RooProdPdf
(
"bkg"
,
"bkg"
, px, py)
35
36
# Construct the composite model sig+bkg
37
f =
ROOT.RooRealVar
(
"f"
,
"f"
, 0.0, 1.0)
38
model =
ROOT.RooAddPdf
(
"model"
,
"model"
, [sig, bkg], [f])
39
40
# Sample 10000 events in (x,y) from the model
41
modelData =
model.generate
({x, y}, 10000)
42
43
# Define signal and sideband regions
44
# -------------------------------------------------------------------
45
46
# Construct the SideBand1,SideBand2, regions
47
#
48
# |
49
# +-------------+-----------+
50
# | | |
51
# | Side | Sig |
52
# | Band1 | nal |
53
# | | |
54
# --+-------------+-----------+--
55
# | |
56
# | Side |
57
# | Band2 |
58
# | |
59
# +-------------+-----------+
60
# |
61
62
x.setRange
(
"SB1"
, -10, +10)
63
y.setRange
(
"SB1"
, -10, 0)
64
65
x.setRange
(
"SB2"
, -10, 0)
66
y.setRange
(
"SB2"
, 0, +10)
67
68
x.setRange
(
"SIG"
, 0, +10)
69
y.setRange
(
"SIG"
, 0, +10)
70
71
x.setRange
(
"FULL"
, -10, +10)
72
y.setRange
(
"FULL"
, -10, +10)
73
74
# Perform fits in individual sideband regions
75
# -------------------------------------------------------------------------------------
76
77
# Perform fit in SideBand1 region (ROOT.RooAddPdf coefficients will be
78
# interpreted in full range)
79
r_sb1 =
model.fitTo
(modelData, Range=
"SB1"
, Save=
True
, PrintLevel=-1)
80
81
# Perform fit in SideBand2 region (ROOT.RooAddPdf coefficients will be
82
# interpreted in full range)
83
r_sb2 =
model.fitTo
(modelData, Range=
"SB2"
, Save=
True
, PrintLevel=-1)
84
85
# Perform fits in joint sideband regions
86
# -----------------------------------------------------------------------------
87
88
# Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
89
# (ROOT.RooAddPdf coefficients will be interpreted in full range)
90
r_sb12 =
model.fitTo
(modelData, Range=
"SB1,SB2"
, Save=
True
, PrintLevel=-1)
91
92
# Print results for comparison
93
r_sb1.Print
()
94
r_sb2.Print
()
95
r_sb12.Print
()
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
rf312_multirangefit.py
ROOT v6-32 - Reference Guide Generated on Thu Mar 13 2025 14:14:36 (GVA Time) using Doxygen 1.10.0