ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf208_convolution.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// 'ADDITION AND CONVOLUTION' RooFit tutorial macro #208
4 ///
5 /// One-dimensional numeric convolution
6 /// (require ROOT to be compiled with --enable-fftw3)
7 ///
8 /// pdf = landau(t) (x) gauss(t)
9 ///
10 ///
11 /// \macro_code
12 /// \author 07/2008 - Wouter Verkerke
13 
14 
15 #ifndef __CINT__
16 #include "RooGlobalFunc.h"
17 #endif
18 #include "RooRealVar.h"
19 #include "RooDataSet.h"
20 #include "RooGaussian.h"
21 #include "RooLandau.h"
22 #include "RooFFTConvPdf.h"
23 #include "RooPlot.h"
24 #include "TCanvas.h"
25 #include "TAxis.h"
26 #include "TH1.h"
27 using namespace RooFit ;
28 
29 
30 
31 void rf208_convolution()
32 {
33  // S e t u p c o m p o n e n t p d f s
34  // ---------------------------------------
35 
36  // Construct observable
37  RooRealVar t("t","t",-10,30) ;
38 
39  // Construct landau(t,ml,sl) ;
40  RooRealVar ml("ml","mean landau",5.,-20,20) ;
41  RooRealVar sl("sl","sigma landau",1,0.1,10) ;
42  RooLandau landau("lx","lx",t,ml,sl) ;
43 
44  // Construct gauss(t,mg,sg)
45  RooRealVar mg("mg","mg",0) ;
46  RooRealVar sg("sg","sg",2,0.1,10) ;
47  RooGaussian gauss("gauss","gauss",t,mg,sg) ;
48 
49 
50  // C o n s t r u c t c o n v o l u t i o n p d f
51  // ---------------------------------------
52 
53  // Set #bins to be used for FFT sampling to 10000
54  t.setBins(10000,"cache") ;
55 
56  // Construct landau (x) gauss
57  RooFFTConvPdf lxg("lxg","landau (X) gauss",t,landau,gauss) ;
58 
59 
60 
61  // 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
62  // ----------------------------------------------------------------------
63 
64  // Sample 1000 events in x from gxlx
65  RooDataSet* data = lxg.generate(t,10000) ;
66 
67  // Fit gxlx to data
68  lxg.fitTo(*data) ;
69 
70  // Plot data, landau pdf, landau (X) gauss pdf
71  RooPlot* frame = t.frame(Title("landau (x) gauss convolution")) ;
72  data->plotOn(frame) ;
73  lxg.plotOn(frame) ;
74  landau.plotOn(frame,LineStyle(kDashed)) ;
75 
76 
77  // Draw frame on canvas
78  new TCanvas("rf208_convolution","rf208_convolution",600,600) ;
79  gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
80 
81 }
82 
83 
84 
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title Offset is a correction factor with respect to the "s...
Definition: TAttAxis.cxx:245
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1118
Landau Distribution p.d.f.
Definition: RooLandau.h:24
RooCmdArg Title(const char *name)
RooCmdArg LineStyle(Style_t style)
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
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
Plot dataset on specified frame.
Definition: RooAbsData.cxx:626
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
TThread * t[5]
Definition: threadsh1.C:13
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
#define gPad
Definition: TVirtualPad.h:288
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:559