Logo ROOT   6.14/05
Reference Guide
rf208_convolution.C File Reference

Detailed Description

View in nbviewer Open in SWAN '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)

pict1_rf208_convolution.C.png
Processing /mnt/build/workspace/root-makedoc-v614/rootspi/rdoc/src/v6-14-00-patches/tutorials/roofit/rf208_convolution.C...
RooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:Eval -- RooRealVar::setRange(t) new range named 'refrange_fft_lxg' created with bounds [-10,30]
[#1] INFO:NumericIntegration -- RooRealIntegral::init(lx_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x14a3e60 with pdf lx_CONV_gauss_CACHE_Obs[t] for nset (t) with code 0
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 ml 5.00000e+00 4.00000e+00 -2.00000e+01 2.00000e+01
2 sg 2.00000e+00 9.50000e-01 1.00000e-01 1.00000e+01
3 sl 1.00000e+00 4.50000e-01 1.00000e-01 1.00000e+01
**********
** 3 **SET ERR 0.5
**********
**********
** 4 **SET PRINT 1
**********
**********
** 5 **SET STR 1
**********
NOW USING STRATEGY 1: TRY TO BALANCE SPEED AGAINST RELIABILITY
**********
** 6 **MIGRAD 1500 1
**********
FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=28285.6 FROM MIGRAD STATUS=INITIATE 12 CALLS 13 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 ml 5.00000e+00 4.00000e+00 2.08372e-01 1.28490e+03
2 sg 2.00000e+00 9.50000e-01 2.51381e-01 1.33500e+02
3 sl 1.00000e+00 4.50000e-01 1.63378e-01 8.59464e+01
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=28280.5 FROM MIGRAD STATUS=CONVERGED 63 CALLS 64 TOTAL
EDM=7.08918e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 ml 4.88711e+00 3.70488e-02 1.95153e-04 -4.21500e-01
2 sg 1.89676e+00 4.61915e-02 9.28724e-04 4.64118e-02
3 sl 1.03032e+00 2.92154e-02 8.41380e-04 5.88413e-02
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 3 ERR DEF=0.5
1.373e-03 8.103e-04 -3.221e-04
8.103e-04 2.134e-03 -9.406e-04
-3.221e-04 -9.406e-04 8.536e-04
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2 3
1 0.47561 1.000 0.473 -0.298
2 0.75061 0.473 1.000 -0.697
3 0.69792 -0.298 -0.697 1.000
**********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 1500
**********
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=28280.5 FROM HESSE STATUS=OK 16 CALLS 80 TOTAL
EDM=7.12199e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 ml 4.88711e+00 3.71053e-02 3.90306e-05 2.46855e-01
2 sg 1.89676e+00 4.64029e-02 3.71490e-05 -6.90624e-01
3 sl 1.03032e+00 2.93441e-02 3.36552e-05 -9.47669e-01
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 3 ERR DEF=0.5
1.377e-03 8.197e-04 -3.287e-04
8.197e-04 2.153e-03 -9.535e-04
-3.287e-04 -9.535e-04 8.611e-04
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2 3
1 0.47808 1.000 0.476 -0.302
2 0.75324 0.476 1.000 -0.700
3 0.70113 -0.302 -0.700 1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(lx_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lxg) creating new cache 0x14a5d90 with pdf lx_CONV_gauss_CACHE_Obs[t] for nset (t) with code 0
[#1] INFO:NumericIntegration -- RooRealIntegral::init(lx_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooLandau.h"
#include "RooFFTConvPdf.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
using namespace RooFit ;
void rf208_convolution()
{
// S e t u p c o m p o n e n t p d f s
// ---------------------------------------
// Construct observable
RooRealVar t("t","t",-10,30) ;
// Construct landau(t,ml,sl) ;
RooRealVar ml("ml","mean landau",5.,-20,20) ;
RooRealVar sl("sl","sigma landau",1,0.1,10) ;
RooLandau landau("lx","lx",t,ml,sl) ;
// Construct gauss(t,mg,sg)
RooRealVar mg("mg","mg",0) ;
RooRealVar sg("sg","sg",2,0.1,10) ;
RooGaussian gauss("gauss","gauss",t,mg,sg) ;
// C o n s t r u c t c o n v o l u t i o n p d f
// ---------------------------------------
// Set #bins to be used for FFT sampling to 10000
t.setBins(10000,"cache") ;
// Construct landau (x) gauss
RooFFTConvPdf lxg("lxg","landau (X) gauss",t,landau,gauss) ;
// S a m p l e , f i t a n d p l o t c o n v o l u t e d p d f
// ----------------------------------------------------------------------
// Sample 1000 events in x from gxlx
RooDataSet* data = lxg.generate(t,10000) ;
// Fit gxlx to data
lxg.fitTo(*data) ;
// Plot data, landau pdf, landau (X) gauss pdf
RooPlot* frame = t.frame(Title("landau (x) gauss convolution")) ;
data->plotOn(frame) ;
lxg.plotOn(frame) ;
landau.plotOn(frame,LineStyle(kDashed)) ;
// Draw frame on canvas
new TCanvas("rf208_convolution","rf208_convolution",600,600) ;
gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
}
Author
07/2008 - Wouter Verkerke

Definition in file rf208_convolution.C.