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