Logo ROOT  
Reference Guide
rf210_angularconv.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Addition and convolution: convolution in cyclical angular observables theta
5///
6/// and construction of p.d.f in terms of transformed angular coordinates, e.g. cos(theta),
7/// where the convolution is performed in theta rather than cos(theta)
8///
9/// pdf(theta) = T(theta) (x) gauss(theta)
10/// pdf(cosTheta) = T(acos(cosTheta)) (x) gauss(acos(cosTheta))
11///
12/// This tutorial requires FFT3 to be enabled.
13///
14/// \macro_image
15/// \macro_output
16/// \macro_code
17/// \author 04/2009 - Wouter Verkerke
18
19#include "RooRealVar.h"
20#include "RooDataSet.h"
21#include "RooGaussian.h"
22#include "RooGenericPdf.h"
23#include "RooFormulaVar.h"
24#include "RooFFTConvPdf.h"
25#include "RooPlot.h"
26#include "TCanvas.h"
27#include "TAxis.h"
28#include "TH1.h"
29using namespace RooFit;
30
31void rf210_angularconv()
32{
33 // S e t u p c o m p o n e n t p d f s
34 // ---------------------------------------
35
36 // Define angle psi
37 RooRealVar psi("psi", "psi", 0, 3.14159268);
38
39 // Define physics p.d.f T(psi)
40 RooGenericPdf Tpsi("Tpsi", "1+sin(2*@0)", psi);
41
42 // Define resolution R(psi)
43 RooRealVar gbias("gbias", "gbias", 0.2, 0., 1);
44 RooRealVar greso("greso", "greso", 0.3, 0.1, 1.0);
45 RooGaussian Rpsi("Rpsi", "Rpsi", psi, gbias, greso);
46
47 // Define cos(psi) and function psif that calculates psi from cos(psi)
48 RooRealVar cpsi("cpsi", "cos(psi)", -1, 1);
49 RooFormulaVar psif("psif", "acos(cpsi)", cpsi);
50
51 // Define physics p.d.f. also as function of cos(psi): T(psif(cpsi)) = T(cpsi) ;
52 RooGenericPdf Tcpsi("T", "1+sin(2*@0)", psif);
53
54 // C o n s t r u c t c o n v o l u t i o n p d f i n p s i
55 // --------------------------------------------------------------
56
57 // Define convoluted p.d.f. as function of psi: M=[T(x)R](psi) = M(psi)
58 RooFFTConvPdf Mpsi("Mf", "Mf", psi, Tpsi, Rpsi);
59
60 // Set the buffer fraction to zero to obtain a true cyclical convolution
61 Mpsi.setBufferFraction(0);
62
63 // 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 ( p s i )
64 // --------------------------------------------------------------------------------
65
66 // Generate some events in observable psi
67 RooDataSet *data_psi = Mpsi.generate(psi, 10000);
68
69 // Fit convoluted model as function of angle psi
70 Mpsi.fitTo(*data_psi);
71
72 // Plot cos(psi) frame with Mf(cpsi)
73 RooPlot *frame1 = psi.frame(Title("Cyclical convolution in angle psi"));
74 data_psi->plotOn(frame1);
75 Mpsi.plotOn(frame1);
76
77 // Overlay comparison to unsmeared physics p.d.f T(psi)
78 Tpsi.plotOn(frame1, LineColor(kRed));
79
80 // C o n s t r u c t c o n v o l u t i o n p d f i n c o s ( p s i )
81 // --------------------------------------------------------------------------
82
83 // Define convoluted p.d.f. as function of cos(psi): M=[T(x)R](psif(cpsi)) = M(cpsi)
84 //
85 // Need to give both observable psi here (for definition of convolution)
86 // and function psif here (for definition of observables, ultimately in cpsi)
87 RooFFTConvPdf Mcpsi("Mf", "Mf", psif, psi, Tpsi, Rpsi);
88
89 // Set the buffer fraction to zero to obtain a true cyclical convolution
90 Mcpsi.setBufferFraction(0);
91
92 // 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 ( c o s p s i )
93 // --------------------------------------------------------------------------------
94
95 // Generate some events
96 RooDataSet *data_cpsi = Mcpsi.generate(cpsi, 10000);
97
98 // set psi constant to exclude to be a parameter of the fit
99 psi.setConstant(true);
100
101 // Fit convoluted model as function of cos(psi)
102 Mcpsi.fitTo(*data_cpsi);
103
104 // Plot cos(psi) frame with Mf(cpsi)
105 RooPlot *frame2 = cpsi.frame(Title("Same convolution in psi, expressed in cos(psi)"));
106 data_cpsi->plotOn(frame2);
107 Mcpsi.plotOn(frame2);
108
109 // Overlay comparison to unsmeared physics p.d.f Tf(cpsi)
110 Tcpsi.plotOn(frame2, LineColor(kRed));
111
112 // Draw frame on canvas
113 TCanvas *c = new TCanvas("rf210_angularconv", "rf210_angularconv", 800, 400);
114 c->Divide(2);
115 c->cd(1);
116 gPad->SetLeftMargin(0.15);
117 frame1->GetYaxis()->SetTitleOffset(1.4);
118 frame1->Draw();
119 c->cd(2);
120 gPad->SetLeftMargin(0.15);
121 frame2->GetYaxis()->SetTitleOffset(1.4);
122 frame2->Draw();
123}
#define c(i)
Definition: RSha256.hxx:101
@ kRed
Definition: Rtypes.h:64
#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
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Definition: RooFormulaVar.h:29
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooGenericPdf is a concrete implementation of a probability density function, which takes a RooArgLis...
Definition: RooGenericPdf.h:25
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 LineColor(Color_t color)
const char * Title
Definition: TXMLSetup.cxx:67