Logo ROOT  
Reference Guide
rf305_condcorrprod.py
Go to the documentation of this file.
1 ## \file
2 ## \ingroup tutorial_roofit
3 ## \notebook
4 ##
5 ## Multidimensional models: multi-dimensional p.d.f.s with conditional p.d.fs in product
6 ##
7 ## `pdf = gauss(x,f(y),sx | y ) * gauss(y,ms,sx)` with `f(y) = a0 + a1*y`
8 ##
9 ## \macro_code
10 ##
11 ## \date February 2018
12 ## \authors Clemens Lange, Wouter Verkerke (C++ version)
13 
14 import ROOT
15 
16 # Create conditional pdf gx(x|y)
17 # -----------------------------------------------------------
18 
19 # Create observables
20 x = ROOT.RooRealVar("x", "x", -5, 5)
21 y = ROOT.RooRealVar("y", "y", -5, 5)
22 
23 # Create function f(y) = a0 + a1*y
24 a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
25 a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
26 fy = ROOT.RooPolyVar("fy", "fy", y, ROOT.RooArgList(a0, a1))
27 
28 # Create gaussx(x,f(y),sx)
29 sigmax = ROOT.RooRealVar("sigma", "width of gaussian", 0.5)
30 gaussx = ROOT.RooGaussian(
31  "gaussx", "Gaussian in x with shifting mean in y", x, fy, sigmax)
32 
33 # Create pdf gy(y)
34 # -----------------------------------------------------------
35 
36 # Create gaussy(y,0,5)
37 gaussy = ROOT.RooGaussian(
38  "gaussy",
39  "Gaussian in y",
40  y,
41  ROOT.RooFit.RooConst(0),
42  ROOT.RooFit.RooConst(3))
43 
44 # Create product gx(x|y)*gy(y)
45 # -------------------------------------------------------
46 
47 # Create gaussx(x,sx|y) * gaussy(y)
48 model = ROOT.RooProdPdf(
49  "model",
50  "gaussx(x|y)*gaussy(y)",
51  ROOT.RooArgSet(gaussy),
52  ROOT.RooFit.Conditional(
53  ROOT.RooArgSet(gaussx),
54  ROOT.RooArgSet(x)))
55 
56 # Sample, fit and plot product pdf
57 # ---------------------------------------------------------------
58 
59 # Generate 1000 events in x and y from model
60 data = model.generate(ROOT.RooArgSet(x, y), 10000)
61 
62 # Plot x distribution of data and projection of model x = Int(dy)
63 # model(x,y)
64 xframe = x.frame()
65 data.plotOn(xframe)
66 model.plotOn(xframe)
67 
68 # Plot x distribution of data and projection of model y = Int(dx)
69 # model(x,y)
70 yframe = y.frame()
71 data.plotOn(yframe)
72 model.plotOn(yframe)
73 
74 # Make two-dimensional plot in x vs y
75 hh_model = model.createHistogram("hh_model", x, ROOT.RooFit.Binning(
76  50), ROOT.RooFit.YVar(y, ROOT.RooFit.Binning(50)))
77 hh_model.SetLineColor(ROOT.kBlue)
78 
79 # Make canvas and draw ROOT.RooPlots
80 c = ROOT.TCanvas("rf305_condcorrprod", "rf05_condcorrprod", 1200, 400)
81 c.Divide(3)
82 c.cd(1)
83 ROOT.gPad.SetLeftMargin(0.15)
84 xframe.GetYaxis().SetTitleOffset(1.6)
85 xframe.Draw()
86 c.cd(2)
87 ROOT.gPad.SetLeftMargin(0.15)
88 yframe.GetYaxis().SetTitleOffset(1.6)
89 yframe.Draw()
90 c.cd(3)
91 ROOT.gPad.SetLeftMargin(0.20)
92 hh_model.GetZaxis().SetTitleOffset(2.5)
93 hh_model.Draw("surf")
94 
95 c.SaveAs("rf305_condcorrprod.png")