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_code
10##
11## \date February 2018
12## \authors Clemens Lange, Wouter Verkerke (C++ version)
13
14import ROOT
15
16# Setup composed model gauss(x, m(y), s)
17# -----------------------------------------------------------------------
18
19# Create observables
20x = ROOT.RooRealVar("x", "x", -5, 5)
21y = ROOT.RooRealVar("y", "y", -5, 5)
22
23# Create function f(y) = a0 + a1*y
24a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
25a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
26fy = ROOT.RooPolyVar("fy", "fy", y, [a0, a1])
27
28# Creat gauss(x,f(y),s)
29sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
30model = ROOT.RooGaussian("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({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, Binning=50, YVar=dict(var=y, Binning=50))
52hh_model.SetLineColor(ROOT.kBlue)
53
54# Make canvas and draw ROOT.RooPlots
55c = ROOT.TCanvas("rf301_composition", "rf301_composition", 1200, 400)
56c.Divide(3)
57c.cd(1)
58ROOT.gPad.SetLeftMargin(0.15)
59xframe.GetYaxis().SetTitleOffset(1.4)
60xframe.Draw()
61c.cd(2)
62ROOT.gPad.SetLeftMargin(0.15)
63yframe.GetYaxis().SetTitleOffset(1.4)
64yframe.Draw()
65c.cd(3)
66ROOT.gPad.SetLeftMargin(0.20)
67hh_model.GetZaxis().SetTitleOffset(2.5)
68hh_model.Draw("surf")
69
70c.SaveAs("rf301_composition.png")