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