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