Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf303_conditional.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #303
5## Use of tailored p.d.f as conditional p.d.fs.s
6##
7## pdf = gauss(x,f(y),sx | y ) with f(y) = a0 + a1*y
8##
9## \macro_image
10## \macro_code
11## \macro_output
12##
13## \date February 2018
14## \authors Clemens Lange, Wouter Verkerke (C version)
15
16import ROOT
17
18
19def makeFakeDataXY():
20 x = ROOT.RooRealVar("x", "x", -10, 10)
21 y = ROOT.RooRealVar("y", "y", -10, 10)
22 coord = {x, y}
23
24 d = ROOT.RooDataSet("d", "d", coord)
25
26 for i in range(10000):
27 tmpy = ROOT.gRandom.Gaus(0, 10)
28 tmpx = ROOT.gRandom.Gaus(0.5 * tmpy, 1)
29 if (abs(tmpy) < 10) and (abs(tmpx) < 10):
30 x.setVal(tmpx)
31 y.setVal(tmpy)
32 d.add(coord)
33
34 return d
35
36
37# Set up composed model gauss(x, m(y), s)
38# -----------------------------------------------------------------------
39
40# Create observables
41x = ROOT.RooRealVar("x", "x", -10, 10)
42y = ROOT.RooRealVar("y", "y", -10, 10)
43
44# Create function f(y) = a0 + a1*y
45a0 = ROOT.RooRealVar("a0", "a0", -0.5, -5, 5)
46a1 = ROOT.RooRealVar("a1", "a1", -0.5, -1, 1)
47fy = ROOT.RooPolyVar("fy", "fy", y, [a0, a1])
48
49# Creat gauss(x,f(y),s)
50sigma = ROOT.RooRealVar("sigma", "width of gaussian", 0.5, 0.1, 2.0)
51model = ROOT.RooGaussian("model", "Gaussian with shifting mean", x, fy, sigma)
52
53# Obtain fake external experimental dataset with values for x and y
54expDataXY = makeFakeDataXY()
55
56# Generate data from conditional p.d.f. model(x|y)
57# ---------------------------------------------------------------------------------------------
58
59# Make subset of experimental data with only y values
60expDataY = expDataXY.reduce({y})
61
62# Generate 10000 events in x obtained from _conditional_ model(x|y) with y
63# values taken from experimental data
64data = model.generate({x}, ProtoData=expDataY)
65data.Print()
66
67# Fit conditional p.d.f model(x|y) to data
68# ---------------------------------------------------------------------------------------------
69
70model.fitTo(expDataXY, ConditionalObservables={y}, PrintLevel=-1)
71
72# Project conditional p.d.f on x and y dimensions
73# ---------------------------------------------------------------------------------------------
74
75# Plot x distribution of data and projection of model x = 1/Ndata
76# sum(data(y_i)) model(x;y_i)
77xframe = x.frame()
78expDataXY.plotOn(xframe)
79model.plotOn(xframe, ProjWData=expDataY)
80
81# Speed up (and approximate) projection by using binned clone of data for
82# projection
83binnedDataY = expDataY.binnedClone()
84model.plotOn(xframe, ProjWData=binnedDataY, LineColor="c", LineStyle=":")
85
86# Show effect of projection with too coarse binning
87(expDataY.get().find("y")).setBins(5)
88binnedDataY2 = expDataY.binnedClone()
89model.plotOn(xframe, ProjWData=binnedDataY2, LineColor="r")
90
91# Make canvas and draw ROOT.RooPlots
92c = ROOT.TCanvas("rf303_conditional", "rf303_conditional", 600, 460)
93ROOT.gPad.SetLeftMargin(0.15)
94xframe.GetYaxis().SetTitleOffset(1.2)
95xframe.Draw()
96
97c.SaveAs("rf303_conditional.png")