Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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"
29 using namespace RooFit;
30 
31 void 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 }