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.0)
22
23mixState = ROOT.RooCategory("mixState", "B0/B0bar mixing state", {"mixed": -1, "unmixed": 1})
24tagFlav = ROOT.RooCategory("tagFlav", "Flavour of the tagged B0", {"B0": 1, "B0bar": -1})
25
26# Build a gaussian resolution model
27dterr = ROOT.RooRealVar("dterr", "dterr", 0.1, 1.0)
28bias1 = ROOT.RooRealVar("bias1", "bias1", 0)
29sigma1 = ROOT.RooRealVar("sigma1", "sigma1", 0.1)
30gm1 = ROOT.RooGaussModel("gm1", "gauss model 1", dt, bias1, sigma1)
31
32# Construct Bdecay (x) gauss
33bmix = ROOT.RooBMixDecay("bmix", "decay", dt, mixState, tagFlav, tau, dm, w, dw, gm1, type="DoubleSided")
34
35# Sample data from model
36# --------------------------------------------
37
38# Sample 2000 events in (dt,mixState,tagFlav) from bmix
39data = bmix.generate({dt, mixState, tagFlav}, 2000)
40
41# Show dt distribution with custom binning
42# -------------------------------------------------------------------------------
43
44# Make plot of dt distribution of data in range (-15,15) with fine binning
45# for dt>0 and coarse binning for dt<0
46
47# Create binning object with range (-15,15)
48tbins = ROOT.RooBinning(-15, 15)
49
50# Add 60 bins with uniform spacing in range (-15,0)
51tbins.addUniform(60, -15, 0)
52
53# Add 15 bins with uniform spacing in range (0,15)
54tbins.addUniform(15, 0, 15)
55
56# Make plot with specified binning
57dtframe = dt.frame(Range=(-15, 15), Title="dt distribution with custom binning")
58data.plotOn(dtframe, Binning=tbins)
59bmix.plotOn(dtframe)
60
61# NB: Note that bin density for each bin is adjusted to that of default frame binning as shown
62# in Y axis label (100 bins -. Events/0.4*Xaxis-dim) so that all bins
63# represent a consistent density distribution
64
65# Show mixstate asymmetry with custom binning
66# ------------------------------------------------------------------------------------
67
68# Make plot of dt distribution of data asymmetry in 'mixState' with
69# variable binning
70
71# Create binning object with range (-10,10)
72abins = ROOT.RooBinning(-10, 10)
73
74# Add boundaries at 0, (-1,1), (-2,2), (-3,3), (-4,4) and (-6,6)
75abins.addBoundary(0)
76abins.addBoundaryPair(1)
77abins.addBoundaryPair(2)
78abins.addBoundaryPair(3)
79abins.addBoundaryPair(4)
80abins.addBoundaryPair(6)
81
82# Create plot frame in dt
83aframe = dt.frame(Range=(-10, 10), Title="mixState asymmetry distribution with custom binning")
84
85# Plot mixState asymmetry of data with specified customg binning
86data.plotOn(aframe, Asymmetry=mixState, Binning=abins)
87
88# Plot corresponding property of pdf
89bmix.plotOn(aframe, Asymmetry=mixState)
90
91# Adjust vertical range of plot to sensible values for an asymmetry
92aframe.SetMinimum(-1.1)
93aframe.SetMaximum(1.1)
94
95# NB: For asymmetry distributions no density corrects are needed (and are
96# thus not applied)
97
98# Draw plots on canvas
99c = ROOT.TCanvas("rf108_plotbinning", "rf108_plotbinning", 800, 400)
100c.Divide(2)
101c.cd(1)
102ROOT.gPad.SetLeftMargin(0.15)
103dtframe.GetYaxis().SetTitleOffset(1.6)
104dtframe.Draw()
105c.cd(2)
106ROOT.gPad.SetLeftMargin(0.15)
107aframe.GetYaxis().SetTitleOffset(1.6)
108aframe.Draw()
109
110c.SaveAs("rf108_plotbinning.png")