Logo ROOT  
Reference Guide
rf301_composition.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ##
5 ## \brief 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 
15 import ROOT
16 
17 # Setup composed model gauss(x, m(y), s)
18 # -----------------------------------------------------------------------
19 
20 # Create observables
21 x = ROOT.RooRealVar("x", "x", -5, 5)
22 y = ROOT.RooRealVar("y", "y", -5, 5)
23 
24 # Create function f(y) = a0 + a1*y
25 a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
26 a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
27 fy = ROOT.RooPolyVar("fy", "fy", y, ROOT.RooArgList(a0, a1))
28 
29 # Creat gauss(x,f(y),s)
30 sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
31 model = 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
38 data = 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)
42 xframe = x.frame()
43 data.plotOn(xframe)
44 model.plotOn(xframe)
45 
46 # Plot x distribution of data and projection of model y = Int(dx)
47 # model(x,y)
48 yframe = y.frame()
49 data.plotOn(yframe)
50 model.plotOn(yframe)
51 
52 # Make two-dimensional plot in x vs y
53 hh_model = model.createHistogram("hh_model", x, ROOT.RooFit.Binning(
54  50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
55 hh_model.SetLineColor(ROOT.kBlue)
56 
57 # Make canvas and draw ROOT.RooPlots
58 c = ROOT.TCanvas("rf301_composition", "rf301_composition", 1200, 400)
59 c.Divide(3)
60 c.cd(1)
61 ROOT.gPad.SetLeftMargin(0.15)
62 xframe.GetYaxis().SetTitleOffset(1.4)
63 xframe.Draw()
64 c.cd(2)
65 ROOT.gPad.SetLeftMargin(0.15)
66 yframe.GetYaxis().SetTitleOffset(1.4)
67 yframe.Draw()
68 c.cd(3)
69 ROOT.gPad.SetLeftMargin(0.20)
70 hh_model.GetZaxis().SetTitleOffset(2.5)
71 hh_model.Draw("surf")
72 
73 c.SaveAs("rf301_composition.png")