Logo ROOT  
Reference Guide
rf110_normintegration.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## Basic functionality: examples on normalization and integration of p.d.fs, construction
6## of cumulative distribution functions from monodimensional p.d.f.s
7##
8## \macro_code
9##
10## \date February 2018
11## \authors Clemens Lange, Wouter Verkerke (C++ version)
12
13from __future__ import print_function
14import ROOT
15
16# Set up model
17# ---------------------
18
19# Create observables x,y
20x = ROOT.RooRealVar("x", "x", -10, 10)
21
22# Create p.d.f. gaussx(x,-2,3)
23gx = ROOT.RooGaussian(
24 "gx", "gx", x, ROOT.RooFit.RooConst(-2), ROOT.RooFit.RooConst(3))
25
26# Retrieve raw & normalized values of RooFit p.d.f.s
27# --------------------------------------------------------------------------------------------------
28
29# Return 'raw' unnormalized value of gx
30print("gx = ", gx.getVal())
31
32# Return value of gx normalized over x in range [-10,10]
33nset = ROOT.RooArgSet(x)
34print("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]
38igx = gx.createIntegral(ROOT.RooArgSet(x))
39print("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
45x.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 p.d.f. gx_Norm[x] which is in the
49# range named "signal"
50xset = ROOT.RooArgSet(x)
51igx_sig = gx.createIntegral(xset, ROOT.RooFit.NormSet(xset), ROOT.RooFit.Range("signal"))
52print("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'
59gx_cdf = gx.createCdf(ROOT.RooArgSet(x))
60
61# Plot cdf of gx versus x
62frame = x.frame(ROOT.RooFit.Title("c.d.f of Gaussian p.d.f"))
63gx_cdf.plotOn(frame)
64
65# Draw plot on canvas
66c = ROOT.TCanvas("rf110_normintegration",
67 "rf110_normintegration", 600, 600)
68ROOT.gPad.SetLeftMargin(0.15)
69frame.GetYaxis().SetTitleOffset(1.6)
70frame.Draw()
71
72c.SaveAs("rf110_normintegration.png")