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