Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf316_llratioplot.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: using the likelihood ratio techique to construct a signal
5## enhanced one-dimensional projection of a multi-dimensional pdf
6##
7## \macro_code
8##
9## \date February 2018
10## \authors Clemens Lange, Wouter Verkerke (C++ version)
11
12import ROOT
13
14
15# Create 3D pdf and data
16# -------------------------------------------
17
18# Create observables
19x = ROOT.RooRealVar("x", "x", -5, 5)
20y = ROOT.RooRealVar("y", "y", -5, 5)
21z = ROOT.RooRealVar("z", "z", -5, 5)
22
23# Create signal pdf gauss(x)*gauss(y)*gauss(z)
24gx = ROOT.RooGaussian("gx", "gx", x, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
25gy = ROOT.RooGaussian("gy", "gy", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
26gz = ROOT.RooGaussian("gz", "gz", z, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(1))
27sig = ROOT.RooProdPdf("sig", "sig", [gx, gy, gz])
28
29# Create background pdf poly(x)*poly(y)*poly(z)
30px = ROOT.RooPolynomial("px", "px", x, [-0.1, 0.004])
31py = ROOT.RooPolynomial("py", "py", y, [0.1, -0.004])
32pz = ROOT.RooPolynomial("pz", "pz", z)
33bkg = ROOT.RooProdPdf("bkg", "bkg", [px, py, pz])
34
35# Create composite pdf sig+bkg
36fsig = ROOT.RooRealVar("fsig", "signal fraction", 0.1, 0.0, 1.0)
37model = ROOT.RooAddPdf("model", "model", [sig, bkg], [fsig])
38
39data = model.generate({x, y, z}, 20000)
40
41# Project pdf and data on x
42# -------------------------------------------------
43
44# Make plain projection of data and pdf on x observable
45frame = x.frame(Title="Projection of 3D data and pdf on X", Bins=40)
46data.plotOn(frame)
47model.plotOn(frame)
48
49# Define projected signal likelihood ratio
50# ----------------------------------------------------------------------------------
51
52# Calculate projection of signal and total likelihood on (y,z) observables
53# i.e. integrate signal and composite model over x
54sigyz = sig.createProjection({x})
55totyz = model.createProjection({x})
56
57# Construct the log of the signal / signal+background probability
58llratio_func = ROOT.RooFormulaVar("llratio", "log10(@0)-log10(@1)", [sigyz, totyz])
59
60# Plot data with a LL ratio cut
61# -------------------------------------------------------
62
63# Calculate the llratio value for each event in the dataset
64data.addColumn(llratio_func)
65
66# Extract the subset of data with large signal likelihood
67dataSel = data.reduce(Cut="llratio>0.7")
68
69# Make plot frame
70frame2 = x.frame(Title="Same projection on X with LLratio(y,z)>0.7", Bins=40)
71
72# Plot select data on frame
73dataSel.plotOn(frame2)
74
75# Make MC projection of pdf with same LL ratio cut
76# ---------------------------------------------------------------------------------------------
77
78# Generate large number of events for MC integration of pdf projection
79mcprojData = model.generate({x, y, z}, 10000)
80
81# Calculate LL ratio for each generated event and select MC events with
82# llratio)0.7
83mcprojData.addColumn(llratio_func)
84mcprojDataSel = mcprojData.reduce(Cut="llratio>0.7")
85
86# Project model on x, projected observables (y,z) with Monte Carlo technique
87# on set of events with the same llratio cut as was applied to data
88model.plotOn(frame2, ProjWData=mcprojDataSel)
89
90c = ROOT.TCanvas("rf316_llratioplot", "rf316_llratioplot", 800, 400)
91c.Divide(2)
92c.cd(1)
93ROOT.gPad.SetLeftMargin(0.15)
94frame.GetYaxis().SetTitleOffset(1.4)
95frame.Draw()
96c.cd(2)
97ROOT.gPad.SetLeftMargin(0.15)
98frame2.GetYaxis().SetTitleOffset(1.4)
99frame2.Draw()
100c.SaveAs("rf316_llratioplot.png")