Logo ROOT  
Reference Guide
 
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 transformed 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_image
15## \macro_code
16## \macro_output
17##
18## \date February 2018
19## \authors Clemens Lange, Wouter Verkerke (C version)
20
21import ROOT
22
23# Set up component pdfs
24# ---------------------------------------
25
26# Define angle psi
27psi = ROOT.RooRealVar("psi", "psi", 0, 3.14159268)
28
29# Define physics p.d.f T(psi)
30Tpsi = ROOT.RooGenericPdf("Tpsi", "1+sin(2*@0)", [psi])
31
32# Define resolution R(psi)
33gbias = ROOT.RooRealVar("gbias", "gbias", 0.2, 0.0, 1)
34greso = ROOT.RooRealVar("greso", "greso", 0.3, 0.1, 1.0)
35Rpsi = ROOT.RooGaussian("Rpsi", "Rpsi", psi, gbias, greso)
36
37# Define cos(psi) and function psif that calculates psi from cos(psi)
38cpsi = ROOT.RooRealVar("cpsi", "cos(psi)", -1, 1)
39psif = ROOT.RooFormulaVar("psif", "acos(cpsi)", [cpsi])
40
41# Define physics p.d.f. also as function of cos(psi): T(psif(cpsi)) = T(cpsi)
42Tcpsi = ROOT.RooGenericPdf("T", "1+sin(2*@0)", [psif])
43
44# Construct convolution pdf in psi
45# --------------------------------------------------------------
46
47# Define convoluted p.d.f. as function of psi: M=[T(x)R](psi) = M(psi)
48Mpsi = ROOT.RooFFTConvPdf("Mf", "Mf", psi, Tpsi, Rpsi)
49
50# Set the buffer fraction to zero to obtain a ROOT.True cyclical
51# convolution
52Mpsi.setBufferFraction(0)
53
54# Sample, fit and plot convoluted pdf (psi)
55# --------------------------------------------------------------------------------
56
57# Generate some events in observable psi
58data_psi = Mpsi.generate({psi}, 10000)
59
60# Fit convoluted model as function of angle psi
61Mpsi.fitTo(data_psi, PrintLevel=-1)
62
63# Plot cos(psi) frame with Mf(cpsi)
64frame1 = psi.frame(Title="Cyclical convolution in angle psi")
65data_psi.plotOn(frame1)
66Mpsi.plotOn(frame1)
67
68# Overlay comparison to unsmeared physics p.d.f ROOT.T(psi)
69Tpsi.plotOn(frame1, LineColor="r")
70
71# Construct convolution pdf in cos(psi)
72# --------------------------------------------------------------------------
73
74# Define convoluted p.d.f. as function of cos(psi): M=[T(x)R](psif cpsi)) = M(cpsi:
75#
76# Need to give both observable psi here (for definition of convolution)
77# and function psif here (for definition of observables, in cpsi)
78Mcpsi = ROOT.RooFFTConvPdf("Mf", "Mf", psif, psi, Tpsi, Rpsi)
79
80# Set the buffer fraction to zero to obtain a ROOT.True cyclical
81# convolution
82Mcpsi.setBufferFraction(0)
83
84# Sample, fit and plot convoluted pdf (cospsi)
85# --------------------------------------------------------------------------------
86
87# Generate some events
88data_cpsi = Mcpsi.generate({cpsi}, 10000)
89
90# set psi constant to exclude to be a parameter of the fit
91psi.setConstant(True)
92
93# Fit convoluted model as function of cos(psi)
94Mcpsi.fitTo(data_cpsi, PrintLevel=-1)
95
96# Plot cos(psi) frame with Mf(cpsi)
97frame2 = cpsi.frame(Title="Same convolution in psi, in cos(psi)")
98data_cpsi.plotOn(frame2)
99Mcpsi.plotOn(frame2)
100
101# Overlay comparison to unsmeared physics p.d.f ROOT.Tf(cpsi)
102Tcpsi.plotOn(frame2, LineColor="r")
103
104# Draw frame on canvas
105c = ROOT.TCanvas("rf210_angularconv", "rf210_angularconv", 800, 400)
106c.Divide(2)
107c.cd(1)
108ROOT.gPad.SetLeftMargin(0.15)
109frame1.GetYaxis().SetTitleOffset(1.4)
110frame1.Draw()
111c.cd(2)
112ROOT.gPad.SetLeftMargin(0.15)
113frame2.GetYaxis().SetTitleOffset(1.4)
114frame2.Draw()
115
116c.SaveAs("rf210_angularconv.png")