Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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("eff", "0.5*(TMath::Erf((t-1)/0.5)+1)", [t])
29
30# Define decay pdf with efficiency
31# ---------------------------------------------------------------
32
33# Multiply pdf(t) with efficiency in t
34modelEff = ROOT.RooEffProd("modelEff", "model with efficiency", model, eff)
35
36# Plot efficiency, pdf
37# ----------------------------------------
38
39frame1 = t.frame(Title="Efficiency")
40eff.plotOn(frame1, LineColor="r")
41
42frame2 = t.frame(Title="Pdf with and without efficiency")
43
44model.plotOn(frame2, LineStyle="--")
45modelEff.plotOn(frame2)
46
47# Generate toy data, fit model eff to data
48# ------------------------------------------------------------------------------
49
50# Generate events. If the input pdf has an internal generator, internal generator
51# is used and an accept/reject sampling on the efficiency is applied.
52data = modelEff.generate({t}, 10000)
53
54# Fit pdf. The normalization integral is calculated numerically.
55modelEff.fitTo(data)
56
57# Plot generated data and overlay fitted pdf
58frame3 = t.frame(Title="Fitted pdf with efficiency")
59data.plotOn(frame3)
60modelEff.plotOn(frame3)
61
62c = ROOT.TCanvas("rf703_effpdfprod", "rf703_effpdfprod", 1200, 400)
63c.Divide(3)
64c.cd(1)
65ROOT.gPad.SetLeftMargin(0.15)
66frame1.GetYaxis().SetTitleOffset(1.4)
67frame1.Draw()
68c.cd(2)
69ROOT.gPad.SetLeftMargin(0.15)
70frame2.GetYaxis().SetTitleOffset(1.6)
71frame2.Draw()
72c.cd(3)
73ROOT.gPad.SetLeftMargin(0.15)
74frame3.GetYaxis().SetTitleOffset(1.6)
75frame3.Draw()
76
77c.SaveAs("rf703_effpdfprod.png")