Logo ROOT  
Reference Guide
rf208_convolution.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Addition and convolution: one-dimensional numeric convolution
5///
6/// pdf = landau(t) (x) gauss(t)
7///
8/// This tutorial requires FFT3 to be enabled.
9///
10/// \macro_image
11/// \macro_output
12/// \macro_code
13/// \author 07/2008 - Wouter Verkerke
14
15#include "RooRealVar.h"
16#include "RooDataSet.h"
17#include "RooGaussian.h"
18#include "RooLandau.h"
19#include "RooFFTConvPdf.h"
20#include "RooPlot.h"
21#include "TCanvas.h"
22#include "TAxis.h"
23#include "TH1.h"
24using namespace RooFit;
25
26void rf208_convolution()
27{
28 // S e t u p c o m p o n e n t p d f s
29 // ---------------------------------------
30
31 // Construct observable
32 RooRealVar t("t", "t", -10, 30);
33
34 // Construct landau(t,ml,sl) ;
35 RooRealVar ml("ml", "mean landau", 5., -20, 20);
36 RooRealVar sl("sl", "sigma landau", 1, 0.1, 10);
37 RooLandau landau("lx", "lx", t, ml, sl);
38
39 // Construct gauss(t,mg,sg)
40 RooRealVar mg("mg", "mg", 0);
41 RooRealVar sg("sg", "sg", 2, 0.1, 10);
42 RooGaussian gauss("gauss", "gauss", t, mg, sg);
43
44 // C o n s t r u c t c o n v o l u t i o n p d f
45 // ---------------------------------------
46
47 // Set #bins to be used for FFT sampling to 10000
48 t.setBins(10000, "cache");
49
50 // Construct landau (x) gauss
51 RooFFTConvPdf lxg("lxg", "landau (X) gauss", t, landau, gauss);
52
53 // 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
54 // ----------------------------------------------------------------------
55
56 // Sample 1000 events in x from gxlx
57 RooDataSet *data = lxg.generate(t, 10000);
58
59 // Fit gxlx to data
60 lxg.fitTo(*data);
61
62 // Plot data, landau pdf, landau (X) gauss pdf
63 RooPlot *frame = t.frame(Title("landau (x) gauss convolution"));
64 data->plotOn(frame);
65 lxg.plotOn(frame);
66 landau.plotOn(frame, LineStyle(kDashed));
67
68 // Draw frame on canvas
69 new TCanvas("rf208_convolution", "rf208_convolution", 600, 600);
70 gPad->SetLeftMargin(0.15);
71 frame->GetYaxis()->SetTitleOffset(1.4);
72 frame->Draw();
73}
@ kDashed
Definition: TAttLine.h:48
#define gPad
Definition: TVirtualPad.h:286
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:550
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
PDF for the numerical (FFT) convolution of two PDFs.
Definition: RooFFTConvPdf.h:26
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
Landau distribution p.d.f.
Definition: RooLandau.h:24
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition: RooPlot.cxx:277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
The Canvas class.
Definition: TCanvas.h:31
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooCmdArg LineStyle(Style_t style)
static constexpr double gauss
static constexpr double mg
const char * Title
Definition: TXMLSetup.cxx:67