'ADDITION AND CONVOLUTION' RooFit tutorial macro #208 One-dimensional numeric convolution (require ROOT to be compiled with –enable-fftw3)
pdf = landau(t) (x) gauss(t)
import ROOT
t = ROOT.RooRealVar("t", "t", -10, 30)
ml = ROOT.RooRealVar("ml", "mean landau", 5.0, -20, 20)
sl = ROOT.RooRealVar("sl", "sigma landau", 1, 0.1, 10)
landau = ROOT.RooLandau("lx", "lx", t, ml, sl)
mg = ROOT.RooRealVar("mg", "mg", 0)
sg = ROOT.RooRealVar("sg", "sg", 2, 0.1, 10)
gauss = ROOT.RooGaussian("gauss", "gauss", t, mg, sg)
t.setBins(10000, "cache")
lxg = ROOT.RooFFTConvPdf("lxg", "landau (X) gauss", t, landau, gauss)
data = lxg.generate({t}, 10000)
lxg.fitTo(data, PrintLevel=-1)
frame = t.frame(Title="landau (x) gauss convolution")
data.plotOn(frame)
lxg.plotOn(frame)
landau.plotOn(frame, LineStyle="--")
c = ROOT.TCanvas("rf208_convolution", "rf208_convolution", 600, 600)
ROOT.gPad.SetLeftMargin(0.15)
frame.GetYaxis().SetTitleOffset(1.4)
frame.Draw()
c.SaveAs("rf208_convolution.png")
[#1] INFO:Eval -- RooRealVar::setRange(t) new range named 'refrange_fft_lxg' created with bounds [-10,30]
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x85ab260 with pdf lx_CONV_gauss_CACHE_Obs[t]_NORM_t for nset (t) with code 0
[#1] INFO:Fitting -- RooAbsPdf::fitTo(lxg_over_lxg_Int[t]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_lxg_over_lxg_Int[t]_lxgData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x88084e0 with pdf lx_CONV_gauss_CACHE_Obs[t] for nset () with code 1 from preexisting content.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x82935f0 with pdf lx_CONV_gauss_CACHE_Obs[t]_NORM_t for nset (t) with code 0
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C version)
Definition in file rf208_convolution.py.