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 MeV.
auto filename =
gROOT->GetTutorialDir() +
"/analysis/dataframe/df017_vecOpsHEP.root";
auto treename = "myDataset";
void WithTTreeReader()
{
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]));
}
}
}
void WithRDataFrame()
{
for (
auto i=0U;i < px.
size(); ++i) {
if (E[i] > 100) {
v.emplace_back(
sqrt(px[i]*px[i] + py[i]*py[i]));
}
}
};
f.Define(
"pt", CalcPt, {
"px",
"py",
"E"})
.Histo1D<RVecD>({"pt", "pt", 16, 0, 4}, "pt")->DrawCopy();
}
void WithRDataFrameVecOps()
{
auto pt =
sqrt(px*px + py*py);
};
f.Define(
"good_pt", CalcPt, {
"px",
"py",
"E"})
.Histo1D<RVecD>({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
}
void WithRDataFrameVecOpsJit()
{
f.Define(
"good_pt",
"sqrt(px*px + py*py)[E>100]")
.Histo1D({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
}
void df017_vecOpsHEP()
{
WithTTreeReader();
WithRDataFrame();
WithRDataFrameVecOps();
WithRDataFrameVecOpsJit();
}
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
A file, usually with extension .root, that stores data and code in the form of serialized objects in ...
1-D histogram with a float per channel (see TH1 documentation)
An interface for reading collections stored in ROOT columnar datasets.
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
ROOT::VecOps::RVec< double > RVecD
constexpr Double_t E()
Base of natural log: .