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)?)
21# Define state labels of discrete observables
22mixState = ROOT.RooCategory("mixState", "B0/B0bar mixing state", {"mixed": -1, "unmixed": 1})
23tagFlav = ROOT.RooCategory("tagFlav", "Flavour of the tagged B0", {"B0": 1, "B0bar": -1})
24
25# Model parameters
26dm = ROOT.RooRealVar("dm", "delta m(B)", 0.472, 0.0, 1.0)
27tau = ROOT.RooRealVar("tau", "B0 decay time", 1.547, 1.0, 2.0)
28w = ROOT.RooRealVar("w", "Flavor Mistag rate", 0.03, 0.0, 1.0)
29dw = ROOT.RooRealVar("dw", "Flavor Mistag rate difference between B0 and B0bar", 0.01)
30
31# Build a gaussian resolution model
32bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
33sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 0.01)
34gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
35
36# Construct a decay pdf, with single gaussian resolution model
37bmix_gm1 = ROOT.RooBMixDecay("bmix", "decay", dt, mixState, tagFlav, tau, dm, w, dw, gm1, type="DoubleSided")
38
39# Generate BMixing data with above set of event errors
40data = bmix_gm1.generate({dt, tagFlav, mixState}, 20000)
41
42# Plot full decay distribution
43# ----------------------------------------------------------
44
45# Create frame, data and pdf projection (integrated over tagFlav and
46# mixState)
47frame = dt.frame(Title="Inclusive decay distribution")
48data.plotOn(frame)
49bmix_gm1.plotOn(frame)
50
51# Plot decay distribution for mixed and unmixed slice of mixState
52# -------------------------------------------------------------------------------------------
53
54# Create frame, data (mixed only)
55frame2 = dt.frame(Title="Decay distribution of mixed events")
56data.plotOn(frame2, Cut="mixState==mixState::mixed")
57
58# Position slice in mixState at "mixed" and plot slice of pdf in mixstate
59# over data (integrated over tagFlav)
60bmix_gm1.plotOn(frame2, Slice=(mixState, "mixed"))
61
62# Create frame, data (unmixed only)
63frame3 = dt.frame(Title="Decay distribution of unmixed events")
64data.plotOn(frame3, Cut="mixState==mixState::unmixed")
65
66# Position slice in mixState at "unmixed" and plot slice of pdf in
67# mixstate over data (integrated over tagFlav)
68bmix_gm1.plotOn(frame3, Slice=(mixState, "unmixed"))
69
70c = ROOT.TCanvas("rf310_sliceplot", "rf310_sliceplot", 1200, 400)
71c.Divide(3)
72c.cd(1)
73ROOT.gPad.SetLeftMargin(0.15)
74frame.GetYaxis().SetTitleOffset(1.4)
75ROOT.gPad.SetLogy()
76frame.Draw()
77c.cd(2)
78ROOT.gPad.SetLeftMargin(0.15)
79frame2.GetYaxis().SetTitleOffset(1.4)
80ROOT.gPad.SetLogy()
81frame2.Draw()
82c.cd(3)
83ROOT.gPad.SetLeftMargin(0.15)
84frame3.GetYaxis().SetTitleOffset(1.4)
85ROOT.gPad.SetLogy()
86frame3.Draw()
87
88c.SaveAs("rf310_sliceplot.png")