Logo ROOT  
Reference Guide
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 
13 import 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
21 x = ROOT.RooRealVar("x", "x", -10, 10)
22 
23 # Efficiency function eff(x;a,b)
24 a = ROOT.RooRealVar("a", "a", 0.4, 0, 1)
25 b = ROOT.RooRealVar("b", "b", 5)
26 c = ROOT.RooRealVar("c", "c", -1, -10, 10)
27 effFunc = ROOT.RooFormulaVar(
28  "effFunc", "(1-a)+a*cos((x-c)/b)", ROOT.RooArgList(a, b, c, x))
29 
30 # Construct conditional efficiency pdf E(cut|x)
31 # ------------------------------------------------------------------------------------------
32 
33 # Acceptance state cut (1 or 0)
34 cut = ROOT.RooCategory("cut", "cutr")
35 cut.defineType("accept", 1)
36 cut.defineType("reject", 0)
37 
38 # Construct efficiency pdf eff(cut|x)
39 effPdf = ROOT.RooEfficiency("effPdf", "effPdf", effFunc, cut, "accept")
40 
41 # Generate data (x, cut) from a toy model
42 # -----------------------------------------------------------------------------
43 
44 # Construct global shape pdf shape(x) and product model(x,cut) = eff(cut|x)*shape(x)
45 # (These are _only_ needed to generate some toy MC here to be used later)
46 shapePdf = ROOT.RooPolynomial(
47  "shapePdf", "shapePdf", x, ROOT.RooArgList(ROOT.RooFit.RooConst(-0.095)))
48 model = ROOT.RooProdPdf(
49  "model",
50  "model",
51  ROOT.RooArgSet(shapePdf),
52  ROOT.RooFit.Conditional(
53  ROOT.RooArgSet(effPdf),
54  ROOT.RooArgSet(cut)))
55 
56 # Generate some toy data from model
57 data = model.generate(ROOT.RooArgSet(x, cut), 10000)
58 
59 # Fit conditional efficiency pdf to data
60 # --------------------------------------------------------------------------
61 
62 # Fit conditional efficiency pdf to data
63 effPdf.fitTo(data, ROOT.RooFit.ConditionalObservables(ROOT.RooArgSet(x)))
64 
65 # Plot fitted, data efficiency
66 # --------------------------------------------------------
67 
68 # Plot distribution of all events and accepted fraction of events on frame
69 frame1 = x.frame(ROOT.RooFit.Bins(
70  20), ROOT.RooFit.Title("Data (all, accepted)"))
71 data.plotOn(frame1)
72 data.plotOn(
73  frame1,
74  ROOT.RooFit.Cut("cut==cut::accept"),
75  ROOT.RooFit.MarkerColor(
76  ROOT.kRed),
77  ROOT.RooFit.LineColor(
78  ROOT.kRed))
79 
80 # Plot accept/reject efficiency on data overlay fitted efficiency curve
81 frame2 = x.frame(ROOT.RooFit.Bins(
82  20), ROOT.RooFit.Title("Fitted efficiency"))
83 data.plotOn(frame2, ROOT.RooFit.Efficiency(cut)) # needs ROOT version >= 5.21
84 effFunc.plotOn(frame2, ROOT.RooFit.LineColor(ROOT.kRed))
85 
86 # Draw all frames on a canvas
87 ca = ROOT.TCanvas("rf701_efficiency", "rf701_efficiency", 800, 400)
88 ca.Divide(2)
89 ca.cd(1)
90 ROOT.gPad.SetLeftMargin(0.15)
91 frame1.GetYaxis().SetTitleOffset(1.6)
92 frame1.Draw()
93 ca.cd(2)
94 ROOT.gPad.SetLeftMargin(0.15)
95 frame2.GetYaxis().SetTitleOffset(1.4)
96 frame2.Draw()
97 
98 ca.SaveAs("rf701_efficiencyfit.png")