Logo ROOT   6.16/01
Reference Guide
rf209_anaconv.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'ADDITION AND CONVOLUTION' RooFit tutorial macro #209
5## Decay function p.d.fs with optional B physics
6## effects (mixing and CP violation) that can be
7## analytically convolved with e.g. Gaussian resolution
8## functions
9##
10## pdf1 = decay(t,tau) (x) delta(t)
11## pdf2 = decay(t,tau) (x) gauss(t,m,s)
12## pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
13##
14## \macro_code
15##
16## \date February 2018
17## \author Clemens Lange
18## \author Wouter Verkerke (C version)
19
20import ROOT
21
22# B-physics pdf with truth resolution
23# ---------------------------------------------------------------------
24
25# Variables of decay p.d.f.
26dt = ROOT.RooRealVar("dt", "dt", -10, 10)
27tau = ROOT.RooRealVar("tau", "tau", 1.548)
28
29# Build a truth resolution model (delta function)
30tm = ROOT.RooTruthModel("tm", "truth model", dt)
31
32# Construct decay(t) (x) delta(t)
33decay_tm = ROOT.RooDecay("decay_tm", "decay", dt,
34 tau, tm, ROOT.RooDecay.DoubleSided)
35
36# Plot p.d.f. (dashed)
37frame = dt.frame(ROOT.RooFit.Title("Bdecay (x) resolution"))
38decay_tm.plotOn(frame, ROOT.RooFit.LineStyle(ROOT.kDashed))
39
40# B-physics pdf with Gaussian resolution
41# ----------------------------------------------------------------------------
42
43# Build a gaussian resolution model
44bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
45sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 1)
46gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
47
48# Construct decay(t) (x) gauss1(t)
49decay_gm1 = ROOT.RooDecay("decay_gm1", "decay",
50 dt, tau, gm1, ROOT.RooDecay.DoubleSided)
51
52# Plot p.d.f.
53decay_gm1.plotOn(frame)
54
55# B-physics pdf with double Gaussian resolution
56# ------------------------------------------------------------------------------------------
57
58# Build another gaussian resolution model
59bias2 = ROOT.RooRealVar("bias2", "bias2", 0)
60sigma2 = ROOT.RooRealVar("sigma2", "sigma2", 5)
61gm2 = ROOT.RooGaussModel("gm2", "gauss model 2", dt, bias2, sigma2)
62
63# Build a composite resolution model f*gm1+(1-f)*gm2
64gm1frac = ROOT.RooRealVar("gm1frac", "fraction of gm1", 0.5)
65gmsum = ROOT.RooAddModel(
66 "gmsum", "sum of gm1 and gm2", ROOT.RooArgList(gm1, gm2), ROOT.RooArgList(gm1frac))
67
68# Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
69decay_gmsum = ROOT.RooDecay(
70 "decay_gmsum", "decay", dt, tau, gmsum, ROOT.RooDecay.DoubleSided)
71
72# Plot p.d.f. (red)
73decay_gmsum.plotOn(frame, ROOT.RooFit.LineColor(ROOT.kRed))
74
75# Draw all frames on canvas
76c = ROOT.TCanvas("rf209_anaconv", "rf209_anaconv", 600, 600)
77ROOT.gPad.SetLeftMargin(0.15)
78frame.GetYaxis().SetTitleOffset(1.6)
79frame.Draw()
80
81c.SaveAs("rf209_anaconv.png")