-
Notifications
You must be signed in to change notification settings - Fork 18
/
rf312_multirangefit.py
106 lines (80 loc) · 3.28 KB
/
rf312_multirangefit.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
#####################################
#
# 'MULTIDIMENSIONAL MODELS' ROOT.RooFit tutorial macro #312
#
# Performing fits in multiple (disjoint) ranges in one or more dimensions
#
#
#
# 07/2008 - Wouter Verkerke
#
# /
import ROOT
def rf312_multirangefit():
# C r e a t e 2 D p d f a n d d a t a
# -------------------------------------------
# Define observables x,y
x = ROOT.RooRealVar("x", "x", -10, 10)
y = ROOT.RooRealVar("y", "y", -10, 10)
# Construct the signal pdf gauss(x)*gauss(y)
mx = ROOT.RooRealVar("mx", "mx", 1, -10, 10)
my = ROOT.RooRealVar("my", "my", 1, -10, 10)
gx = ROOT.RooGaussian("gx", "gx", x, mx, ROOT.RooFit.RooConst(1))
gy = ROOT.RooGaussian("gy", "gy", y, my, ROOT.RooFit.RooConst(1))
sig = ROOT.RooProdPdf("sig", "sig", gx, gy)
# Construct the background pdf (flat in x,y)
px = ROOT.RooPolynomial("px", "px", x)
py = ROOT.RooPolynomial("py", "py", y)
bkg = ROOT.RooProdPdf("bkg", "bkg", px, py)
# Construct the composite model sig+bkg
f = ROOT.RooRealVar("f", "f", 0., 1.)
model = ROOT.RooAddPdf("model", "model", ROOT.RooArgList(sig, bkg), ROOT.RooArgList(f))
# Sample 10000 events in (x,y) from the model
modelData = model.generate(ROOT.RooArgSet(x, y), 10000)
# D e f i n e s i g n a l a n d s i d e b a n d r e g i o n s
# -------------------------------------------------------------------
# Construct the SideBand1,SideBand2, regions
#
# |
# +-------------+-----------+
# | | |
# | Side | Sig |
# | Band1 | nal |
# | | |
# --+-------------+-----------+--
# | |
# | Side |
# | Band2 |
# | |
# +-------------+-----------+
# |
x.setRange("SB1", -10, +10)
y.setRange("SB1", -10, 0)
x.setRange("SB2", -10, 0)
y.setRange("SB2", 0, +10)
x.setRange("SIG", 0, +10)
y.setRange("SIG", 0, +10)
x.setRange("FULL", -10, +10)
y.setRange("FULL", -10, +10)
# P e r f o r m f i t s i n i n d i v i d u a l s i d e b a n d r e g i o n s
# -------------------------------------------------------------------------------------
# Perform fit in SideBand1 region (ROOT.RooAddPdf coefficients will be
# interpreted in full range)
r_sb1 = model.fitTo(modelData, ROOT.RooFit.Range(
"SB1"), ROOT.RooFit.Save())
# Perform fit in SideBand2 region (ROOT.RooAddPdf coefficients will be
# interpreted in full range)
r_sb2 = model.fitTo(modelData, ROOT.RooFit.Range(
"SB2"), ROOT.RooFit.Save())
# P e r f o r m f i t s i n j o i n t s i d e b a n d r e g i o n s
# -----------------------------------------------------------------------------
# Now perform fit to joint 'L-shaped' sideband region 'SB1|SB2'
# (ROOT.RooAddPdf coefficients will be interpreted in full range)
r_sb12 = model.fitTo(modelData, ROOT.RooFit.Range(
"SB1,SB2"), ROOT.RooFit.Save())
# Print results for comparison
r_sb1.Print()
r_sb2.Print()
r_sb12.Print()
if __name__ == "__main__":
rf312_multirangefit()