Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf703_effpdfprod.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Special pdf's: using a product of an (acceptance) efficiency and a pdf as pdf
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13
14# Define observables and decay pdf
15# ---------------------------------------------------------------
16
17# Declare observables
18t = ROOT.RooRealVar("t", "t", 0, 5)
19
20# Make pdf
21tau = ROOT.RooRealVar("tau", "tau", -1.54, -4, -0.1)
22model = ROOT.RooExponential("model", "model", t, tau)
23
24# Define efficiency function
25# ---------------------------------------------------
26
27# Use error function to simulate turn-on slope
28eff = ROOT.RooFormulaVar(
29 "eff",
30 "0.5*(TMath::Erf((t-1)/0.5)+1)",
31 ROOT.RooArgList(t))
32
33# Define decay pdf with efficiency
34# ---------------------------------------------------------------
35
36# Multiply pdf(t) with efficiency in t
37modelEff = ROOT.RooEffProd("modelEff", "model with efficiency", model, eff)
38
39# Plot efficiency, pdf
40# ----------------------------------------
41
42frame1 = t.frame(ROOT.RooFit.Title("Efficiency"))
43eff.plotOn(frame1, ROOT.RooFit.LineColor(ROOT.kRed))
44
45frame2 = t.frame(ROOT.RooFit.Title("Pdf with and without efficiency"))
46
47model.plotOn(frame2, ROOT.RooFit.LineStyle(ROOT.kDashed))
48modelEff.plotOn(frame2)
49
50# Generate toy data, fit model eff to data
51# ------------------------------------------------------------------------------
52
53# Generate events. If the input pdf has an internal generator, internal generator
54# is used and an accept/reject sampling on the efficiency is applied.
55data = modelEff.generate(ROOT.RooArgSet(t), 10000)
56
57# Fit pdf. The normalization integral is calculated numerically.
58modelEff.fitTo(data)
59
60# Plot generated data and overlay fitted pdf
61frame3 = t.frame(ROOT.RooFit.Title("Fitted pdf with efficiency"))
62data.plotOn(frame3)
63modelEff.plotOn(frame3)
64
65c = ROOT.TCanvas("rf703_effpdfprod", "rf703_effpdfprod", 1200, 400)
66c.Divide(3)
67c.cd(1)
68ROOT.gPad.SetLeftMargin(0.15)
69frame1.GetYaxis().SetTitleOffset(1.4)
70frame1.Draw()
71c.cd(2)
72ROOT.gPad.SetLeftMargin(0.15)
73frame2.GetYaxis().SetTitleOffset(1.6)
74frame2.Draw()
75c.cd(3)
76ROOT.gPad.SetLeftMargin(0.15)
77frame3.GetYaxis().SetTitleOffset(1.6)
78frame3.Draw()
79
80c.SaveAs("rf703_effpdfprod.png")