Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf704_amplitudefit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Special pdf's: using a pdf defined by a sum of real-valued amplitude components
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13
14# Setup 2D amplitude functions
15# -------------------------------------------------------
16
17# Observables
18t = ROOT.RooRealVar("t", "time", -1.0, 15.0)
19cosa = ROOT.RooRealVar("cosa", "cos(alpha)", -1.0, 1.0)
20
21# Use ROOT.RooTruthModel to obtain compiled implementation of sinh/cosh
22# modulated decay functions
23tau = ROOT.RooRealVar("tau", "#tau", 1.5)
24deltaGamma = ROOT.RooRealVar("deltaGamma", "deltaGamma", 0.3)
25tm = ROOT.RooTruthModel("tm", "tm", t)
26coshGBasis = ROOT.RooFormulaVar("coshGBasis", "exp(-@0/ @1)*cosh(@0*@2/2)", [t, tau, deltaGamma])
27sinhGBasis = ROOT.RooFormulaVar("sinhGBasis", "exp(-@0/ @1)*sinh(@0*@2/2)", [t, tau, deltaGamma])
28coshGConv = tm.convolution(coshGBasis, t)
29sinhGConv = tm.convolution(sinhGBasis, t)
30
31# Construct polynomial amplitudes in cos(a)
32poly1 = ROOT.RooPolyVar("poly1", "poly1", cosa, [0.5, 0.2, 0.2], 0)
33poly2 = ROOT.RooPolyVar("poly2", "poly2", cosa, [1.0, -0.2, 3.0], 0)
34
35# Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
36ampl1 = ROOT.RooProduct("ampl1", "amplitude 1", [poly1, coshGConv])
37ampl2 = ROOT.RooProduct("ampl2", "amplitude 2", [poly2, sinhGConv])
38
39# Construct amplitude sum pdf
40# -----------------------------------------------------
41
42# Amplitude strengths
43f1 = ROOT.RooRealVar("f1", "f1", 1, 0, 2)
44f2 = ROOT.RooRealVar("f2", "f2", 0.5, 0, 2)
45
46# Construct pdf
47pdf = ROOT.RooRealSumPdf("pdf", "pdf", [ampl1, ampl2], [f1, f2])
48
49# Generate some toy data from pdf
50data = pdf.generate({t, cosa}, 10000)
51
52# Fit pdf to toy data with only amplitude strength floating
53pdf.fitTo(data)
54
55# Plot amplitude sum pdf
56# -------------------------------------------
57
58# Make 2D plots of amplitudes
59hh_cos = ampl1.createHistogram("hh_cos", t, Binning=50, YVar=dict(var=cosa, Binning=50))
60hh_sin = ampl2.createHistogram("hh_sin", t, Binning=50, YVar=dict(var=cosa, Binning=50))
61hh_cos.SetLineColor(ROOT.kBlue)
62hh_sin.SetLineColor(ROOT.kRed)
63
64# Make projection on t, data, and its components
65# Note component projections may be larger than sum because amplitudes can
66# be negative
67frame1 = t.frame()
68data.plotOn(frame1)
69pdf.plotOn(frame1)
70pdf.plotOn(frame1, Components=ampl1, LineStyle="--")
71pdf.plotOn(frame1, Components=ampl2, LineStyle="--", LineColor="r")
72
73# Make projection on cosa, data, and its components
74# Note that components projection may be larger than sum because
75# amplitudes can be negative
76frame2 = cosa.frame()
77data.plotOn(frame2)
78pdf.plotOn(frame2)
79pdf.plotOn(frame2, Components=ampl1, LineStyle="--")
80pdf.plotOn(frame2, Components=ampl2, LineStyle="--", LineColor="r")
81
82c = ROOT.TCanvas("rf704_amplitudefit", "rf704_amplitudefit", 800, 800)
83c.Divide(2, 2)
84c.cd(1)
85ROOT.gPad.SetLeftMargin(0.15)
86frame1.GetYaxis().SetTitleOffset(1.4)
87frame1.Draw()
88c.cd(2)
89ROOT.gPad.SetLeftMargin(0.15)
90frame2.GetYaxis().SetTitleOffset(1.4)
91frame2.Draw()
92c.cd(3)
93ROOT.gPad.SetLeftMargin(0.20)
94hh_cos.GetZaxis().SetTitleOffset(2.3)
95hh_cos.Draw("surf")
96c.cd(4)
97ROOT.gPad.SetLeftMargin(0.20)
98hh_sin.GetZaxis().SetTitleOffset(2.3)
99hh_sin.Draw("surf")
100
101c.SaveAs("rf704_amplitudefit.png")