Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
df018_customActions.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// \notebook
4 /// This tutorial shows how to implement a custom action.
5 /// As an example, we build a helper for filling THns.
6 ///
7 /// \macro_code
8 /// \macro_output
9 ///
10 /// \date April 2018
11 /// \author Enrico Guiraud, Danilo Piparo
12 
13 // This is a custom action which respects a well defined interface. It supports parallelism,
14 // in the sense that it behaves correctly if implicit multi threading is enabled.
15 // We template it on:
16 // - The type of the internal THnT(s)
17 // - The dimension of the internal THnT(s)
18 // Note the plural: in presence of a MT execution, internally more than a single THnT is created.
19 template <typename T, unsigned int NDIM>
20 class THnHelper : public ROOT::Detail::RDF::RActionImpl<THnHelper<T, NDIM>> {
21 public:
22  /// This is a handy, expressive shortcut.
23  using THn_t = THnT<T>;
24  /// This type is a requirement for every helper.
25  using Result_t = THn_t;
26 
27 private:
28  std::vector<std::shared_ptr<THn_t>> fHistos; // one per data processing slot
29 
30 public:
31  /// This constructor takes all the parameters necessary to build the THnTs. In addition, it requires the names of
32  /// the columns which will be used.
33  THnHelper(std::string_view name, std::string_view title, std::array<int, NDIM> nbins, std::array<double, NDIM> xmins,
34  std::array<double, NDIM> xmax)
35  {
36  const auto nSlots = ROOT::IsImplicitMTEnabled() ? ROOT::GetImplicitMTPoolSize() : 1;
37  for (auto i : ROOT::TSeqU(nSlots)) {
38  fHistos.emplace_back(std::make_shared<THn_t>(std::string(name).c_str(), std::string(title).c_str(),
39  NDIM, nbins.data(), xmins.data(), xmax.data()));
40  (void)i;
41  }
42  }
43  THnHelper(THnHelper &&) = default;
44  THnHelper(const THnHelper &) = delete;
45  std::shared_ptr<THn_t> GetResultPtr() const { return fHistos[0]; }
46  void Initialize() {}
47  void InitTask(TTreeReader *, unsigned int) {}
48  /// This is a method executed at every entry
49  template <typename... ColumnTypes>
50  void Exec(unsigned int slot, ColumnTypes... values)
51  {
52  // Since THnT<T>::Fill expects a double*, we build it passing through a std::array.
53  std::array<double, sizeof...(ColumnTypes)> valuesArr{static_cast<double>(values)...};
54  fHistos[slot]->Fill(valuesArr.data());
55  }
56  /// This method is called at the end of the event loop. It is used to merge all the internal THnTs which
57  /// were used in each of the data processing slots.
58  void Finalize()
59  {
60  auto &res = fHistos[0];
61  for (auto slot : ROOT::TSeqU(1, fHistos.size())) {
62  res->Add(fHistos[slot].get());
63  }
64  }
65 
66  std::string GetActionName(){
67  return "THnHelper";
68  }
69 };
70 
71 void df018_customActions()
72 {
73  // We enable implicit parallelism
74  ROOT::EnableImplicitMT();
75 
76  // We create an empty RDataFrame which contains 4 columns filled with random numbers.
77  // The type of the numbers held by the columns are: double, double, float, int.
78  ROOT::RDataFrame d(128);
79  auto genD = []() { return gRandom->Uniform(-5, 5); };
80  auto genF = [&genD]() { return (float)genD(); };
81  auto genI = [&genD]() { return (int)genD(); };
82  auto dd = d.Define("x0", genD).Define("x1", genD).Define("x2", genF).Define("x3", genI);
83 
84  // Our Helper type: templated on the internal THnT type, the size, the types of the columns
85  // we'll use to fill.
86  using Helper_t = THnHelper<float, 4>;
87 
88  Helper_t helper{"myThN", // Name
89  "A THn with 4 dimensions", // Title
90  {4, 4, 8, 2}, // NBins
91  {-10., -10, -4., -6.}, // Axes min values
92  {10., 10, 5., 7.}}; // Axes max values
93 
94  // We book the action: it will be treated during the event loop.
95  auto myTHnT = dd.Book<double, double, float, int>(std::move(helper), {"x0", "x1", "x2", "x3"});
96 
97  myTHnT->Print();
98 }