Logo ROOT   6.16/01
Reference Guide
rf108_plotbinning.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## 'BASIC FUNCTIONALITY' RooFit tutorial macro #108
5## Plotting unbinned data with alternate and variable binnings
6##
7## \macro_code
8##
9## \date February 2018
10## \author Clemens Lange
11## \author Wouter Verkerke (C version)
12
13import ROOT
14
15# Set up model
16# ---------------------
17
18# Build a B decay p.d.f with mixing
19dt = ROOT.RooRealVar("dt", "dt", -20, 20)
20dm = ROOT.RooRealVar("dm", "dm", 0.472)
21tau = ROOT.RooRealVar("tau", "tau", 1.547)
22w = ROOT.RooRealVar("w", "mistag rate", 0.1)
23dw = ROOT.RooRealVar("dw", "delta mistag rate", 0.)
24
25mixState = ROOT.RooCategory("mixState", "B0/B0bar mixing state")
26mixState.defineType("mixed", -1)
27mixState.defineType("unmixed", 1)
28tagFlav = ROOT.RooCategory("tagFlav", "Flavour of the tagged B0")
29tagFlav.defineType("B0", 1)
30tagFlav.defineType("B0bar", -1)
31
32# Build a gaussian resolution model
33dterr = ROOT.RooRealVar("dterr", "dterr", 0.1, 1.0)
34bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
35sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 0.1)
36gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
37
38# Construct Bdecay (x) gauss
39bmix = ROOT.RooBMixDecay("bmix", "decay", dt, mixState, tagFlav,
40 tau, dm, w, dw, gm1, ROOT.RooBMixDecay.DoubleSided)
41
42# Sample data from model
43# --------------------------------------------
44
45# Sample 2000 events in (dt,mixState,tagFlav) from bmix
46data = bmix.generate(ROOT.RooArgSet(dt, mixState, tagFlav), 2000)
47
48# Show dt distribution with custom binning
49# -------------------------------------------------------------------------------
50
51# Make plot of dt distribution of data in range (-15,15) with fine binning
52# for dt>0 and coarse binning for dt<0
53
54# Create binning object with range (-15,15)
55tbins = ROOT.RooBinning(-15, 15)
56
57# Add 60 bins with uniform spacing in range (-15,0)
58tbins.addUniform(60, -15, 0)
59
60# Add 15 bins with uniform spacing in range (0,15)
61tbins.addUniform(15, 0, 15)
62
63# Make plot with specified binning
64dtframe = dt.frame(ROOT.RooFit.Range(-15, 15),
65 ROOT.RooFit.Title("dt distribution with custom binning"))
66data.plotOn(dtframe, ROOT.RooFit.Binning(tbins))
67bmix.plotOn(dtframe)
68
69# NB: Note that bin density for each bin is adjusted to that of default frame binning as shown
70# in Y axis label (100 bins -. Events/0.4*Xaxis-dim) so that all bins
71# represent a consistent density distribution
72
73# Show mixstate asymmetry with custom binning
74# ------------------------------------------------------------------------------------
75
76# Make plot of dt distribution of data asymmetry in 'mixState' with
77# variable binning
78
79# Create binning object with range (-10,10)
80abins = ROOT.RooBinning(-10, 10)
81
82# Add boundaries at 0, (-1,1), (-2,2), (-3,3), (-4,4) and (-6,6)
83abins.addBoundary(0)
84abins.addBoundaryPair(1)
85abins.addBoundaryPair(2)
86abins.addBoundaryPair(3)
87abins.addBoundaryPair(4)
88abins.addBoundaryPair(6)
89
90# Create plot frame in dt
91aframe = dt.frame(ROOT.RooFit.Range(-10, 10), ROOT.RooFit.Title(
92 "mixState asymmetry distribution with custom binning"))
93
94# Plot mixState asymmetry of data with specified customg binning
95data.plotOn(aframe, ROOT.RooFit.Asymmetry(
96 mixState), ROOT.RooFit.Binning(abins))
97
98# Plot corresponding property of p.d.f
99bmix.plotOn(aframe, ROOT.RooFit.Asymmetry(mixState))
100
101# Adjust vertical range of plot to sensible values for an asymmetry
102aframe.SetMinimum(-1.1)
103aframe.SetMaximum(1.1)
104
105# NB: For asymmetry distributions no density corrects are needed (and are
106# thus not applied)
107
108# Draw plots on canvas
109c = ROOT.TCanvas("rf108_plotbinning", "rf108_plotbinning", 800, 400)
110c.Divide(2)
111c.cd(1)
112ROOT.gPad.SetLeftMargin(0.15)
113dtframe.GetYaxis().SetTitleOffset(1.6)
114dtframe.Draw()
115c.cd(2)
116ROOT.gPad.SetLeftMargin(0.15)
117aframe.GetYaxis().SetTitleOffset(1.6)
118aframe.Draw()
119
120c.SaveAs("rf108_plotbinning.png")