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