Logo ROOT  
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 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_code
14 ##
15 ## \date February 2018
16 ## \authors Clemens Lange, Wouter Verkerke (C++ version)
17 
18 import ROOT
19 
20 # B-physics pdf with truth resolution
21 # ---------------------------------------------------------------------
22 
23 # Variables of decay pdf
24 dt = ROOT.RooRealVar("dt", "dt", -10, 10)
25 tau = ROOT.RooRealVar("tau", "tau", 1.548)
26 
27 # Build a truth resolution model (delta function)
28 tm = ROOT.RooTruthModel("tm", "truth model", dt)
29 
30 # Construct decay(t) (x) delta(t)
31 decay_tm = ROOT.RooDecay("decay_tm", "decay", dt,
32  tau, tm, ROOT.RooDecay.DoubleSided)
33 
34 # Plot pdf (dashed)
35 frame = dt.frame(ROOT.RooFit.Title("Bdecay (x) resolution"))
36 decay_tm.plotOn(frame, ROOT.RooFit.LineStyle(ROOT.kDashed))
37 
38 # B-physics pdf with Gaussian resolution
39 # ----------------------------------------------------------------------------
40 
41 # Build a gaussian resolution model
42 bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
43 sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 1)
44 gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
45 
46 # Construct decay(t) (x) gauss1(t)
47 decay_gm1 = ROOT.RooDecay("decay_gm1", "decay",
48  dt, tau, gm1, ROOT.RooDecay.DoubleSided)
49 
50 # Plot pdf
51 decay_gm1.plotOn(frame)
52 
53 # B-physics pdf with double Gaussian resolution
54 # ------------------------------------------------------------------------------------------
55 
56 # Build another gaussian resolution model
57 bias2 = ROOT.RooRealVar("bias2", "bias2", 0)
58 sigma2 = ROOT.RooRealVar("sigma2", "sigma2", 5)
59 gm2 = ROOT.RooGaussModel("gm2", "gauss model 2", dt, bias2, sigma2)
60 
61 # Build a composite resolution model f*gm1+(1-f)*gm2
62 gm1frac = ROOT.RooRealVar("gm1frac", "fraction of gm1", 0.5)
63 gmsum = ROOT.RooAddModel(
64  "gmsum",
65  "sum of gm1 and gm2",
66  ROOT.RooArgList(
67  gm1,
68  gm2),
69  ROOT.RooArgList(gm1frac))
70 
71 # Construct decay(t) (x) (f*gm1 + (1-f)*gm2)
72 decay_gmsum = ROOT.RooDecay(
73  "decay_gmsum", "decay", dt, tau, gmsum, ROOT.RooDecay.DoubleSided)
74 
75 # Plot pdf (red)
76 decay_gmsum.plotOn(frame, ROOT.RooFit.LineColor(ROOT.kRed))
77 
78 # Draw all frames on canvas
79 c = ROOT.TCanvas("rf209_anaconv", "rf209_anaconv", 600, 600)
80 ROOT.gPad.SetLeftMargin(0.15)
81 frame.GetYaxis().SetTitleOffset(1.6)
82 frame.Draw()
83 
84 c.SaveAs("rf209_anaconv.png")