Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
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## \author Clemens Lange
14## \author Wouter Verkerke (C version)
15
16import ROOT
17
18# Set up component pdfs
19# ---------------------------------------
20
21# Gaussian g(x ; mean,sigma)
22x = ROOT.RooRealVar("x", "x", -10, 10)
23mean = ROOT.RooRealVar("mean", "mean", -3, 3)
24sigma = ROOT.RooRealVar("sigma", "sigma", 0.5, 0.1, 10)
25modelx = ROOT.RooGaussian("gx", "gx", x, mean, sigma)
26
27# Block function in mean
28a = ROOT.RooRealVar("a", "a", 2, 1, 10)
29model_mean = ROOT.RooGenericPdf("model_mean", "abs(mean)<a", [mean, a])
30
31# Convolution in mean model = g(x,mean,sigma) (x) block(mean)
32x.setBins(1000, "cache")
33mean.setBins(50, "cache")
34model = ROOT.RooFFTConvPdf("model", "model", mean, modelx, model_mean)
35
36# Configure convolution to construct a 2-D cache in (x,mean)
37# rather than a 1-d cache in mean that needs to be recalculated
38# for each value of x
39model.setCacheObservables({x})
40model.setBufferFraction(1.0)
41
42# Integrate model over projModel = Int model dmean
43projModel = model.createProjection({mean})
44
45# Generate 1000 toy events
46d = projModel.generateBinned({x}, 1000)
47
48# Fit p.d.f. to toy data
49projModel.fitTo(d, Verbose=True)
50
51# Plot data and fitted p.d.f.
52frame = x.frame(Bins=25)
53d.plotOn(frame)
54projModel.plotOn(frame)
55
56# Make 2d histogram of model(x;mean)
57hh = model.createHistogram(
58 "hh",
59 x,
60 Binning=50,
61 YVar=dict(var=mean, Binning=50),
62 ConditionalObservables={mean},
63)
64hh.SetTitle("histogram of model(x|mean)")
65hh.SetLineColor(ROOT.kBlue)
66
67# Draw frame on canvas
68c = ROOT.TCanvas("rf211_paramconv", "rf211_paramconv", 800, 400)
69c.Divide(2)
70c.cd(1)
71ROOT.gPad.SetLeftMargin(0.15)
72frame.GetYaxis().SetTitleOffset(1.4)
73frame.Draw()
74c.cd(2)
75ROOT.gPad.SetLeftMargin(0.20)
76hh.GetZaxis().SetTitleOffset(2.5)
77hh.Draw("surf")
78
79c.SaveAs("rf211_paramconv.png")