Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
df002_dataModel.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// \notebook -draw
4 /// This tutorial shows the possibility to use data models which are more
5 /// complex than flat ntuples with RDataFrame
6 ///
7 /// \macro_code
8 /// \macro_image
9 ///
10 /// \date December 2016
11 /// \author Danilo Piparo
12 using FourVector = ROOT::Math::XYZTVector;
13 using FourVectorVec = std::vector<FourVector>;
14 using FourVectorRVec = ROOT::VecOps::RVec<FourVector>;
15 using CylFourVector = ROOT::Math::RhoEtaPhiVector;
16 
17 // A simple helper function to fill a test tree: this makes the example
18 // stand-alone.
19 void fill_tree(const char *filename, const char *treeName)
20 {
21  const double M = 0.13957; // set pi+ mass
22  TRandom3 R(1);
23 
24  auto genTracks = [&](){
25  FourVectorVec tracks;
26  const auto nPart = R.Poisson(15);
27  tracks.reserve(nPart);
28  for (int j = 0; j < nPart; ++j) {
29  const auto px = R.Gaus(0, 10);
30  const auto py = R.Gaus(0, 10);
31  const auto pt = sqrt(px * px + py * py);
32  const auto eta = R.Uniform(-3, 3);
33  const auto phi = R.Uniform(0.0, 2 * TMath::Pi());
34  CylFourVector vcyl(pt, eta, phi);
35  // set energy
36  auto E = sqrt(vcyl.R() * vcyl.R() + M * M);
37  // fill track vector
38  tracks.emplace_back(vcyl.X(), vcyl.Y(), vcyl.Z(), E);
39  }
40  return tracks;
41  };
42 
43  ROOT::RDataFrame d(64);
44  d.Define("tracks", genTracks).Snapshot<FourVectorVec>(treeName, filename, {"tracks"});
45 }
46 
47 int df002_dataModel()
48 {
49 
50  // We prepare an input tree to run on
51  auto fileName = "df002_dataModel.root";
52  auto treeName = "myTree";
53  fill_tree(fileName, treeName);
54 
55  // We read the tree from the file and create a RDataFrame, a class that
56  // allows us to interact with the data contained in the tree.
57  ROOT::RDataFrame d(treeName, fileName, {"tracks"});
58 
59  // ## Operating on branches which are collection of objects
60  // Here we deal with the simplest of the cuts: we decide to accept the event
61  // only if the number of tracks is greater than 5.
62  auto n_cut = [](const FourVectorRVec &tracks) { return tracks.size() > 8; };
63  auto nentries = d.Filter(n_cut, {"tracks"}).Count();
64 
65  std::cout << *nentries << " passed all filters" << std::endl;
66 
67  // Another possibility consists in creating a new column containing the
68  // quantity we are interested in.
69  // In this example, we will cut on the number of tracks and plot their
70  // transverse momentum.
71  auto getPt = [](const FourVectorRVec &tracks) {
72  return ROOT::VecOps::Map(tracks, [](const FourVector& v){return v.Pt();});
73  };
74 
75  // We do the same for the weights.
76  auto getPtWeights = [](const FourVectorRVec &tracks) {
77  return ROOT::VecOps::Map(tracks, [](const FourVector& v){ return 1. / v.Pt();});
78  };
79 
80  auto augmented_d = d.Define("tracks_n", [](const FourVectorRVec &tracks) { return (int)tracks.size(); })
81  .Filter([](int tracks_n) { return tracks_n > 2; }, {"tracks_n"})
82  .Define("tracks_pts", getPt)
83  .Define("tracks_pts_weights", getPtWeights);
84 
85  auto trN = augmented_d.Histo1D({"", "", 40, -.5, 39.5}, "tracks_n");
86  auto trPts = augmented_d.Histo1D("tracks_pts");
87  auto trWPts = augmented_d.Histo1D("tracks_pts", "tracks_pts_weights");
88 
89  auto c1 = new TCanvas();
90  trN->DrawCopy();
91 
92  auto c2 = new TCanvas();
93  trPts->DrawCopy();
94 
95  auto c3 = new TCanvas();
96  trWPts->DrawCopy();
97 
98  return 0;
99 }