Loading [MathJax]/extensions/tex2jax.js
Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
rf210_angularconv.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Convolution in cyclical angular observables theta, and
5## construction of p.d.f in terms of trasnformed angular
6## coordinates, e.g. cos(theta), the convolution
7## is performed in theta rather than cos(theta)
8##
9## (require ROOT to be compiled with --enable-fftw3)
10##
11## pdf(theta) = ROOT.T(theta) (x) gauss(theta)
12## pdf(cosTheta) = ROOT.T(acos(cosTheta)) (x) gauss(acos(cosTheta))
13##
14## \macro_code
15##
16## \date February 2018
17## \author Clemens Lange
18## \author Wouter Verkerke (C version)
19
20import ROOT
21
22# Set up component pdfs
23# ---------------------------------------
24
25# Define angle psi
26psi = ROOT.RooRealVar("psi", "psi", 0, 3.14159268)
27
28# Define physics p.d.f T(psi)
29Tpsi = ROOT.RooGenericPdf("Tpsi", "1+sin(2*@0)", [psi])
30
31# Define resolution R(psi)
32gbias = ROOT.RooRealVar("gbias", "gbias", 0.2, 0.0, 1)
33greso = ROOT.RooRealVar("greso", "greso", 0.3, 0.1, 1.0)
34Rpsi = ROOT.RooGaussian("Rpsi", "Rpsi", psi, gbias, greso)
35
36# Define cos(psi) and function psif that calculates psi from cos(psi)
37cpsi = ROOT.RooRealVar("cpsi", "cos(psi)", -1, 1)
38psif = ROOT.RooFormulaVar("psif", "acos(cpsi)", [cpsi])
39
40# Define physics p.d.f. also as function of cos(psi): T(psif(cpsi)) = T(cpsi)
41Tcpsi = ROOT.RooGenericPdf("T", "1+sin(2*@0)", [psif])
42
43# Construct convolution pdf in psi
44# --------------------------------------------------------------
45
46# Define convoluted p.d.f. as function of psi: M=[T(x)R](psi) = M(psi)
47Mpsi = ROOT.RooFFTConvPdf("Mf", "Mf", psi, Tpsi, Rpsi)
48
49# Set the buffer fraction to zero to obtain a ROOT.True cyclical
50# convolution
51Mpsi.setBufferFraction(0)
52
53# Sample, fit and plot convoluted pdf (psi)
54# --------------------------------------------------------------------------------
55
56# Generate some events in observable psi
57data_psi = Mpsi.generate({psi}, 10000)
58
59# Fit convoluted model as function of angle psi
60Mpsi.fitTo(data_psi)
61
62# Plot cos(psi) frame with Mf(cpsi)
63frame1 = psi.frame(Title="Cyclical convolution in angle psi")
64data_psi.plotOn(frame1)
65Mpsi.plotOn(frame1)
66
67# Overlay comparison to unsmeared physics p.d.f ROOT.T(psi)
68Tpsi.plotOn(frame1, LineColor="r")
69
70# Construct convolution pdf in cos(psi)
71# --------------------------------------------------------------------------
72
73# Define convoluted p.d.f. as function of cos(psi): M=[T(x)R](psif cpsi)) = M(cpsi:
74#
75# Need to give both observable psi here (for definition of convolution)
76# and function psif here (for definition of observables, in cpsi)
77Mcpsi = ROOT.RooFFTConvPdf("Mf", "Mf", psif, psi, Tpsi, Rpsi)
78
79# Set the buffer fraction to zero to obtain a ROOT.True cyclical
80# convolution
81Mcpsi.setBufferFraction(0)
82
83# Sample, fit and plot convoluted pdf (cospsi)
84# --------------------------------------------------------------------------------
85
86# Generate some events
87data_cpsi = Mcpsi.generate({cpsi}, 10000)
88
89# set psi constant to exclude to be a parameter of the fit
90psi.setConstant(True)
91
92# Fit convoluted model as function of cos(psi)
93Mcpsi.fitTo(data_cpsi)
94
95# Plot cos(psi) frame with Mf(cpsi)
96frame2 = cpsi.frame(Title="Same convolution in psi, in cos(psi)")
97data_cpsi.plotOn(frame2)
98Mcpsi.plotOn(frame2)
99
100# Overlay comparison to unsmeared physics p.d.f ROOT.Tf(cpsi)
101Tcpsi.plotOn(frame2, LineColor="r")
102
103# Draw frame on canvas
104c = ROOT.TCanvas("rf210_angularconv", "rf210_angularconv", 800, 400)
105c.Divide(2)
106c.cd(1)
107ROOT.gPad.SetLeftMargin(0.15)
108frame1.GetYaxis().SetTitleOffset(1.4)
109frame1.Draw()
110c.cd(2)
111ROOT.gPad.SetLeftMargin(0.15)
112frame2.GetYaxis().SetTitleOffset(1.4)
113frame2.Draw()
114
115c.SaveAs("rf210_angularconv.png")