Logo ROOT   6.16/01
Reference Guide
rf405_realtocatfuncs.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4##
5## 'DATA AND CATEGORIES' RooFit tutorial macro #405
6##
7## Demonstration of real-.discrete mapping functions
8##
9## \macro_code
10##
11## \date February 2018
12## \author Clemens Lange
13## \author Wouter Verkerke (C version)
14
15
16import ROOT
17
18
19# Define pdf in x, sample dataset in x
20# ------------------------------------------------------------------------
21
22# Define a dummy PDF in x
23x = ROOT.RooRealVar("x", "x", 0, 10)
24a = ROOT.RooArgusBG("a", "argus(x)", x, ROOT.RooFit.RooConst(
25 10), ROOT.RooFit.RooConst(-1))
26
27# Generate a dummy dataset
28data = a.generate(ROOT.RooArgSet(x), 10000)
29
30# Create a threshold real -> cat function
31# --------------------------------------------------------------------------
32
33# A RooThresholdCategory is a category function that maps regions in a real-valued
34# input observable observables to state names. At construction time a 'default'
35# state name must be specified to which all values of x are mapped that are not
36# otherwise assigned
37xRegion = ROOT.RooThresholdCategory(
38 "xRegion", "region of x", x, "Background")
39
40# Specify thresholds and state assignments one-by-one.
41# Each statement specifies that all values _below_ the given value
42# (and above any lower specified threshold) are mapped to the
43# category state with the given name
44#
45# Background | SideBand | Signal | SideBand | Background
46# 4.23 5.23 8.23 9.23
47xRegion.addThreshold(4.23, "Background")
48xRegion.addThreshold(5.23, "SideBand")
49xRegion.addThreshold(8.23, "Signal")
50xRegion.addThreshold(9.23, "SideBand")
51
52# Use threshold function to plot data regions
53# ----------------------------------------------
54
55# Add values of threshold function to dataset so that it can be used as
56# observable
57data.addColumn(xRegion)
58
59# Make plot of data in x
60xframe = x.frame(ROOT.RooFit.Title(
61 "Demo of threshold and binning mapping functions"))
62data.plotOn(xframe)
63
64# Use calculated category to select sideband data
65data.plotOn(xframe, ROOT.RooFit.Cut("xRegion==xRegion::SideBand"),
66 ROOT.RooFit.MarkerColor(ROOT.kRed), ROOT.RooFit.LineColor(ROOT.kRed))
67
68# Create a binning real -> cat function
69# ----------------------------------------------------------------------
70
71# A RooBinningCategory is a category function that maps bins of a (named) binning definition
72# in a real-valued input observable observables to state names. The state names are automatically
73# constructed from the variable name, binning name and the bin number. If no binning name
74# is specified the default binning is mapped
75
76x.setBins(10, "coarse")
77xBins = ROOT.RooBinningCategory("xBins", "coarse bins in x", x, "coarse")
78
79# Use binning function for tabulation and plotting
80# -----------------------------------------------------------------------------------------------
81
82# Print table of xBins state multiplicity. Note that xBins does not need to be an observable in data
83# it can be a function of observables in data as well
84xbtable = data.table(xBins)
85xbtable.Print("v")
86
87# Add values of xBins function to dataset so that it can be used as
88# observable
89xb = data.addColumn(xBins)
90
91# Define range "alt" as including bins 1,3,5,7,9
92xb.setRange(
93 "alt", "x_coarse_bin1,x_coarse_bin3,x_coarse_bin5,x_coarse_bin7,x_coarse_bin9")
94
95# Construct subset of data matching range "alt" but only for the first
96# 5000 events and plot it on the frame
97dataSel = data.reduce(ROOT.RooFit.CutRange(
98 "alt"), ROOT.RooFit.EventRange(0, 5000))
99dataSel.plotOn(xframe, ROOT.RooFit.MarkerColor(ROOT.kGreen),
100 ROOT.RooFit.LineColor(ROOT.kGreen))
101
102c = ROOT.TCanvas("rf405_realtocatfuncs", "rf405_realtocatfuncs", 600, 600)
103xframe.SetMinimum(0.01)
104ROOT.gPad.SetLeftMargin(0.15)
105xframe.GetYaxis().SetTitleOffset(1.4)
106xframe.Draw()
107
108c.SaveAs("rf405_realtocatfuncs.png")