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