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, ROOT.RooArgList(a0, a1))
27
28# Creat gauss(x,f(y),s)
29sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
30model = ROOT.RooGaussian(
31 "model", "Gaussian with shifting mean", x, fy, sigma)
32
33# Sample data, plot data and pdf on x and y
34# ---------------------------------------------------------------------------------
35
36# Generate 10000 events in x and y from model
37data = model.generate(ROOT.RooArgSet(x, y), 10000)
38
39# Plot x distribution of data and projection of model x = Int(dy)
40# model(x,y)
41xframe = x.frame()
42data.plotOn(xframe)
43model.plotOn(xframe)
44
45# Plot x distribution of data and projection of model y = Int(dx)
46# model(x,y)
47yframe = y.frame()
48data.plotOn(yframe)
49model.plotOn(yframe)
50
51# Make two-dimensional plot in x vs y
52hh_model = model.createHistogram("hh_model", x, ROOT.RooFit.Binning(
53 50), ROOT.RooFit.YVar(y, ROOT.RooFit.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")