Logo ROOT   6.16/01
Reference Guide
rf315_projectpdf.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #315
6##
7## Marginizalization of multi-dimensional p.d.f.s through integration
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15
16import ROOT
17
18
19# Create pdf m(x,y) = gx(x|y) * g(y)
20# --------------------------------------------------------------
21
22# Increase default precision of numeric integration
23# as self exercise has high sensitivity to numeric integration precision
24ROOT.RooAbsPdf.defaultIntegratorConfig().setEpsRel(1e-8)
25ROOT.RooAbsPdf.defaultIntegratorConfig().setEpsAbs(1e-8)
26
27# Create observables
28x = ROOT.RooRealVar("x", "x", -5, 5)
29y = ROOT.RooRealVar("y", "y", -2, 2)
30
31# Create function f(y) = a0 + a1*y
32a0 = ROOT.RooRealVar("a0", "a0", 0)
33a1 = ROOT.RooRealVar("a1", "a1", -1.5, -3, 1)
34fy = ROOT.RooPolyVar("fy", "fy", y, ROOT.RooArgList(a0, a1))
35
36# Create gaussx(x,f(y),sx)
37sigmax = ROOT.RooRealVar("sigmax", "width of gaussian", 0.5)
38gaussx = ROOT.RooGaussian(
39 "gaussx", "Gaussian in x with shifting mean in y", x, fy, sigmax)
40
41# Create gaussy(y,0,2)
42gaussy = ROOT.RooGaussian(
43 "gaussy", "Gaussian in y", y, ROOT.RooFit.RooConst(0), ROOT.RooFit.RooConst(2))
44
45# Create gaussx(x,sx|y) * gaussy(y)
46model = ROOT.RooProdPdf("model", "gaussx(x|y)*gaussy(y)", ROOT.RooArgSet(
47 gaussy), ROOT.RooFit.Conditional(ROOT.RooArgSet(gaussx), ROOT.RooArgSet(x)))
48
49# Marginalize m(x,y) to m(x)
50# ----------------------------------------------------
51
52# modelx(x) = Int model(x,y) dy
53modelx = model.createProjection(ROOT.RooArgSet(y))
54
55# Use marginalized p.d.f. as regular 1D p.d.f.
56# -----------------------------------------------
57
58# Sample 1000 events from modelx
59data = modelx.generateBinned(ROOT.RooArgSet(x), 1000)
60
61# Fit modelx to toy data
62modelx.fitTo(data, ROOT.RooFit.Verbose())
63
64# Plot modelx over data
65frame = x.frame(40)
66data.plotOn(frame)
67modelx.plotOn(frame)
68
69# Make 2D histogram of model(x,y)
70hh = model.createHistogram("x,y")
71hh.SetLineColor(ROOT.kBlue)
72
73c = ROOT.TCanvas("rf315_projectpdf", "rf315_projectpdf", 800, 400)
74c.Divide(2)
75c.cd(1)
76ROOT.gPad.SetLeftMargin(0.15)
77frame.GetYaxis().SetTitleOffset(1.4)
78frame.Draw()
79c.cd(2)
80ROOT.gPad.SetLeftMargin(0.20)
81hh.GetZaxis().SetTitleOffset(2.5)
82hh.Draw("surf")
83c.SaveAs("rf315_projectpdf.png")