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