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