Logo ROOT   6.18/05
Reference Guide
rf308_normintegration2d.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## Multidimensional models: normalization and integration of p.d.fs, construction of cumulative distribution functions from p.d.f.s in two dimensions
6##
7## \macro_code
8##
9## \date February 2018
10## \author 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)
20y = ROOT.RooRealVar("y", "y", -10, 10)
21
22# Create p.d.f. gaussx(x,-2,3), gaussy(y,2,2)
23gx = ROOT.RooGaussian(
24 "gx", "gx", x, ROOT.RooFit.RooConst(-2), ROOT.RooFit.RooConst(3))
25gy = ROOT.RooGaussian(
26 "gy", "gy", y, ROOT.RooFit.RooConst(+2), ROOT.RooFit.RooConst(2))
27
28# gxy = gx(x)*gy(y)
29gxy = ROOT.RooProdPdf("gxy", "gxy", ROOT.RooArgList(gx, gy))
30
31# Retrieve raw & normalized values of RooFit p.d.f.s
32# --------------------------------------------------------------------------------------------------
33
34# Return 'raw' unnormalized value of gx
35print("gxy = ", gxy.getVal())
36
37# Return value of gxy normalized over x _and_ y in range [-10,10]
38nset_xy = ROOT.RooArgSet(x, y)
39print("gx_Norm[x,y] = ", gxy.getVal(nset_xy))
40
41# Create object representing integral over gx
42# which is used to calculate gx_Norm[x,y] == gx / gx_Int[x,y]
43x_and_y = ROOT.RooArgSet(x, y)
44igxy = gxy.createIntegral(x_and_y)
45print("gx_Int[x,y] = ", igxy.getVal())
46
47# NB: it is also possible to do the following
48
49# Return value of gxy normalized over x in range [-10,10] (i.e. treating y
50# as parameter)
51nset_x = ROOT.RooArgSet(x)
52print("gx_Norm[x] = ", gxy.getVal(nset_x))
53
54# Return value of gxy normalized over y in range [-10,10] (i.e. treating x
55# as parameter)
56nset_y = ROOT.RooArgSet(y)
57print("gx_Norm[y] = ", gxy.getVal(nset_y))
58
59# Integarte normalizes pdf over subrange
60# ----------------------------------------------------------------------------
61
62# Define a range named "signal" in x from -5,5
63x.setRange("signal", -5, 5)
64y.setRange("signal", -3, 3)
65
66# Create an integral of gxy_Norm[x,y] over x and y in range "signal"
67# ROOT.This is the fraction of of p.d.f. gxy_Norm[x,y] which is in the
68# range named "signal"
69
70igxy_sig = gxy.createIntegral(x_and_y, ROOT.RooFit.NormSet(
71 x_and_y), ROOT.RooFit.Range("signal"))
72print("gx_Int[x,y|signal]_Norm[x,y] = ", igxy_sig.getVal())
73
74# Construct cumulative distribution function from pdf
75# -----------------------------------------------------------------------------------------------------
76
77# Create the cumulative distribution function of gx
78# i.e. calculate Int[-10,x] gx(x') dx'
79gxy_cdf = gxy.createCdf(ROOT.RooArgSet(x, y))
80
81# Plot cdf of gx versus x
82hh_cdf = gxy_cdf.createHistogram("hh_cdf", x, ROOT.RooFit.Binning(
83 40), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(40)))
84hh_cdf.SetLineColor(ROOT.kBlue)
85
86c = ROOT.TCanvas("rf308_normintegration2d",
87 "rf308_normintegration2d", 600, 600)
88ROOT.gPad.SetLeftMargin(0.15)
89hh_cdf.GetZaxis().SetTitleOffset(1.8)
90hh_cdf.Draw("surf")
91
92c.SaveAs("rf308_normintegration2d.png")