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(
23 "gx", "gx", x, ROOT.RooFit.RooConst(-2), ROOT.RooFit.RooConst(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 = ROOT.RooArgSet(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(ROOT.RooArgSet(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 = ROOT.RooArgSet(x)
50igx_sig = gx.createIntegral(xset, ROOT.RooFit.NormSet(xset), ROOT.RooFit.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(ROOT.RooArgSet(x))
59
60# Plot cdf of gx versus x
61frame = x.frame(ROOT.RooFit.Title("cdf of Gaussian pdf"))
62gx_cdf.plotOn(frame)
63
64# Draw plot on canvas
65c = ROOT.TCanvas("rf110_normintegration",
66 "rf110_normintegration", 600, 600)
67ROOT.gPad.SetLeftMargin(0.15)
68frame.GetYaxis().SetTitleOffset(1.6)
69frame.Draw()
70
71c.SaveAs("rf110_normintegration.png")