Logo ROOT  
Reference Guide
rf211_paramconv.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #211
5## Working a with a p.d.f. with a convolution operator in terms
6## of a parameter
7##
8## (require ROOT to be compiled with --enable-fftw3)
9##
10## \macro_code
11##
12## \date February 2018
13## \authors Clemens Lange, Wouter Verkerke (C version)
14
15import ROOT
16
17# Set up component pdfs
18# ---------------------------------------
19
20# Gaussian g(x ; mean,sigma)
21x = ROOT.RooRealVar("x", "x", -10, 10)
22mean = ROOT.RooRealVar("mean", "mean", -3, 3)
23sigma = ROOT.RooRealVar("sigma", "sigma", 0.5, 0.1, 10)
24modelx = ROOT.RooGaussian("gx", "gx", x, mean, sigma)
25
26# Block function in mean
27a = ROOT.RooRealVar("a", "a", 2, 1, 10)
28model_mean = ROOT.RooGenericPdf("model_mean", "abs(mean)<a", [mean, a])
29
30# Convolution in mean model = g(x,mean,sigma) (x) block(mean)
31x.setBins(1000, "cache")
32mean.setBins(50, "cache")
33model = ROOT.RooFFTConvPdf("model", "model", mean, modelx, model_mean)
34
35# Configure convolution to construct a 2-D cache in (x,mean)
36# rather than a 1-d cache in mean that needs to be recalculated
37# for each value of x
38model.setCacheObservables({x})
39model.setBufferFraction(1.0)
40
41# Integrate model over projModel = Int model dmean
42projModel = model.createProjection({mean})
43
44# Generate 1000 toy events
45d = projModel.generateBinned({x}, 1000)
46
47# Fit p.d.f. to toy data
48projModel.fitTo(d, Verbose=True)
49
50# Plot data and fitted p.d.f.
51frame = x.frame(Bins=25)
52d.plotOn(frame)
53projModel.plotOn(frame)
54
55# Make 2d histogram of model(x;mean)
56hh = model.createHistogram(
57 "hh",
58 x,
59 Binning=50,
60 YVar=dict(var=mean, Binning=50),
61 ConditionalObservables={mean},
62)
63hh.SetTitle("histogram of model(x|mean)")
64hh.SetLineColor(ROOT.kBlue)
65
66# Draw frame on canvas
67c = ROOT.TCanvas("rf211_paramconv", "rf211_paramconv", 800, 400)
68c.Divide(2)
69c.cd(1)
70ROOT.gPad.SetLeftMargin(0.15)
71frame.GetYaxis().SetTitleOffset(1.4)
72frame.Draw()
73c.cd(2)
74ROOT.gPad.SetLeftMargin(0.20)
75hh.GetZaxis().SetTitleOffset(2.5)
76hh.Draw("surf")
77
78c.SaveAs("rf211_paramconv.png")