Logo ROOT   6.18/05
Reference Guide
rf209_anaconv.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Addition and convolution: decay function p.d.fs with optional B physics effects (mixing and CP violation) that can be analytically convolved with e.g. Gaussian resolution functions
5##
6## pdf1 = decay(t,tau) (x) delta(t)
7## pdf2 = decay(t,tau) (x) gauss(t,m,s)
8## pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
9##
10## \macro_code
11##
12## \date February 2018
13## \author Clemens Lange, Wouter Verkerke (C++ version)
14
15import ROOT
16
17# B-physics pdf with truth resolution
18# ---------------------------------------------------------------------
19
20# Variables of decay p.d.f.
21dt = ROOT.RooRealVar("dt", "dt", -10, 10)
22tau = ROOT.RooRealVar("tau", "tau", 1.548)
23
24# Build a truth resolution model (delta function)
25tm = ROOT.RooTruthModel("tm", "truth model", dt)
26
27# Construct decay(t) (x) delta(t)
28decay_tm = ROOT.RooDecay("decay_tm", "decay", dt,
29 tau, tm, ROOT.RooDecay.DoubleSided)
30
31# Plot p.d.f. (dashed)
32frame = dt.frame(ROOT.RooFit.Title("Bdecay (x) resolution"))
33decay_tm.plotOn(frame, ROOT.RooFit.LineStyle(ROOT.kDashed))
34
35# B-physics pdf with Gaussian resolution
36# ----------------------------------------------------------------------------
37
38# Build a gaussian resolution model
39bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
40sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 1)
41gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
42
43# Construct decay(t) (x) gauss1(t)
44decay_gm1 = ROOT.RooDecay("decay_gm1", "decay",
45 dt, tau, gm1, ROOT.RooDecay.DoubleSided)
46
47# Plot p.d.f.
48decay_gm1.plotOn(frame)
49
50# B-physics pdf with double Gaussian resolution
51# ------------------------------------------------------------------------------------------
52
53# Build another gaussian resolution model
54bias2 = ROOT.RooRealVar("bias2", "bias2", 0)
55sigma2 = ROOT.RooRealVar("sigma2", "sigma2", 5)
56gm2 = ROOT.RooGaussModel("gm2", "gauss model 2", dt, bias2, sigma2)
57
58# Build a composite resolution model f*gm1+(1-f)*gm2
59gm1frac = ROOT.RooRealVar("gm1frac", "fraction of gm1", 0.5)
60gmsum = ROOT.RooAddModel(
61 "gmsum",
62 "sum of gm1 and gm2",
63 ROOT.RooArgList(
64 gm1,
65 gm2),
66 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")