-
Notifications
You must be signed in to change notification settings - Fork 18
/
Copy pathrf309_ndimplot.py
116 lines (92 loc) · 4.12 KB
/
rf309_ndimplot.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
107
108
109
110
111
112
113
114
115
116
#####################################
#
# 'MULTIDIMENSIONAL MODELS' ROOT.RooFit tutorial macro #308
#
# Making 2/3 dimensional plots of p.d.f.s and datasets
#
#
#
# 07/2008 - Wouter Verkerke
#
# /
import ROOT
def rf309_ndimplot():
# C r e a t e 2 D m o d e l a n d d a t a s e t
# -----------------------------------------------------
# Create observables
x = ROOT.RooRealVar("x", "x", -5, 5)
y = ROOT.RooRealVar("y", "y", -5, 5)
# Create parameters
a0 = ROOT.RooRealVar("a0", "a0", -3.5, -5, 5)
a1 = ROOT.RooRealVar("a1", "a1", -1.5, -1, 1)
sigma = ROOT.RooRealVar("sigma", "width of gaussian", 1.5)
# Create interpreted function f(y) = a0 - a1*sqrt(10*abs(y))
fy = ROOT.RooFormulaVar("fy", "a0-a1*sqrt(10*abs(y))",
ROOT.RooArgList(y, a0, a1))
# Create gauss(x,f(y),s)
model = ROOT.RooGaussian(
"model", "Gaussian with shifting mean", x, fy, sigma)
# Sample dataset from gauss(x,y)
data = model.generate(ROOT.RooArgSet(x, y), 10000)
# M a k e 2 D p l o t s o f d a t a a n d m o d e l
# -------------------------------------------------------------
# Create and fill ROOT 2D histogram (20x20 bins) with contents of dataset
# hh_data = data.createHistogram("hh_data",x, ROOT.RooFit.Binning(20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
# hh_data = data.createHistogram("x,y", 20, 20) # does not work, see
# https://root.cern.ch/phpBB3/viewtopic.php?t=16648
hh_data = ROOT.RooAbsData.createHistogram(data, "x,y", x, ROOT.RooFit.Binning(
20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)))
# Create and fill ROOT 2D histogram (50x50 bins) with sampling of pdf
# hh_pdf = model.createHistogram("hh_model",x, ROOT.RooFit.Binning(50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
hh_pdf = model.createHistogram("x,y", 50, 50)
hh_pdf.SetLineColor(ROOT.kBlue)
# C r e a t e 3 D m o d e l a n d d a t a s e t
# -----------------------------------------------------
# Create observables
z = ROOT.RooRealVar("z", "z", -5, 5)
gz = ROOT.RooGaussian(
"gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(2))
model3 = ROOT.RooProdPdf("model3", "model3", ROOT.RooArgList(model, gz))
data3 = model3.generate(ROOT.RooArgSet(x, y, z), 10000)
# M a k e 3 D p l o t s o f d a t a a n d m o d e l
# -------------------------------------------------------------
# Create and fill ROOT 2D histogram (8x8x8 bins) with contents of dataset
# hh_data3 = data3.createHistogram("hh_data3", x, ROOT.RooFit.Binning(8), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(8)), ROOT.RooFit.ZVar(z, ROOT.RooFit.Binning(8)))
hh_data3 = ROOT.RooAbsData.createHistogram(data3, "hh_data3", x, ROOT.RooFit.Binning(
8), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(8)), ROOT.RooFit.ZVar(z, ROOT.RooFit.Binning(8)))
# Create and fill ROOT 2D histogram (20x20x20 bins) with sampling of pdf
hh_pdf3 = model3.createHistogram("hh_model3", x, ROOT.RooFit.Binning(
20), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(20)), ROOT.RooFit.ZVar(z, ROOT.RooFit.Binning(20)))
hh_pdf3.SetFillColor(ROOT.kBlue)
c1 = ROOT.TCanvas("rf309_2dimplot", "rf309_2dimplot", 800, 800)
c1.Divide(2, 2)
c1.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
hh_data.GetZaxis().SetTitleOffset(1.4)
hh_data.Draw("lego")
c1.cd(2)
ROOT.gPad.SetLeftMargin(0.20)
hh_pdf.GetZaxis().SetTitleOffset(2.5)
hh_pdf.Draw("surf")
c1.cd(3)
ROOT.gPad.SetLeftMargin(0.15)
hh_data.GetZaxis().SetTitleOffset(1.4)
hh_data.Draw("box")
c1.cd(4)
ROOT.gPad.SetLeftMargin(0.15)
hh_pdf.GetZaxis().SetTitleOffset(2.5)
hh_pdf.Draw("cont3")
c1.SaveAs("rf309_2dimplot.png")
c2 = ROOT.TCanvas("rf309_3dimplot", "rf309_3dimplot", 800, 400)
c2.Divide(2)
c2.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
hh_data3.GetZaxis().SetTitleOffset(1.4)
hh_data3.Draw("lego")
c2.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
hh_pdf3.GetZaxis().SetTitleOffset(1.4)
hh_pdf3.Draw("iso")
c2.SaveAs("rf309_3dimplot.png")
if __name__ == "__main__":
rf309_ndimplot()