Logo ROOT  
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
14import ROOT
15
16# Set up model
17# ---------------------
18
19# Create observables x,y
20x = ROOT.RooRealVar("x", "x", -10, 10)
21
22# Create pdf gaussx(x,-2,3)
23gx = ROOT.RooGaussian("gx", "gx", x, -2, 3)
24
25# Retrieve raw & normalized values of RooFit pdfs
26# --------------------------------------------------------------------------------------------------
27
28# Return 'raw' unnormalized value of gx
29print("gx = ", gx.getVal())
30
31# Return value of gx normalized over x in range [-10,10]
32nset = {x}
33print("gx_Norm[x] = ", gx.getVal(nset))
34
35# Create object representing integral over gx
36# which is used to calculate gx_Norm[x] == gx / gx_Int[x]
37igx = gx.createIntegral({x})
38print("gx_Int[x] = ", igx.getVal())
39
40# Integrate normalized pdf over subrange
41# ----------------------------------------------------------------------------
42
43# Define a range named "signal" in x from -5,5
44x.setRange("signal", -5, 5)
45
46# Create an integral of gx_Norm[x] over x in range "signal"
47# ROOT.This is the fraction of of pdf gx_Norm[x] which is in the
48# range named "signal"
49xset = {x}
50igx_sig = gx.createIntegral(xset, NormSet=xset, Range="signal")
51print("gx_Int[x|signal]_Norm[x] = ", igx_sig.getVal())
52
53# Construct cumulative distribution function from pdf
54# -----------------------------------------------------------------------------------------------------
55
56# Create the cumulative distribution function of gx
57# i.e. calculate Int[-10,x] gx(x') dx'
58gx_cdf = gx.createCdf({x})
59
60# Plot cdf of gx versus x
61frame = x.frame(Title="cdf of Gaussian pdf")
62gx_cdf.plotOn(frame)
63
64# Draw plot on canvas
65c = ROOT.TCanvas("rf110_normintegration", "rf110_normintegration", 600, 600)
66ROOT.gPad.SetLeftMargin(0.15)
67frame.GetYaxis().SetTitleOffset(1.6)
68frame.Draw()
69
70c.SaveAs("rf110_normintegration.png")