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