ROOT
v6-32
Reference Guide
Loading...
Searching...
No Matches
rf110_normintegration.py
Go to the documentation of this file.
1
## \file
2
## \ingroup tutorial_roofit
3
## \notebook
4
## Basic functionality: examples on normalization and integration of pdfs, construction
5
## of cumulative distribution functions from monodimensional pdfs
6
##
7
## \macro_image
8
## \macro_code
9
## \macro_output
10
##
11
## \date February 2018
12
## \authors Clemens Lange, Wouter Verkerke (C++ version)
13
14
from
__future__
import
print_function
15
import
ROOT
16
17
# Set up model
18
# ---------------------
19
20
# Create observables x,y
21
x =
ROOT.RooRealVar
(
"x"
,
"x"
, -10, 10)
22
23
# Create pdf gaussx(x,-2,3)
24
gx =
ROOT.RooGaussian
(
"gx"
,
"gx"
, x, -2, 3)
25
26
# Retrieve raw & normalized values of RooFit pdfs
27
# --------------------------------------------------------------------------------------------------
28
29
# Return 'raw' unnormalized value of gx
30
print(
"gx = "
,
gx.getVal
())
31
32
# Return value of gx normalized over x in range [-10,10]
33
nset = {x}
34
print(
"gx_Norm[x] = "
,
gx.getVal
(nset))
35
36
# Create object representing integral over gx
37
# which is used to calculate gx_Norm[x] == gx / gx_Int[x]
38
igx =
gx.createIntegral
({x})
39
print(
"gx_Int[x] = "
,
igx.getVal
())
40
41
# Integrate normalized pdf over subrange
42
# ----------------------------------------------------------------------------
43
44
# Define a range named "signal" in x from -5,5
45
x.setRange
(
"signal"
, -5, 5)
46
47
# Create an integral of gx_Norm[x] over x in range "signal"
48
# ROOT.This is the fraction of of pdf gx_Norm[x] which is in the
49
# range named "signal"
50
xset = {x}
51
igx_sig =
gx.createIntegral
(xset, NormSet=xset, Range=
"signal"
)
52
print(
"gx_Int[x|signal]_Norm[x] = "
,
igx_sig.getVal
())
53
54
# Construct cumulative distribution function from pdf
55
# -----------------------------------------------------------------------------------------------------
56
57
# Create the cumulative distribution function of gx
58
# i.e. calculate Int[-10,x] gx(x') dx'
59
gx_cdf =
gx.createCdf
({x})
60
61
# Plot cdf of gx versus x
62
frame =
x.frame
(Title=
"cdf of Gaussian pdf"
)
63
gx_cdf.plotOn
(frame)
64
65
# Draw plot on canvas
66
c =
ROOT.TCanvas
(
"rf110_normintegration"
,
"rf110_normintegration"
, 600, 600)
67
ROOT.gPad.SetLeftMargin
(0.15)
68
frame.GetYaxis
().SetTitleOffset(1.6)
69
frame.Draw
()
70
71
c.SaveAs
(
"rf110_normintegration.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
rf110_normintegration.py
ROOT v6-32 - Reference Guide Generated on Thu Mar 13 2025 14:14:36 (GVA Time) using Doxygen 1.10.0