Logo ROOT   6.16/01
Reference Guide
rf110_normintegration.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'BASIC FUNCTIONALITY' RooFit tutorial macro #110
5## Examples on normalization of p.d.f.s,
6## integration of p.d.fs, construction
7## of cumulative distribution functions from p.d.f.s
8## in one dimension
9##
10## \macro_code
11##
12## \date February 2018
13## \author Clemens Lange
14## \author Wouter Verkerke (C version)
15
16from __future__ import print_function
17import ROOT
18
19# Set up model
20# ---------------------
21
22# Create observables x,y
23x = ROOT.RooRealVar("x", "x", -10, 10)
24
25# Create p.d.f. gaussx(x,-2,3)
26gx = ROOT.RooGaussian(
27 "gx", "gx", x, ROOT.RooFit.RooConst(-2), ROOT.RooFit.RooConst(3))
28
29# Retrieve raw & normalized values of RooFit p.d.f.s
30# --------------------------------------------------------------------------------------------------
31
32# Return 'raw' unnormalized value of gx
33print("gx = ", gx.getVal())
34
35# Return value of gx normalized over x in range [-10,10]
36nset = ROOT.RooArgSet(x)
37print("gx_Norm[x] = ", gx.getVal(nset))
38
39# Create object representing integral over gx
40# which is used to calculate gx_Norm[x] == gx / gx_Int[x]
41igx = gx.createIntegral(ROOT.RooArgSet(x))
42print("gx_Int[x] = ", igx.getVal())
43
44# Integrate normalized pdf over subrange
45# ----------------------------------------------------------------------------
46
47# Define a range named "signal" in x from -5,5
48x.setRange("signal", -5, 5)
49
50# Create an integral of gx_Norm[x] over x in range "signal"
51# ROOT.This is the fraction of of p.d.f. gx_Norm[x] which is in the
52# range named "signal"
53igx_sig = gx.createIntegral(ROOT.RooArgSet(x), ROOT.RooFit.NormSet(
54 ROOT.RooArgSet(x)), ROOT.RooFit.Range("signal"))
55print("gx_Int[x|signal]_Norm[x] = ", igx_sig.getVal())
56
57# Construct cumulative distribution function from pdf
58# -----------------------------------------------------------------------------------------------------
59
60# Create the cumulative distribution function of gx
61# i.e. calculate Int[-10,x] gx(x') dx'
62gx_cdf = gx.createCdf(ROOT.RooArgSet(x))
63
64# Plot cdf of gx versus x
65frame = x.frame(ROOT.RooFit.Title("c.d.f of Gaussian p.d.f"))
66gx_cdf.plotOn(frame)
67
68# Draw plot on canvas
69c = ROOT.TCanvas("rf110_normintegration",
70 "rf110_normintegration", 600, 600)
71ROOT.gPad.SetLeftMargin(0.15)
72frame.GetYaxis().SetTitleOffset(1.6)
73frame.Draw()
74
75c.SaveAs("rf110_normintegration.png")