Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf301_composition.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Multidimensional models: multi-dimensional pdfs through composition, e.g. substituting
5## a pdf parameter with a function that depends on other observables
6##
7## `pdf = gauss(x,f(y),s)` with `f(y) = a0 + a1*y`
8##
9## \macro_image
10## \macro_code
11## \macro_output
12##
13## \date February 2018
14## \authors Clemens Lange, Wouter Verkerke (C++ version)
15
16import ROOT
17
18# Setup composed model gauss(x, m(y), s)
19# -----------------------------------------------------------------------
20
21# Create observables
22x = ROOT.RooRealVar("x", "x", -5, 5)
23y = ROOT.RooRealVar("y", "y", -5, 5)
24
25# Create function f(y) = a0 + a1*y
26a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
27a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
28fy = ROOT.RooPolyVar("fy", "fy", y, [a0, a1])
29
30# Creat gauss(x,f(y),s)
31sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
32model = ROOT.RooGaussian("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({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, Binning=50, YVar=dict(var=y, Binning=50))
54hh_model.SetLineColor(ROOT.kBlue)
55
56# Make canvas and draw ROOT.RooPlots
57c = ROOT.TCanvas("rf301_composition", "rf301_composition", 1200, 400)
58c.Divide(3)
59c.cd(1)
60ROOT.gPad.SetLeftMargin(0.15)
61xframe.GetYaxis().SetTitleOffset(1.4)
62xframe.Draw()
63c.cd(2)
64ROOT.gPad.SetLeftMargin(0.15)
65yframe.GetYaxis().SetTitleOffset(1.4)
66yframe.Draw()
67c.cd(3)
68ROOT.gPad.SetLeftMargin(0.20)
69hh_model.GetZaxis().SetTitleOffset(2.5)
70hh_model.Draw("surf")
71
72c.SaveAs("rf301_composition.png")