Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf701_efficiencyfit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Special pdf's: unbinned maximum likelihood fit of an efficiency eff(x) function to a
5## dataset D(x,cut), cut is a category encoding a selection, which the efficiency as function
6## of x should be described by eff(x)
7##
8## \macro_code
9##
10## \date February 2018
11## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13import ROOT
14
15
16# Construct efficiency function e(x)
17# -------------------------------------------------------------------
18
19# Declare variables x,mean, with associated name, title, value and allowed
20# range
21x = ROOT.RooRealVar("x", "x", -10, 10)
22
23# Efficiency function eff(x;a,b)
24a = ROOT.RooRealVar("a", "a", 0.4, 0, 1)
25b = ROOT.RooRealVar("b", "b", 5)
26c = ROOT.RooRealVar("c", "c", -1, -10, 10)
27effFunc = ROOT.RooFormulaVar("effFunc", "(1-a)+a*cos((x-c)/b)", [a, b, c, x])
28
29# Construct conditional efficiency pdf E(cut|x)
30# ------------------------------------------------------------------------------------------
31
32# Acceptance state cut (1 or 0)
33cut = ROOT.RooCategory("cut", "cutr", {"accept": 1, "reject": 0})
34
35# Construct efficiency pdf eff(cut|x)
36effPdf = ROOT.RooEfficiency("effPdf", "effPdf", effFunc, cut, "accept")
37
38# Generate data (x, cut) from a toy model
39# -----------------------------------------------------------------------------
40
41# Construct global shape pdf shape(x) and product model(x,cut) = eff(cut|x)*shape(x)
42# (These are _only_ needed to generate some toy MC here to be used later)
43shapePdf = ROOT.RooPolynomial("shapePdf", "shapePdf", x, [-0.095])
44model = ROOT.RooProdPdf("model", "model", {shapePdf}, Conditional=({effPdf}, {cut}))
45
46# Generate some toy data from model
47data = model.generate({x, cut}, 10000)
48
49# Fit conditional efficiency pdf to data
50# --------------------------------------------------------------------------
51
52# Fit conditional efficiency pdf to data
53effPdf.fitTo(data, ConditionalObservables={x})
54
55# Plot fitted, data efficiency
56# --------------------------------------------------------
57
58# Plot distribution of all events and accepted fraction of events on frame
59frame1 = x.frame(Bins=20, Title="Data (all, accepted)")
60data.plotOn(frame1)
61data.plotOn(frame1, Cut="cut==cut::accept", MarkerColor="r", LineColor="r")
62
63# Plot accept/reject efficiency on data overlay fitted efficiency curve
64frame2 = x.frame(Bins=20, Title="Fitted efficiency")
65data.plotOn(frame2, Efficiency=cut) # needs ROOT version >= 5.21
66effFunc.plotOn(frame2, LineColor="r")
67
68# Draw all frames on a canvas
69ca = ROOT.TCanvas("rf701_efficiency", "rf701_efficiency", 800, 400)
70ca.Divide(2)
71ca.cd(1)
72ROOT.gPad.SetLeftMargin(0.15)
73frame1.GetYaxis().SetTitleOffset(1.6)
74frame1.Draw()
75ca.cd(2)
76ROOT.gPad.SetLeftMargin(0.15)
77frame2.GetYaxis().SetTitleOffset(1.4)
78frame2.Draw()
79
80ca.SaveAs("rf701_efficiencyfit.png")