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