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