Logo ROOT   6.14/05
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 ///
13 /// \date March 2018
14 /// \author Danilo Piparo, Andre Vieira Silva
15 
16 auto filename = gROOT->GetTutorialDir() + "/dataframe/df017_vecOpsHEP.root";
17 auto treename = "myDataset";
18 using doubles = ROOT::VecOps::RVec<double>;
19 using RDF = ROOT::RDataFrame;
20 
21 void WithTTreeReader()
22 {
23  TFile f(filename);
24  TTreeReader tr(treename, &f);
25  TTreeReaderArray<double> px(tr, "px");
26  TTreeReaderArray<double> py(tr, "py");
27  TTreeReaderArray<double> E(tr, "E");
28 
29  TH1F h("pt", "pt", 16, 0, 4);
30 
31  while (tr.Next()) {
32  for (auto i=0U;i < px.GetSize(); ++i) {
33  if (E[i] > 100) h.Fill(sqrt(px[i]*px[i] + py[i]*py[i]));
34  }
35  }
36  h.DrawCopy();
37 }
38 
39 void WithRDataFrame()
40 {
41  RDF f(treename, filename.Data());
42  auto CalcPt = [](doubles &px, doubles &py, doubles &E) {
43  doubles v;
44  for (auto i=0U;i < px.size(); ++i) {
45  if (E[i] > 100) {
46  v.emplace_back(sqrt(px[i]*px[i] + py[i]*py[i]));
47  }
48  }
49  return v;
50  };
51  f.Define("pt", CalcPt, {"px", "py", "E"})
52  .Histo1D<doubles>({"pt", "pt", 16, 0, 4}, "pt")->DrawCopy();
53 }
54 
55 void WithRDataFrameVecOps()
56 {
57  RDF f(treename, filename.Data());
58  auto CalcPt = [](doubles &px, doubles &py, doubles &E) {
59  auto pt = sqrt(px*px + py*py);
60  return pt[E>100];
61  };
62  f.Define("good_pt", CalcPt, {"px", "py", "E"})
63  .Histo1D<doubles>({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
64 }
65 
66 void WithRDataFrameVecOpsJit()
67 {
68  RDF f(treename, filename.Data());
69  f.Define("good_pt", "sqrt(px*px + py*py)[E>100]")
70  .Histo1D({"pt", "pt", 16, 0, 4}, "good_pt")->DrawCopy();
71 }
72 
73 void df017_vecOpsHEP()
74 {
75  // We plot four times the same quantity, the key is to look into the implementation
76  // of the functions above
77  auto c = new TCanvas();
78  c->Divide(2,2);
79  c->cd(1);
80  WithTTreeReader();
81  c->cd(2);
82  WithRDataFrame();
83  c->cd(3);
84  WithRDataFrameVecOps();
85  c->cd(4);
86  WithRDataFrameVecOpsJit();
87 }
TTreeReader is a simple, robust and fast interface to read values from a TTree, TChain or TNtuple...
Definition: TTreeReader.h:43
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:47
#define gROOT
Definition: TROOT.h:410
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:567
#define f(i)
Definition: RSha256.hxx:104
double sqrt(double)
A "std::vector"-like collection of values implementing handy operation to analyse them...
Definition: RVec.hxx:146
TPaveText * pt
SVector< double, 2 > v
Definition: Dict.h:5
constexpr Double_t E()
Base of natural log: .
Definition: TMath.h:97
#define h(i)
Definition: RSha256.hxx:106
ROOT&#39;s RDataFrame offers a high level interface for analyses of data stored in TTrees, CSV&#39;s and other data formats.
Definition: RDataFrame.hxx:42
The Canvas class.
Definition: TCanvas.h:31
Extracts array data from a TTree.
#define c(i)
Definition: RSha256.hxx:101