Logo ROOT   6.16/01
Reference Guide
rf310_sliceplot.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'MULTIDIMENSIONAL MODELS' RooFit tutorial macro #310
6##
7## Projecting p.d.f and data slices in discrete observables
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15import ROOT
16
17# Create B decay pdf with mixing
18# ----------------------------------------------------------
19
20# Decay time observables
21dt = ROOT.RooRealVar("dt", "dt", -20, 20)
22
23# Discrete observables mixState (B0tag==B0reco?) and tagFlav
24# (B0tag==B0(bar)?)
25mixState = ROOT.RooCategory("mixState", "B0/B0bar mixing state")
26tagFlav = ROOT.RooCategory("tagFlav", "Flavour of the tagged B0")
27
28# Define state labels of discrete observables
29mixState.defineType("mixed", -1)
30mixState.defineType("unmixed", 1)
31tagFlav.defineType("B0", 1)
32tagFlav.defineType("B0bar", -1)
33
34# Model parameters
35dm = ROOT.RooRealVar("dm", "delta m(B)", 0.472, 0., 1.0)
36tau = ROOT.RooRealVar("tau", "B0 decay time", 1.547, 1.0, 2.0)
37w = ROOT.RooRealVar("w", "Flavor Mistag rate", 0.03, 0.0, 1.0)
38dw = ROOT.RooRealVar(
39 "dw", "Flavor Mistag rate difference between B0 and B0bar", 0.01)
40
41# Build a gaussian resolution model
42bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
43sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 0.01)
44gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
45
46# Construct a decay pdf, with single gaussian resolution model
47bmix_gm1 = ROOT.RooBMixDecay("bmix", "decay", dt, mixState,
48 tagFlav, tau, dm, w, dw, gm1, ROOT.RooBMixDecay.DoubleSided)
49
50# Generate BMixing data with above set of event errors
51data = bmix_gm1.generate(ROOT.RooArgSet(dt, tagFlav, mixState), 20000)
52
53# Plot full decay distribution
54# ----------------------------------------------------------
55
56# Create frame, data and pdf projection (integrated over tagFlav and
57# mixState)
58frame = dt.frame(ROOT.RooFit.Title("Inclusive decay distribution"))
59data.plotOn(frame)
60bmix_gm1.plotOn(frame)
61
62# Plot decay distribution for mixed and unmixed slice of mixState
63# -------------------------------------------------------------------------------------------
64
65# Create frame, data (mixed only)
66frame2 = dt.frame(ROOT.RooFit.Title("Decay distribution of mixed events"))
67data.plotOn(frame2, ROOT.RooFit.Cut("mixState==mixState::mixed"))
68
69# Position slice in mixState at "mixed" and plot slice of pdf in mixstate
70# over data (integrated over tagFlav)
71bmix_gm1.plotOn(frame2, ROOT.RooFit.Slice(mixState, "mixed"))
72
73# Create frame, data (unmixed only)
74frame3 = dt.frame(ROOT.RooFit.Title(
75 "Decay distribution of unmixed events"))
76data.plotOn(frame3, ROOT.RooFit.Cut("mixState==mixState::unmixed"))
77
78# Position slice in mixState at "unmixed" and plot slice of pdf in
79# mixstate over data (integrated over tagFlav)
80bmix_gm1.plotOn(frame3, ROOT.RooFit.Slice(mixState, "unmixed"))
81
82c = ROOT.TCanvas("rf310_sliceplot", "rf310_sliceplot", 1200, 400)
83c.Divide(3)
84c.cd(1)
85ROOT.gPad.SetLeftMargin(0.15)
86frame.GetYaxis().SetTitleOffset(1.4)
87ROOT.gPad.SetLogy()
88frame.Draw()
89c.cd(2)
90ROOT.gPad.SetLeftMargin(0.15)
91frame2.GetYaxis().SetTitleOffset(1.4)
92ROOT.gPad.SetLogy()
93frame2.Draw()
94c.cd(3)
95ROOT.gPad.SetLeftMargin(0.15)
96frame3.GetYaxis().SetTitleOffset(1.4)
97ROOT.gPad.SetLogy()
98frame3.Draw()
99
100c.SaveAs("rf310_sliceplot.png")