19 template <
typename T,
unsigned int NDIM>
20 class THnHelper :
public ROOT::Detail::RDF::RActionImpl<THnHelper<T, NDIM>> {
23 using THn_t = THnT<T>;
25 using Result_t = THn_t;
28 std::vector<std::shared_ptr<THn_t>> fHistos;
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)
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()));
43 THnHelper(THnHelper &&) =
default;
44 THnHelper(
const THnHelper &) =
delete;
45 std::shared_ptr<THn_t> GetResultPtr()
const {
return fHistos[0]; }
47 void InitTask(TTreeReader *,
unsigned int) {}
49 template <
typename... ColumnTypes>
50 void Exec(
unsigned int slot, ColumnTypes... values)
53 std::array<double,
sizeof...(ColumnTypes)> valuesArr{
static_cast<double>(values)...};
54 fHistos[slot]->Fill(valuesArr.data());
60 auto &res = fHistos[0];
61 for (
auto slot : ROOT::TSeqU(1, fHistos.size())) {
62 res->Add(fHistos[slot].
get());
66 std::string GetActionName(){
71 void df018_customActions()
74 ROOT::EnableImplicitMT();
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);
86 using Helper_t = THnHelper<float, 4>;
88 Helper_t helper{
"myThN",
89 "A THn with 4 dimensions",
91 {-10., -10, -4., -6.},
95 auto myTHnT = dd.Book<double, double, float,
int>(std::move(helper), {
"x0",
"x1",
"x2",
"x3"});