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