Logo ROOT   6.16/01
Reference Guide
rf704_amplitudefit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'SPECIAL PDFS' RooFit tutorial macro #704
6##
7## Using a p.d.f defined by a sum of real-valued amplitude components
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15
16import ROOT
17
18
19# Setup 2D amplitude functions
20# -------------------------------------------------------
21
22# Observables
23t = ROOT.RooRealVar("t", "time", -1., 15.)
24cosa = ROOT.RooRealVar("cosa", "cos(alpha)", -1., 1.)
25
26# Use ROOT.RooTruthModel to obtain compiled implementation of sinh/cosh
27# modulated decay functions
28tau = ROOT.RooRealVar("tau", "#tau", 1.5)
29deltaGamma = ROOT.RooRealVar("deltaGamma", "deltaGamma", 0.3)
30tm = ROOT.RooTruthModel("tm", "tm", t)
31coshGBasis = ROOT.RooFormulaVar(
32 "coshGBasis", "exp(-@0/ @1)*cosh(@0*@2/2)", ROOT.RooArgList(t, tau, deltaGamma))
33sinhGBasis = ROOT.RooFormulaVar(
34 "sinhGBasis", "exp(-@0/ @1)*sinh(@0*@2/2)", ROOT.RooArgList(t, tau, deltaGamma))
35coshGConv = tm.convolution(coshGBasis, t)
36sinhGConv = tm.convolution(sinhGBasis, t)
37
38# Construct polynomial amplitudes in cos(a)
39poly1 = ROOT.RooPolyVar("poly1", "poly1", cosa, ROOT.RooArgList(
40 ROOT.RooFit.RooConst(0.5), ROOT.RooFit.RooConst(0.2), ROOT.RooFit.RooConst(0.2)), 0)
41poly2 = ROOT.RooPolyVar("poly2", "poly2", cosa, ROOT.RooArgList(
42 ROOT.RooFit.RooConst(1), ROOT.RooFit.RooConst(-0.2), ROOT.RooFit.RooConst(3)), 0)
43
44# Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
45ampl1 = ROOT.RooProduct("ampl1", "amplitude 1",
46 ROOT.RooArgList(poly1, coshGConv))
47ampl2 = ROOT.RooProduct("ampl2", "amplitude 2",
48 ROOT.RooArgList(poly2, sinhGConv))
49
50# Construct amplitude sum pdf
51# -----------------------------------------------------
52
53# Amplitude strengths
54f1 = ROOT.RooRealVar("f1", "f1", 1, 0, 2)
55f2 = ROOT.RooRealVar("f2", "f2", 0.5, 0, 2)
56
57# Construct pdf
58pdf = ROOT.RooRealSumPdf("pdf", "pdf", ROOT.RooArgList(
59 ampl1, ampl2), ROOT.RooArgList(f1, f2))
60
61# Generate some toy data from pdf
62data = pdf.generate(ROOT.RooArgSet(t, cosa), 10000)
63
64# Fit pdf to toy data with only amplitude strength floating
65pdf.fitTo(data)
66
67# Plot amplitude sum pdf
68# -------------------------------------------
69
70# Make 2D plots of amplitudes
71hh_cos = ampl1.createHistogram("hh_cos", t, ROOT.RooFit.Binning(
72 50), ROOT.RooFit.YVar(cosa, ROOT.RooFit.Binning(50)))
73hh_sin = ampl2.createHistogram("hh_sin", t, ROOT.RooFit.Binning(
74 50), ROOT.RooFit.YVar(cosa, ROOT.RooFit.Binning(50)))
75hh_cos.SetLineColor(ROOT.kBlue)
76hh_sin.SetLineColor(ROOT.kRed)
77
78# Make projection on t, data, and its components
79# Note component projections may be larger than sum because amplitudes can
80# be negative
81frame1 = t.frame()
82data.plotOn(frame1)
83pdf.plotOn(frame1)
84# workaround, see https://root.cern.ch/phpBB3/viewtopic.php?t=7764
85ras_ampl1 = ROOT.RooArgSet(ampl1)
86pdf.plotOn(frame1, ROOT.RooFit.Components(ras_ampl1),
87 ROOT.RooFit.LineStyle(ROOT.kDashed))
88ras_ampl2 = ROOT.RooArgSet(ampl2)
89pdf.plotOn(frame1, ROOT.RooFit.Components(ras_ampl2), ROOT.RooFit.LineStyle(
90 ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
91
92# Make projection on cosa, data, and its components
93# Note that components projection may be larger than sum because
94# amplitudes can be negative
95frame2 = cosa.frame()
96data.plotOn(frame2)
97pdf.plotOn(frame2)
98pdf.plotOn(frame2, ROOT.RooFit.Components(ras_ampl1),
99 ROOT.RooFit.LineStyle(ROOT.kDashed))
100pdf.plotOn(frame2, ROOT.RooFit.Components(ras_ampl2), ROOT.RooFit.LineStyle(
101 ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
102
103c = ROOT.TCanvas("rf704_amplitudefit", "rf704_amplitudefit", 800, 800)
104c.Divide(2, 2)
105c.cd(1)
106ROOT.gPad.SetLeftMargin(0.15)
107frame1.GetYaxis().SetTitleOffset(1.4)
108frame1.Draw()
109c.cd(2)
110ROOT.gPad.SetLeftMargin(0.15)
111frame2.GetYaxis().SetTitleOffset(1.4)
112frame2.Draw()
113c.cd(3)
114ROOT.gPad.SetLeftMargin(0.20)
115hh_cos.GetZaxis().SetTitleOffset(2.3)
116hh_cos.Draw("surf")
117c.cd(4)
118ROOT.gPad.SetLeftMargin(0.20)
119hh_sin.GetZaxis().SetTitleOffset(2.3)
120hh_sin.Draw("surf")
121
122c.SaveAs("rf704_amplitudefit.png")