Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf404_categories.py File Reference

Namespaces

namespace  rf404_categories
 

Detailed Description

View in nbviewer Open in SWAN
Data and categories: working with ROOT.RooCategory objects to describe discrete variables

from __future__ import print_function
import ROOT
# Construct a category with labels
# --------------------------------------------
# Define a category with labels only
tagCat = ROOT.RooCategory("tagCat", "Tagging category")
tagCat.defineType("Lepton")
tagCat.defineType("Kaon")
tagCat.defineType("NetTagger-1")
tagCat.defineType("NetTagger-2")
tagCat.Print()
# Construct a category with labels and indices
# ------------------------------------------------
# Define a category with explicitly numbered states
b0flav = ROOT.RooCategory("b0flav", "B0 flavour eigenstate", {"B0": -1, "B0bar": 1})
b0flav.Print()
# Generate dummy data for tabulation demo
# ------------------------------------------------
# Generate a dummy dataset
x = ROOT.RooRealVar("x", "x", 0, 10)
data = ROOT.RooPolynomial("p", "p", x).generate({x, b0flav, tagCat}, 10000)
# Print tables of category contents of datasets
# --------------------------------------------------
# Tables are equivalent of plots for categories
btable = data.table(b0flav)
btable.Print()
btable.Print("v")
# Create table for subset of events matching cut expression
ttable = data.table(tagCat, "x>8.23")
ttable.Print()
ttable.Print("v")
# Create table for all (tagCat x b0flav) state combinations
bttable = data.table({tagCat, b0flav})
bttable.Print("v")
# Retrieve number of events from table
# Number can be non-integer if source dataset has weighed events
nb0 = btable.get("B0")
print("Number of events with B0 flavor is ", nb0)
# Retrieve fraction of events with "Lepton" tag
fracLep = ttable.getFrac("Lepton")
print("Fraction of events tagged with Lepton tag is ", fracLep)
# Defining ranges for plotting, fitting on categories
# ------------------------------------------------------------------------------------------------------
# Define named range as comma separated list of labels
tagCat.setRange("good", "Lepton,Kaon")
# Or add state names one by one
tagCat.addToRange("soso", "NetTagger-1")
tagCat.addToRange("soso", "NetTagger-2")
# Use category range in dataset reduction specification
goodData = data.reduce(CutRange="good")
goodData.table(tagCat).Print("v")
RooCategory::tagCat = Lepton(idx = 0)
RooCategory::b0flav = B0(idx = -1)
Roo1DTable::b0flav = (B0=5058,B0bar=4942)
Table b0flav : pData
+-------+------+
| B0 | 5058 |
| B0bar | 4942 |
+-------+------+
[#1] INFO:InputArguments -- The formula cutVar claims to use the variables (x,b0flav,tagCat) but only (x) seem to be in use.
inputs: x>8.23
Roo1DTable::tagCat = (Lepton=454,Kaon=450,NetTagger-1=432,NetTagger-2=429)
Table tagCat : pData(x>8.23)
+-------------+-----+
| Lepton | 454 |
| Kaon | 450 |
| NetTagger-1 | 432 |
| NetTagger-2 | 429 |
+-------------+-----+
Table (tagCat x b0flav) : pData
+---------------------+------+
| {Lepton;B0} | 1302 |
| {Kaon;B0} | 1232 |
| {NetTagger-1;B0} | 1242 |
| {NetTagger-2;B0} | 1282 |
| {Lepton;B0bar} | 1192 |
| {Kaon;B0bar} | 1314 |
| {NetTagger-1;B0bar} | 1208 |
| {NetTagger-2;B0bar} | 1228 |
+---------------------+------+
Table tagCat : pData
+-------------+------+
| Lepton | 2494 |
| Kaon | 2546 |
| NetTagger-1 | 0 |
| NetTagger-2 | 0 |
+-------------+------+
Number of events with B0 flavor is 5058.0
Fraction of events tagged with Lepton tag is 0.2572237960339943
Date
February 2018
Authors
Clemens Lange, Wouter Verkerke (C++ version)

Definition in file rf404_categories.py.