Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
df016_vecOps.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// \notebook -draw
4 /// This tutorial shows the potential of the VecOps approach for treating collections
5 /// stored in datasets, a situation very common in HEP data analysis.
6 ///
7 /// \macro_code
8 /// \macro_image
9 ///
10 /// \date February 2018
11 /// \author Danilo Piparo
12 
13 using ROOT::RDataFrame;
14 using namespace ROOT::VecOps;
15 
16 int df016_vecOps()
17 {
18  // We re-create a set of points in a square.
19  // This is a technical detail, just to create a dataset to play with!
20  auto unifGen = [](double) { return gRandom->Uniform(-1.0, 1.0); };
21  auto vGen = [&](int len) {
22  RVec<double> v(len);
23  std::transform(v.begin(), v.end(), v.begin(), unifGen);
24  return v;
25  };
26  RDataFrame d(1024);
27  auto d0 = d.Define("len", []() { return (int)gRandom->Uniform(0, 16); })
28  .Define("x", vGen, {"len"})
29  .Define("y", vGen, {"len"});
30 
31  // Now we have in hands d, a RDataFrame with two columns, x and y, which
32  // hold collections of coordinates. The size of these collections vary.
33  // Let's now define radii out of x and y. We'll do it treating the collections
34  // stored in the columns without looping on the individual elements.
35  auto d1 = d0.Define("r", "sqrt(x*x + y*y)");
36 
37  // Now we want to plot 2 quarters of a ring with radii .5 and 1
38  // Note how the cuts are performed on RVecs, comparing them with integers and
39  // among themselves
40  auto ring_h = d1.Define("rInFig", "r > .4 && r < .8 && x*y < 0")
41  .Define("yFig", "y[rInFig]")
42  .Define("xFig", "x[rInFig]")
43  .Histo2D({"fig", "Two quarters of a ring", 64, -1, 1, 64, -1, 1}, "xFig", "yFig");
44 
45  auto cring = new TCanvas();
46  ring_h->DrawCopy("Colz");
47 
48  return 0;
49 }