-
Notifications
You must be signed in to change notification settings - Fork 18
/
rf305_condcorrprod.py
93 lines (73 loc) · 2.78 KB
/
rf305_condcorrprod.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
# /
#
# 'MULTIDIMENSIONAL MODELS' ROOT.RooFit tutorial macro #305
#
# Multi-dimensional p.d.f.s with conditional p.d.fs in product
#
# pdf = gauss(x,f(y),sx | y ) * gauss(y,ms,sx) with f(y) = a0 + a1*y
#
#
# 07/2008 - Wouter Verkerke
#
# /
import ROOT
def rf305_condcorrprod():
# C r e a t e c o n d i t i o n a l p d f g x ( x | y )
# -----------------------------------------------------------
# Create observables
x = ROOT.RooRealVar("x", "x", -5, 5)
y = ROOT.RooRealVar("y", "y", -5, 5)
# Create function f(y) = a0 + a1*y
a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
fy = ROOT.RooPolyVar("fy", "fy", y, ROOT.RooArgList(a0, a1))
# Create gaussx(x,f(y),sx)
sigmax = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
gaussx = ROOT.RooGaussian(
"gaussx", "Gaussian in x with shifting mean in y", x, fy, sigmax)
# C r e a t e p d f g y ( y )
# -----------------------------------------------------------
# Create gaussy(y,0,5)
gaussy = ROOT.RooGaussian(
"gaussy", "Gaussian in y", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(3))
# C r e a t e p r o d u c t g x ( x | y ) * g y ( y )
# -------------------------------------------------------
# Create gaussx(x,sx|y) * gaussy(y)
model = ROOT.RooProdPdf(
"model", "gaussx(x|y)*gaussy(y)", ROOT.RooArgSet(gaussy), ROOT.RooFit.Conditional(ROOT.RooArgSet(gaussx), ROOT.RooArgSet(x)))
# S a m p l e , i t a n d p l o t p r o d u c t p d f
# ---------------------------------------------------------------
# Generate 1000 events in x and y from model
data = model.generate(ROOT.RooArgSet(x, y), 10000)
# Plot x distribution of data and projection of model x = Int(dy)
# model(x,y)
xframe = x.frame()
data.plotOn(xframe)
model.plotOn(xframe)
# Plot x distribution of data and projection of model y = Int(dx)
# model(x,y)
yframe = y.frame()
data.plotOn(yframe)
model.plotOn(yframe)
# Make two-dimensional plot in x vs y
hh_model = model.createHistogram("hh_model", x, ROOT.RooFit.Binning(
50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
hh_model.SetLineColor(ROOT.kBlue)
# Make canvas and draw ROOT.RooPlots
c = ROOT.TCanvas("rf305_condcorrprod", "rf05_condcorrprod", 1200, 400)
c.Divide(3)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
xframe.GetYaxis().SetTitleOffset(1.6)
xframe.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
yframe.GetYaxis().SetTitleOffset(1.6)
yframe.Draw()
c.cd(3)
ROOT.gPad.SetLeftMargin(0.20)
hh_model.GetZaxis().SetTitleOffset(2.5)
hh_model.Draw("surf")
c.SaveAs("rf305_condcorrprod.png")
if __name__ == "__main__":
rf305_condcorrprod()