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