Logo ROOT  
Reference Guide
rf301_composition.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## Multidimensional models: multi-dimensional p.d.f.s through composition, e.g. substituting
6## a p.d.f parameter with a function that depends on other observables
7##
8## `pdf = gauss(x,f(y),s)` with `f(y) = a0 + a1*y`
9##
10## \macro_code
11##
12## \date February 2018
13## \authors Clemens Lange, Wouter Verkerke (C++ version)
14
15import ROOT
16
17# Setup composed model gauss(x, m(y), s)
18# -----------------------------------------------------------------------
19
20# Create observables
21x = ROOT.RooRealVar("x", "x", -5, 5)
22y = ROOT.RooRealVar("y", "y", -5, 5)
23
24# Create function f(y) = a0 + a1*y
25a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
26a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
27fy = ROOT.RooPolyVar("fy", "fy", y, ROOT.RooArgList(a0, a1))
28
29# Creat gauss(x,f(y),s)
30sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
31model = ROOT.RooGaussian(
32 "model", "Gaussian with shifting mean", x, fy, sigma)
33
34# Sample data, plot data and pdf on x and y
35# ---------------------------------------------------------------------------------
36
37# Generate 10000 events in x and y from model
38data = model.generate(ROOT.RooArgSet(x, y), 10000)
39
40# Plot x distribution of data and projection of model x = Int(dy)
41# model(x,y)
42xframe = x.frame()
43data.plotOn(xframe)
44model.plotOn(xframe)
45
46# Plot x distribution of data and projection of model y = Int(dx)
47# model(x,y)
48yframe = y.frame()
49data.plotOn(yframe)
50model.plotOn(yframe)
51
52# Make two-dimensional plot in x vs y
53hh_model = model.createHistogram("hh_model", x, ROOT.RooFit.Binning(
54 50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
55hh_model.SetLineColor(ROOT.kBlue)
56
57# Make canvas and draw ROOT.RooPlots
58c = ROOT.TCanvas("rf301_composition", "rf301_composition", 1200, 400)
59c.Divide(3)
60c.cd(1)
61ROOT.gPad.SetLeftMargin(0.15)
62xframe.GetYaxis().SetTitleOffset(1.4)
63xframe.Draw()
64c.cd(2)
65ROOT.gPad.SetLeftMargin(0.15)
66yframe.GetYaxis().SetTitleOffset(1.4)
67yframe.Draw()
68c.cd(3)
69ROOT.gPad.SetLeftMargin(0.20)
70hh_model.GetZaxis().SetTitleOffset(2.5)
71hh_model.Draw("surf")
72
73c.SaveAs("rf301_composition.png")