Logo ROOT   6.16/01
Reference Guide
df017_vecOpsHEP.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -draw
4/// This tutorial shows how VecOps can be used to slim down the programming
5/// model typically adopted in HEP for analysis.
6/// In this case we have a dataset containing the kinematic properties of
7/// particles stored in individual arrays.
8/// We want to plot the transverse momentum of these particles if the energy is
9/// greater than 100.
10///
11/// \macro_code
12/// \macro_image
13///
14/// \date March 2018
15/// \author Danilo Piparo, Andre Vieira Silva
16
17auto filename = gROOT->GetTutorialDir() + "/dataframe/df017_vecOpsHEP.root";
18auto treename = "myDataset";
19using doubles = ROOT::VecOps::RVec<double>;
20using RDF = ROOT::RDataFrame;
21
22void WithTTreeReader()
23{
24 TFile f(filename);
25 TTreeReader tr(treename, &f);
26 TTreeReaderArray<double> px(tr, "px");
27 TTreeReaderArray<double> py(tr, "py");
29
30 TH1F h("pt", "pt", 16, 0, 4);
31
32 while (tr.Next()) {
33 for (auto i=0U;i < px.GetSize(); ++i) {
34 if (E[i] > 100) h.Fill(sqrt(px[i]*px[i] + py[i]*py[i]));
35 }
36 }
37 h.DrawCopy();
38}
39
40void WithRDataFrame()
41{
42 RDF f(treename, filename.Data());
43 auto CalcPt = [](doubles &px, doubles &py, doubles &E) {
44 doubles v;
45 for (auto i=0U;i < px.size(); ++i) {
46 if (E[i] > 100) {
47 v.emplace_back(sqrt(px[i]*px[i] + py[i]*py[i]));
48 }
49 }
50 return v;
51 };
52 f.Define("pt", CalcPt, {"px", "py", "E"})
53 .Histo1D<doubles>({"pt", "pt", 16, 0, 4}, "pt")->DrawCopy();
54}
55
56void WithRDataFrameVecOps()
57{
58 RDF f(treename, filename.Data());
59 auto CalcPt = [](doubles &px, doubles &py, doubles &E) {
60 auto pt = sqrt(px*px + py*py);
61 return pt[E>100];
62 };
63 f.Define("good_pt", CalcPt, {"px", "py", "E"})
64 .Histo1D<doubles>({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
65}
66
67void WithRDataFrameVecOpsJit()
68{
69 RDF f(treename, filename.Data());
70 f.Define("good_pt", "sqrt(px*px + py*py)[E>100]")
71 .Histo1D({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
72}
73
74void df017_vecOpsHEP()
75{
76 // We plot four times the same quantity, the key is to look into the implementation
77 // of the functions above
78 auto c = new TCanvas();
79 c->Divide(2,2);
80 c->cd(1);
81 WithTTreeReader();
82 c->cd(2);
83 WithRDataFrame();
84 c->cd(3);
85 WithRDataFrameVecOps();
86 c->cd(4);
87 WithRDataFrameVecOpsJit();
88}
SVector< double, 2 > v
Definition: Dict.h:5
#define f(i)
Definition: RSha256.hxx:104
#define c(i)
Definition: RSha256.hxx:101
#define h(i)
Definition: RSha256.hxx:106
double sqrt(double)
#define gROOT
Definition: TROOT.h:410
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:41
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition: RVec.hxx:221
The Canvas class.
Definition: TCanvas.h:31
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
An interface for reading collections stored in ROOT columnar datasets.
A simple, robust and fast interface to read values from ROOT colmnar datasets such as TTree,...
Definition: TTreeReader.h:44
TPaveText * pt
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97