Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf208_convolution.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #208
5## One-dimensional numeric convolution
6## (require ROOT to be compiled with --enable-fftw3)
7##
8## pdf = landau(t) (x) gauss(t)
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# Construct observable
22t = ROOT.RooRealVar("t", "t", -10, 30)
23
24# Construct landau(t,ml,sl)
25ml = ROOT.RooRealVar("ml", "mean landau", 5.0, -20, 20)
26sl = ROOT.RooRealVar("sl", "sigma landau", 1, 0.1, 10)
27landau = ROOT.RooLandau("lx", "lx", t, ml, sl)
28
29# Construct gauss(t,mg,sg)
30mg = ROOT.RooRealVar("mg", "mg", 0)
31sg = ROOT.RooRealVar("sg", "sg", 2, 0.1, 10)
32gauss = ROOT.RooGaussian("gauss", "gauss", t, mg, sg)
33
34# Construct convolution pdf
35# ---------------------------------------
36
37# Set #bins to be used for FFT sampling to 10000
38t.setBins(10000, "cache")
39
40# Construct landau (x) gauss
41lxg = ROOT.RooFFTConvPdf("lxg", "landau (X) gauss", t, landau, gauss)
42
43# Sample, fit and plot convoluted pdf
44# ----------------------------------------------------------------------
45
46# Sample 1000 events in x from gxlx
47data = lxg.generate({t}, 10000)
48
49# Fit gxlx to data
50lxg.fitTo(data)
51
52# Plot data, pdf, landau (X) gauss pdf
53frame = t.frame(Title="landau (x) gauss convolution")
54data.plotOn(frame)
55lxg.plotOn(frame)
56landau.plotOn(frame, LineStyle="--")
57
58# Draw frame on canvas
59c = ROOT.TCanvas("rf208_convolution", "rf208_convolution", 600, 600)
60ROOT.gPad.SetLeftMargin(0.15)
61frame.GetYaxis().SetTitleOffset(1.4)
62frame.Draw()
63
64c.SaveAs("rf208_convolution.png")