ROOT
Version v6.32
master
v6.34
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
rf103_interprfuncs.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## Basic functionality: interpreted functions and pdfs
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
# Generic interpreted pdf
16
# ------------------------------
17
18
# Declare observable x
19
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -20, 20)
20
21
# Construct generic pdf from interpreted expression
22
# ------------------------------------------------------
23
24
# ROOT.To construct a proper pdf, the formula expression is explicitly normalized internally by dividing
25
# it by a numeric integral of the expression over x in the range [-20,20]
26
#
27
alpha =
ROOT.RooRealVar
(
"alpha"
,
"alpha"
, 5, 0.1, 10)
28
genpdf =
ROOT.RooGenericPdf
(
"genpdf"
,
"genpdf"
,
"(1+0.1*abs(x)+sin(sqrt(abs(x*alpha+0.1))))"
, [x, alpha])
29
30
# Sample, fit and plot generic pdf
31
# ---------------------------------------------------------------
32
33
# Generate a toy dataset from the interpreted pdf
34
data =
genpdf.generate
({x}, 10000)
35
36
# Fit the interpreted pdf to the generated data
37
genpdf.fitTo
(data, PrintLevel=-1)
38
39
# Make a plot of the data and the pdf overlaid
40
xframe =
x.frame
(Title=
"Interpreted expression pdf"
)
41
data.plotOn
(xframe)
42
genpdf.plotOn
(xframe)
43
44
# Standard pdf adjust with interpreted helper function
45
# ------------------------------------------------------------------------------------------------------------
46
# Make a gauss(x,sqrt(mean2),sigma) from a standard ROOT.RooGaussian #
47
#
48
# Construct standard pdf with formula replacing parameter
49
# ------------------------------------------------------------------------------------------------------------
50
51
# Construct parameter mean2 and sigma
52
mean2 =
ROOT.RooRealVar
(
"mean2"
,
"mean^2"
, 10, 0, 200)
53
sigma =
ROOT.RooRealVar
(
"sigma"
,
"sigma"
, 3, 0.1, 10)
54
55
# Construct interpreted function mean = sqrt(mean^2)
56
mean =
ROOT.RooFormulaVar
(
"mean"
,
"mean"
,
"sqrt(mean2)"
, [mean2])
57
58
# Construct a gaussian g2(x,sqrt(mean2),sigma)
59
g2 =
ROOT.RooGaussian
(
"g2"
,
"h2"
, x, mean, sigma)
60
61
# Generate toy data
62
# ---------------------------------
63
64
# Construct a separate gaussian g1(x,10,3) to generate a toy Gaussian
65
# dataset with mean 10 and width 3
66
g1 =
ROOT.RooGaussian
(
"g1"
,
"g1"
, x, 10, 3)
67
data2 =
g1.generate
({x}, 1000)
68
69
# Fit and plot tailored standard pdf
70
# -------------------------------------------------------------------
71
72
# Fit g2 to data from g1
73
r =
g2.fitTo
(data2, Save=
True
, PrintLevel=-1)
# ROOT.RooFitResult
74
r.Print
()
75
76
# Plot data on frame and overlay projection of g2
77
xframe2 =
x.frame
(Title=
"Tailored Gaussian pdf"
)
78
data2.plotOn
(xframe2)
79
g2.plotOn
(xframe2)
80
81
# Draw all frames on a canvas
82
c =
ROOT.TCanvas
(
"rf103_interprfuncs"
,
"rf103_interprfuncs"
, 800, 400)
83
c.Divide
(2)
84
c.cd
(1)
85
ROOT.gPad.SetLeftMargin
(0.15)
86
xframe.GetYaxis
().SetTitleOffset(1.4)
87
xframe.Draw
()
88
c.cd
(2)
89
ROOT.gPad.SetLeftMargin
(0.15)
90
xframe2.GetYaxis
().SetTitleOffset(1.4)
91
xframe2.Draw
()
92
93
c.SaveAs
(
"rf103_interprfuncs.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
rf103_interprfuncs.py
ROOT v6-32 - Reference Guide Generated on Wed Apr 2 2025 08:24:49 (GVA Time) using Doxygen 1.10.0