Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df004_cutFlowReport.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Display cut/Filter efficiencies with RDataFrame.

This tutorial shows how to get information about the efficiency of the filters applied

using FourVector = ROOT::Math::XYZTVector;
using FourVectors = std::vector<FourVector>;
using CylFourVector = ROOT::Math::RhoEtaPhiVector;
// A simple helper function to fill a test tree: this makes the example
// stand-alone.
void fill_tree(const char *treeName, const char *fileName)
{
int i(0);
d.Define("b1", [&i]() { return (double)i; })
.Define("b2",
[&i]() {
auto j = i * i;
++i;
return j;
})
.Snapshot(treeName, fileName);
}
{
// We prepare an input tree to run on
auto fileName = "df004_cutFlowReport.root";
auto treeName = "myTree";
fill_tree(treeName, fileName);
// We read the tree from the file and create a RDataFrame
ROOT::RDataFrame d(treeName, fileName, {"b1", "b2"});
// ## Define cuts and create the report
// Here we define two simple cuts
auto cut1 = [](double b1) { return b1 > 25.; };
auto cut2 = [](int b2) { return 0 == b2 % 2; };
// An optional string parameter name can be passed to the Filter method to create a named filter.
// Named filters work as usual, but also keep track of how many entries they accept and reject.
auto filtered1 = d.Filter(cut1, {"b1"}, "Cut1");
auto filtered2 = d.Filter(cut2, {"b2"}, "Cut2");
auto augmented1 = filtered2.Define("b3", [](double b1, int b2) { return b1 / b2; });
auto cut3 = [](double x) { return x < .5; };
auto filtered3 = augmented1.Filter(cut3, {"b3"}, "Cut3");
// Statistics are retrieved through a call to the Report method:
// when Report is called on the main RDataFrame object, it retrieves stats
// for all named filters declared up to that point.
// When called on a stored chain state (i.e. a chain/graph node), it
// retrieves stats for all named filters in the section of the chain between
// the main RDataFrame and that node (included).
// Stats are printed in the same order as named filters that have been added to
// the graph, and refer to the latest event-loop that has been running using the
// relevant RDataFrame.
std::cout << "Cut3 stats:" << std::endl;
filtered3.Report()->Print();
// It is not only possible to print the information about cuts, but also to
// retrieve it to then use it programmatically.
std::cout << "All stats:" << std::endl;
auto allCutsReport = d.Report();
allCutsReport->Print();
// We can now loop on the cuts
std::cout << "Name\tAll\tPass\tEfficiency" << std::endl;
for (auto &&cutInfo : allCutsReport) {
std::cout << cutInfo.GetName() << "\t" << cutInfo.GetAll() << "\t" << cutInfo.GetPass() << "\t"
<< cutInfo.GetEff() << " %" << std::endl;
}
// Or get information about them individually
auto cutName = "Cut1";
auto cut = allCutsReport->At("Cut1");
std::cout << cutName << " efficiency is " << cut.GetEff() << " %" << std::endl;
}
#define d(i)
Definition RSha256.hxx:102
RInterface< RDFDetail::RFilter< F, Proxied >, DS_t > Filter(F f, const ColumnNames_t &columns={}, std::string_view name="")
Append a filter to the call graph.
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
Double_t x[n]
Definition legend1.C:17
DisplacementVector3D< CylindricalEta3D< double >, DefaultCoordinateSystemTag > RhoEtaPhiVector
3D Vector based on the eta based cylindrical coordinates rho, eta, phi in double precision.
Definition Vector3Dfwd.h:62
LorentzVector< PxPyPzE4D< double > > XYZTVector
LorentzVector based on x,y,x,t (or px,py,pz,E) coordinates in double precision with metric (-,...
Definition Vector4Dfwd.h:46
Cut3 stats:
Cut2 : pass=25 all=50 -- eff=50.00 % cumulative eff=50.00 %
Cut3 : pass=23 all=25 -- eff=92.00 % cumulative eff=46.00 %
All stats:
Cut1 : pass=24 all=50 -- eff=48.00 % cumulative eff=48.00 %
Cut2 : pass=25 all=50 -- eff=50.00 % cumulative eff=50.00 %
Cut3 : pass=23 all=25 -- eff=92.00 % cumulative eff=46.00 %
Name All Pass Efficiency
Cut1 50 24 48 %
Cut2 50 25 50 %
Cut3 25 23 92 %
Cut1 efficiency is 48 %
Date
December 2016
Author
Danilo Piparo (CERN)

Definition in file df004_cutFlowReport.C.