Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
17 auto filename = gROOT->GetTutorialDir() + "/dataframe/df017_vecOpsHEP.root";
18 auto treename = "myDataset";
19 using doubles = ROOT::VecOps::RVec<double>;
20 using RDF = ROOT::RDataFrame;
21 
22 void WithTTreeReader()
23 {
24  TFile f(filename);
25  TTreeReader tr(treename, &f);
26  TTreeReaderArray<double> px(tr, "px");
27  TTreeReaderArray<double> py(tr, "py");
28  TTreeReaderArray<double> E(tr, "E");
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 
40 void 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 
56 void 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 
67 void 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 
74 void 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 }