Logo ROOT  
Reference Guide
df018_customActions.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_dataframe
3 /// \notebook
4 /// \brief Implement a custom action to fill THns.
5 ///
6 /// This tutorial shows how to implement a custom action.
7 /// As an example, we build a helper for filling THns.
8 ///
9 /// \macro_code
10 /// \macro_output
11 ///
12 /// \date April 2018
13 /// \author Enrico Guiraud, Danilo Piparo
14 
15 // This is a custom action which respects a well defined interface. It supports parallelism,
16 // in the sense that it behaves correctly if implicit multi threading is enabled.
17 // We template it on:
18 // - The type of the internal THnT(s)
19 // - The dimension of the internal THnT(s)
20 // Note the plural: in presence of a MT execution, internally more than a single THnT is created.
21 template <typename T, unsigned int NDIM>
22 class THnHelper : public ROOT::Detail::RDF::RActionImpl<THnHelper<T, NDIM>> {
23 public:
24  /// This is a handy, expressive shortcut.
25  using THn_t = THnT<T>;
26  /// This type is a requirement for every helper.
27  using Result_t = THn_t;
28 
29 private:
30  std::vector<std::shared_ptr<THn_t>> fHistos; // one per data processing slot
31 
32 public:
33  /// This constructor takes all the parameters necessary to build the THnTs. In addition, it requires the names of
34  /// the columns which will be used.
35  THnHelper(std::string_view name, std::string_view title, std::array<int, NDIM> nbins, std::array<double, NDIM> xmins,
36  std::array<double, NDIM> xmax)
37  {
38  const auto nSlots = ROOT::IsImplicitMTEnabled() ? ROOT::GetThreadPoolSize() : 1;
39  for (auto i : ROOT::TSeqU(nSlots)) {
40  fHistos.emplace_back(std::make_shared<THn_t>(std::string(name).c_str(), std::string(title).c_str(),
41  NDIM, nbins.data(), xmins.data(), xmax.data()));
42  (void)i;
43  }
44  }
45  THnHelper(THnHelper &&) = default;
46  THnHelper(const THnHelper &) = delete;
47  std::shared_ptr<THn_t> GetResultPtr() const { return fHistos[0]; }
48  void Initialize() {}
49  void InitTask(TTreeReader *, unsigned int) {}
50  /// This is a method executed at every entry
51  template <typename... ColumnTypes>
52  void Exec(unsigned int slot, ColumnTypes... values)
53  {
54  // Since THnT<T>::Fill expects a double*, we build it passing through a std::array.
55  std::array<double, sizeof...(ColumnTypes)> valuesArr{static_cast<double>(values)...};
56  fHistos[slot]->Fill(valuesArr.data());
57  }
58  /// This method is called at the end of the event loop. It is used to merge all the internal THnTs which
59  /// were used in each of the data processing slots.
60  void Finalize()
61  {
62  auto &res = fHistos[0];
63  for (auto slot : ROOT::TSeqU(1, fHistos.size())) {
64  res->Add(fHistos[slot].get());
65  }
66  }
67 
68  std::string GetActionName(){
69  return "THnHelper";
70  }
71 };
72 
73 void df018_customActions()
74 {
75  // We enable implicit parallelism
77 
78  // We create an empty RDataFrame which contains 4 columns filled with random numbers.
79  // The type of the numbers held by the columns are: double, double, float, int.
80  ROOT::RDataFrame d(128);
81  auto genD = []() { return gRandom->Uniform(-5, 5); };
82  auto genF = [&genD]() { return (float)genD(); };
83  auto genI = [&genD]() { return (int)genD(); };
84  auto dd = d.Define("x0", genD).Define("x1", genD).Define("x2", genF).Define("x3", genI);
85 
86  // Our Helper type: templated on the internal THnT type, the size, the types of the columns
87  // we'll use to fill.
88  using Helper_t = THnHelper<float, 4>;
89 
90  Helper_t helper{"myThN", // Name
91  "A THn with 4 dimensions", // Title
92  {4, 4, 8, 2}, // NBins
93  {-10., -10, -4., -6.}, // Axes min values
94  {10., 10, 5., 7.}}; // Axes max values
95 
96  // We book the action: it will be treated during the event loop.
97  auto myTHnT = dd.Book<double, double, float, int>(std::move(helper), {"x0", "x1", "x2", "x3"});
98 
99  myTHnT->Print();
100 }
xmax
float xmax
Definition: THbookFile.cxx:95
string_view
basic_string_view< char > string_view
Definition: libcpp_string_view.h:785
TRandom::Uniform
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:635
ROOT::RDataFrame
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTrees,...
Definition: RDataFrame.hxx:42
ROOT::EnableImplicitMT
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Definition: TROOT.cxx:525
ROOT::GetThreadPoolSize
UInt_t GetThreadPoolSize()
Returns the size of ROOT's thread pool.
Definition: TROOT.cxx:563
TTreeReader
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition: TTreeReader.h:44
gRandom
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
double
double
Definition: Converters.cxx:921
void
typedef void((*Func_t)())
TMVA::TMVAGlob::Initialize
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition: tmvaglob.cxx:176
ROOT::IsImplicitMTEnabled
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:556
name
char name[80]
Definition: TGX11.cxx:110
d
#define d(i)
Definition: RSha256.hxx:102
ROOT::TSeq
A pseudo container class which is a generator of indices.
Definition: TSeq.hxx:66
CPyCppyy::Exec
CPYCPPYY_EXTERN bool Exec(const std::string &cmd)
Definition: API.cxx:331
THnT
Templated implementation of the abstract base THn.
Definition: THn.h:222