Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf308_normintegration2d.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: normalization and integration of pdfs, construction of
5## cumulative distribution functions from pdfs in two dimensions
6##
7## \macro_image
8## \macro_code
9## \macro_output
10##
11## \date February 2018
12## \authors Clemens Lange, Wouter Verkerke (C++ version)
13
14import ROOT
15
16# Set up model
17# ---------------------
18
19# Create observables x,y
20x = ROOT.RooRealVar("x", "x", -10, 10)
21y = ROOT.RooRealVar("y", "y", -10, 10)
22
23# Create pdf gaussx(x,-2,3), gaussy(y,2,2)
24gx = ROOT.RooGaussian("gx", "gx", x, -2.0, 3.0)
25gy = ROOT.RooGaussian("gy", "gy", y, +2.0, 2.0)
26
27# gxy = gx(x)*gy(y)
28gxy = ROOT.RooProdPdf("gxy", "gxy", [gx, gy])
29
30# Retrieve raw & normalized values of RooFit pdfs
31# --------------------------------------------------------------------------------------------------
32
33# Return 'raw' unnormalized value of gx
34print("gxy = ", gxy.getVal())
35
36# Return value of gxy normalized over x _and_ y in range [-10,10]
37nset_xy = {x, y}
38print("gx_Norm[x,y] = ", gxy.getVal(nset_xy))
39
40# Create object representing integral over gx
41# which is used to calculate gx_Norm[x,y] == gx / gx_Int[x,y]
42x_and_y = {x, y}
43igxy = gxy.createIntegral(x_and_y)
44print("gx_Int[x,y] = ", igxy.getVal())
45
46# NB: it is also possible to do the following
47
48# Return value of gxy normalized over x in range [-10,10] (i.e. treating y
49# as parameter)
50nset_x = {x}
51print("gx_Norm[x] = ", gxy.getVal(nset_x))
52
53# Return value of gxy normalized over y in range [-10,10] (i.e. treating x
54# as parameter)
55nset_y = {y}
56print("gx_Norm[y] = ", gxy.getVal(nset_y))
57
58# Integrate normalized pdf over subrange
59# ----------------------------------------------------------------------------
60
61# Define a range named "signal" in x from -5,5
62x.setRange("signal", -5, 5)
63y.setRange("signal", -3, 3)
64
65# Create an integral of gxy_Norm[x,y] over x and y in range "signal"
66# ROOT.This is the fraction of of pdf gxy_Norm[x,y] which is in the
67# range named "signal"
68
69igxy_sig = gxy.createIntegral(x_and_y, NormSet=x_and_y, Range="signal")
70print("gx_Int[x,y|signal]_Norm[x,y] = ", igxy_sig.getVal())
71
72# Construct cumulative distribution function from pdf
73# -----------------------------------------------------------------------------------------------------
74
75# Create the cumulative distribution function of gx
76# i.e. calculate Int[-10,x] gx(x') dx'
77gxy_cdf = gxy.createCdf({x, y})
78
79# Plot cdf of gx versus x
80hh_cdf = gxy_cdf.createHistogram("hh_cdf", x, Binning=40, YVar=dict(var=y, Binning=40))
81hh_cdf.SetLineColor(ROOT.kBlue)
82
83c = ROOT.TCanvas("rf308_normintegration2d", "rf308_normintegration2d", 600, 600)
84ROOT.gPad.SetLeftMargin(0.15)
85hh_cdf.GetZaxis().SetTitleOffset(1.8)
86hh_cdf.Draw("surf")
87
88c.SaveAs("rf308_normintegration2d.png")