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