ROOT
Version v6.34
master
v6.32
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
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_image
7
## \macro_code
8
## \macro_output
9
##
10
## \date February 2018
11
## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13
import
ROOT
14
15
16
# Define observables and decay pdf
17
# ---------------------------------------------------------------
18
19
# Declare observables
20
t =
ROOT.RooRealVar
(
"t"
,
"t"
, 0, 5)
21
22
# Make pdf
23
tau =
ROOT.RooRealVar
(
"tau"
,
"tau"
, -1.54, -4, -0.1)
24
model =
ROOT.RooExponential
(
"model"
,
"model"
, t, tau)
25
26
# Define efficiency function
27
# ---------------------------------------------------
28
29
# Use error function to simulate turn-on slope
30
eff =
ROOT.RooFormulaVar
(
"eff"
,
"0.5*(TMath::Erf((t-1)/0.5)+1)"
, [t])
31
32
# Define decay pdf with efficiency
33
# ---------------------------------------------------------------
34
35
# Multiply pdf(t) with efficiency in t
36
modelEff =
ROOT.RooEffProd
(
"modelEff"
,
"model with efficiency"
, model, eff)
37
38
# Plot efficiency, pdf
39
# ----------------------------------------
40
41
frame1 =
t.frame
(Title=
"Efficiency"
)
42
eff.plotOn
(frame1, LineColor=
"r"
)
43
44
frame2 =
t.frame
(Title=
"Pdf with and without efficiency"
)
45
46
model.plotOn
(frame2, LineStyle=
"--"
)
47
modelEff.plotOn
(frame2)
48
49
# Generate toy data, fit model eff to data
50
# ------------------------------------------------------------------------------
51
52
# Generate events. If the input pdf has an internal generator, internal generator
53
# is used and an accept/reject sampling on the efficiency is applied.
54
data =
modelEff.generate
({t}, 10000)
55
56
# Fit pdf. The normalization integral is calculated numerically.
57
modelEff.fitTo
(data, PrintLevel=-1)
58
59
# Plot generated data and overlay fitted pdf
60
frame3 =
t.frame
(Title=
"Fitted pdf with efficiency"
)
61
data.plotOn
(frame3)
62
modelEff.plotOn
(frame3)
63
64
c =
ROOT.TCanvas
(
"rf703_effpdfprod"
,
"rf703_effpdfprod"
, 1200, 400)
65
c.Divide
(3)
66
c.cd
(1)
67
ROOT.gPad.SetLeftMargin
(0.15)
68
frame1.GetYaxis
().SetTitleOffset(1.4)
69
frame1.Draw
()
70
c.cd
(2)
71
ROOT.gPad.SetLeftMargin
(0.15)
72
frame2.GetYaxis
().SetTitleOffset(1.6)
73
frame2.Draw
()
74
c.cd
(3)
75
ROOT.gPad.SetLeftMargin
(0.15)
76
frame3.GetYaxis
().SetTitleOffset(1.6)
77
frame3.Draw
()
78
79
c.SaveAs
(
"rf703_effpdfprod.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
rf703_effpdfprod.py
ROOT tags/6-34-04 - Reference Guide Generated on Fri Mar 21 2025 04:40:18 (GVA Time) using Doxygen 1.10.0