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
rf105_funcbinding.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## 'BASIC FUNCTIONALITY' RooFit tutorial macro #105
5
## Demonstration of binding ROOT Math functions as RooFit functions
6
## and pdfs
7
##
8
## \macro_image
9
## \macro_code
10
## \macro_output
11
##
12
## \date February 2018
13
## \authors Clemens Lange, Wouter Verkerke (C version)
14
15
import
ROOT
16
17
# Bind ROOT TMath::Erf C function
18
# ---------------------------------------------------
19
20
# Bind one-dimensional ROOT.TMath.Erf function as ROOT.RooAbsReal function
21
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -3, 3)
22
erf =
ROOT.RooFit.bindFunction
(
"erf"
,
ROOT.TMath.Erf
, x)
23
24
# Print erf definition
25
erf.Print
()
26
27
# Plot erf on frame
28
frame1 =
x.frame
(Title=
"TMath.Erf bound as ROOT.RooFit function"
)
29
erf.plotOn
(frame1)
30
31
# Bind ROOT::Math::beta_pdf C function
32
# -----------------------------------------------------------------------
33
34
# Bind pdf ROOT.Math.Beta with three variables as ROOT.RooAbsPdf function
35
x2 =
ROOT.RooRealVar
(
"x2"
,
"x2"
, 0, 0.999)
36
a =
ROOT.RooRealVar
(
"a"
,
"a"
, 5, 0, 10)
37
b =
ROOT.RooRealVar
(
"b"
,
"b"
, 2, 0, 10)
38
beta =
ROOT.RooFit.bindPdf
(
"beta"
,
ROOT.Math.beta_pdf
, x2, a, b)
39
40
# Perf beta definition
41
beta.Print
()
42
43
# Generate some events and fit
44
data =
beta.generate
({x2}, 10000)
45
beta.fitTo
(data, PrintLevel=-1)
46
47
# Plot data and pdf on frame
48
frame2 =
x2.frame
(Title=
"ROOT.Math.Beta bound as ROOT.RooFit pdf"
)
49
data.plotOn
(frame2)
50
beta.plotOn
(frame2)
51
52
# Bind ROOT TF1 as RooFit function
53
# ---------------------------------------------------------------
54
55
# Create a ROOT TF1 function
56
fa1 =
ROOT.TF1
(
"fa1"
,
"sin(x)/x"
, 0, 10)
57
58
# Create an observable
59
x3 =
ROOT.RooRealVar
(
"x3"
,
"x3"
, 0.01, 20)
60
61
# Create binding of TF1 object to above observable
62
rfa1 =
ROOT.RooFit.bindFunction
(fa1, x3)
63
64
# Print rfa1 definition
65
rfa1.Print
()
66
67
# Make plot frame in observable, TF1 binding function
68
frame3 =
x3.frame
(Title=
"TF1 bound as ROOT.RooFit function"
)
69
rfa1.plotOn
(frame3)
70
71
c =
ROOT.TCanvas
(
"rf105_funcbinding"
,
"rf105_funcbinding"
, 1200, 400)
72
c.Divide
(3)
73
c.cd
(1)
74
ROOT.gPad.SetLeftMargin
(0.15)
75
frame1.GetYaxis
().SetTitleOffset(1.6)
76
frame1.Draw
()
77
c.cd
(2)
78
ROOT.gPad.SetLeftMargin
(0.15)
79
frame2.GetYaxis
().SetTitleOffset(1.6)
80
frame2.Draw
()
81
c.cd
(3)
82
ROOT.gPad.SetLeftMargin
(0.15)
83
frame3.GetYaxis
().SetTitleOffset(1.6)
84
frame3.Draw
()
85
86
c.SaveAs
(
"rf105_funcbinding.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
rf105_funcbinding.py
ROOT tags/6-34-04 - Reference Guide Generated on Wed Mar 19 2025 04:35:07 (GVA Time) using Doxygen 1.10.0