Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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/// ```
7/// pdf = landau(t) (x) gauss(t)
8/// ```
9///
10/// This tutorial requires FFT3 to be enabled.
11///
12/// \macro_image
13/// \macro_code
14/// \macro_output
15///
16/// \date July 2008
17/// \author Wouter Verkerke
18
19#include "RooRealVar.h"
20#include "RooDataSet.h"
21#include "RooGaussian.h"
22#include "RooLandau.h"
23#include "RooFFTConvPdf.h"
24#include "RooPlot.h"
25#include "TCanvas.h"
26#include "TAxis.h"
27#include "TH1.h"
28using namespace RooFit;
29
31{
32 // S e t u p c o m p o n e n t p d f s
33 // ---------------------------------------
34
35 // Construct observable
36 RooRealVar t("t", "t", -10, 30);
37
38 // Construct landau(t,ml,sl) ;
39 RooRealVar ml("ml", "mean landau", 5., -20, 20);
40 RooRealVar sl("sl", "sigma landau", 1, 0.1, 10);
41 RooLandau landau("lx", "lx", t, ml, sl);
42
43 // Construct gauss(t,mg,sg)
44 RooRealVar mg("mg", "mg", 0);
45 RooRealVar sg("sg", "sg", 2, 0.1, 10);
46 RooGaussian gauss("gauss", "gauss", t, mg, sg);
47
48 // C o n s t r u c t c o n v o l u t i o n p d f
49 // ---------------------------------------
50
51 // Set #bins to be used for FFT sampling to 10000
52 t.setBins(10000, "cache");
53
54 // Construct landau (x) gauss
55 RooFFTConvPdf lxg("lxg", "landau (X) gauss", t, landau, gauss);
56
57 // 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
58 // ----------------------------------------------------------------------
59
60 // Sample 1000 events in x from gxlx
61 std::unique_ptr<RooDataSet> data{lxg.generate(t, 10000)};
62
63 // Fit gxlx to data
64 lxg.fitTo(*data, PrintLevel(-1));
65
66 // Plot data, landau pdf, landau (X) gauss pdf
67 RooPlot *frame = t.frame(Title("landau (x) gauss convolution"));
68 data->plotOn(frame);
69 lxg.plotOn(frame);
70 landau.plotOn(frame, LineStyle(kDashed));
71
72 // Draw frame on canvas
73 new TCanvas("rf208_convolution", "rf208_convolution", 600, 600);
74 gPad->SetLeftMargin(0.15);
75 frame->GetYaxis()->SetTitleOffset(1.4);
76 frame->Draw();
77}
@ kDashed
Definition TAttLine.h:48
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
PDF for the numerical (FFT) convolution of two PDFs.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
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:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:239
TAxis * GetYaxis() const
Definition RooPlot.cxx:1279
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:652
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg PrintLevel(Int_t code)
RooCmdArg LineStyle(Style_t style)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68