Logo ROOT  
Reference Guide
df017_vecOpsHEP.C File Reference

Detailed Description

View in nbviewer Open in SWAN This tutorial shows how VecOps can be used to slim down the programming model typically adopted in HEP for analysis. In this case we have a dataset containing the kinematic properties of particles stored in individual arrays. We want to plot the transverse momentum of these particles if the energy is greater than 100.

auto filename = gROOT->GetTutorialDir() + "/dataframe/df017_vecOpsHEP.root";
auto treename = "myDataset";
using RDF = ROOT::RDataFrame;
void WithTTreeReader()
{
TFile f(filename);
TTreeReader tr(treename, &f);
TH1F h("pt", "pt", 16, 0, 4);
while (tr.Next()) {
for (auto i=0U;i < px.GetSize(); ++i) {
if (E[i] > 100) h.Fill(sqrt(px[i]*px[i] + py[i]*py[i]));
}
}
h.DrawCopy();
}
void WithRDataFrame()
{
RDF f(treename, filename.Data());
auto CalcPt = [](doubles &px, doubles &py, doubles &E) {
doubles v;
for (auto i=0U;i < px.size(); ++i) {
if (E[i] > 100) {
v.emplace_back(sqrt(px[i]*px[i] + py[i]*py[i]));
}
}
return v;
};
f.Define("pt", CalcPt, {"px", "py", "E"})
.Histo1D<doubles>({"pt", "pt", 16, 0, 4}, "pt")->DrawCopy();
}
void WithRDataFrameVecOps()
{
RDF f(treename, filename.Data());
auto CalcPt = [](doubles &px, doubles &py, doubles &E) {
auto pt = sqrt(px*px + py*py);
return pt[E>100];
};
f.Define("good_pt", CalcPt, {"px", "py", "E"})
.Histo1D<doubles>({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
}
void WithRDataFrameVecOpsJit()
{
RDF f(treename, filename.Data());
f.Define("good_pt", "sqrt(px*px + py*py)[E>100]")
.Histo1D({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
}
{
// We plot four times the same quantity, the key is to look into the implementation
// of the functions above
auto c = new TCanvas();
c->Divide(2,2);
c->cd(1);
WithTTreeReader();
c->cd(2);
WithRDataFrame();
c->cd(3);
WithRDataFrameVecOps();
c->cd(4);
WithRDataFrameVecOpsJit();
}
Date
March 2018
Author
Danilo Piparo, Andre Vieira Silva

Definition in file df017_vecOpsHEP.C.

c
#define c(i)
Definition: RSha256.hxx:119
df017_vecOpsHEP
Definition: df017_vecOpsHEP.py:1
f
#define f(i)
Definition: RSha256.hxx:122
TTreeReaderArray
An interface for reading collections stored in ROOT columnar datasets.
Definition: TTreeReaderArray.h:75
ROOT::RDataFrame
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:42
v
@ v
Definition: rootcling_impl.cxx:3635
h
#define h(i)
Definition: RSha256.hxx:124
TTreeReader
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition: TTreeReader.h:44
sqrt
double sqrt(double)
TFile
Definition: TFile.h:54
TCanvas
Definition: TCanvas.h:23
TH1F
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:572
pt
TPaveText * pt
Definition: entrylist_figure1.C:7
TMath::E
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:102
ROOT::VecOps::RVec< double >
gROOT
#define gROOT
Definition: TROOT.h:406