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