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