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