Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf209_anaconv.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Addition and convolution: decay function pdfs with optional B physics effects (mixing
5## and CP violation) that can be analytically convolved with e.g. Gaussian resolution functions
6##
7## ```
8## pdf1 = decay(t,tau) (x) delta(t)
9## pdf2 = decay(t,tau) (x) gauss(t,m,s)
10## pdf3 = decay(t,tau) (x) (f*gauss1(t,m1,s1) + (1-f)*gauss2(t,m1,s1))
11## ```
12##
13## \macro_image
14## \macro_code
15## \macro_output
16##
17## \date February 2018
18## \authors Clemens Lange, Wouter Verkerke (C++ version)
19
20import ROOT
21
22# B-physics pdf with truth resolution
23# ---------------------------------------------------------------------
24
25# Variables of decay pdf
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, tau, tm, type="DoubleSided")
34
35# Plot pdf (dashed)
36frame = dt.frame(Title="Bdecay (x) resolution")
37decay_tm.plotOn(frame, LineStyle="--")
38
39# B-physics pdf with Gaussian resolution
40# ----------------------------------------------------------------------------
41
42# Build a gaussian resolution model
43bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
44sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 1)
45gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
46
47# Construct decay(t) (x) gauss1(t)
48decay_gm1 = ROOT.RooDecay("decay_gm1", "decay", dt, tau, gm1, type="DoubleSided")
49
50# Plot pdf
51decay_gm1.plotOn(frame)
52
53# B-physics pdf with double Gaussian resolution
54# ------------------------------------------------------------------------------------------
55
56# Build another gaussian resolution model
57bias2 = ROOT.RooRealVar("bias2", "bias2", 0)
58sigma2 = ROOT.RooRealVar("sigma2", "sigma2", 5)
59gm2 = ROOT.RooGaussModel("gm2", "gauss model 2", dt, bias2, sigma2)
60
61# Build a composite resolution model f*gm1+(1-f)*gm2
62gm1frac = ROOT.RooRealVar("gm1frac", "fraction of gm1", 0.5)
63gmsum = ROOT.RooAddModel("gmsum", "sum of gm1 and gm2", [gm1, gm2], [gm1frac])
64
65# Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
66decay_gmsum = ROOT.RooDecay("decay_gmsum", "decay", dt, tau, gmsum, ROOT.RooDecay.DoubleSided)
67
68# Plot pdf (red)
69decay_gmsum.plotOn(frame, LineColor="r")
70
71# Draw all frames on canvas
72c = ROOT.TCanvas("rf209_anaconv", "rf209_anaconv", 600, 600)
73ROOT.gPad.SetLeftMargin(0.15)
74frame.GetYaxis().SetTitleOffset(1.6)
75frame.Draw()
76
77c.SaveAs("rf209_anaconv.png")