Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df018_customActions.C File Reference

Detailed Description

View in nbviewer Open in SWAN Implement a custom action to fill THns.

This tutorial shows how to implement a custom action. As an example, we build a helper for filling THns.

// This is a custom action which respects a well defined interface. It supports parallelism,
// in the sense that it behaves correctly if implicit multi threading is enabled.
// We template it on:
// - The type of the internal THnT(s)
// - The dimension of the internal THnT(s)
// Note the plural: in presence of a MT execution, internally more than a single THnT is created.
template <typename T, unsigned int NDIM>
class THnHelper : public ROOT::Detail::RDF::RActionImpl<THnHelper<T, NDIM>> {
public:
/// This is a handy, expressive shortcut.
using THn_t = THnT<T>;
/// This type is a requirement for every helper.
using Result_t = THn_t;
private:
std::vector<std::shared_ptr<THn_t>> fHistos; // one per data processing slot
public:
/// This constructor takes all the parameters necessary to build the THnTs. In addition, it requires the names of
/// the columns which will be used.
THnHelper(std::string_view name, std::string_view title, std::array<int, NDIM> nbins, std::array<double, NDIM> xmins,
std::array<double, NDIM> xmax)
{
const auto nSlots = ROOT::IsImplicitMTEnabled() ? ROOT::GetThreadPoolSize() : 1;
for (auto i : ROOT::TSeqU(nSlots)) {
fHistos.emplace_back(std::make_shared<THn_t>(std::string(name).c_str(), std::string(title).c_str(),
NDIM, nbins.data(), xmins.data(), xmax.data()));
(void)i;
}
}
THnHelper(THnHelper &&) = default;
THnHelper(const THnHelper &) = delete;
std::shared_ptr<THn_t> GetResultPtr() const { return fHistos[0]; }
void Initialize() {}
void InitTask(TTreeReader *, unsigned int) {}
/// This is a method executed at every entry
template <typename... ColumnTypes>
void Exec(unsigned int slot, ColumnTypes... values)
{
// Since THnT<T>::Fill expects a double*, we build it passing through a std::array.
std::array<double, sizeof...(ColumnTypes)> valuesArr{static_cast<double>(values)...};
fHistos[slot]->Fill(valuesArr.data());
}
/// This method is called at the end of the event loop. It is used to merge all the internal THnTs which
/// were used in each of the data processing slots.
void Finalize()
{
auto &res = fHistos[0];
for (auto slot : ROOT::TSeqU(1, fHistos.size())) {
res->Add(fHistos[slot].get());
}
}
std::string GetActionName(){
return "THnHelper";
}
};
void df018_customActions()
{
// We enable implicit parallelism
// We create an empty RDataFrame which contains 4 columns filled with random numbers.
// The type of the numbers held by the columns are: double, double, float, int.
auto genD = []() { return gRandom->Uniform(-5, 5); };
auto genF = [&genD]() { return (float)genD(); };
auto genI = [&genD]() { return (int)genD(); };
auto dd = d.Define("x0", genD).Define("x1", genD).Define("x2", genF).Define("x3", genI);
// Our Helper type: templated on the internal THnT type, the size, the types of the columns
// we'll use to fill.
using Helper_t = THnHelper<float, 4>;
Helper_t helper{"myThN", // Name
"A THn with 4 dimensions", // Title
{4, 4, 8, 2}, // NBins
{-10., -10, -4., -6.}, // Axes min values
{10., 10, 5., 7.}}; // Axes max values
// We book the action: it will be treated during the event loop.
auto myTHnT = dd.Book<double, double, float, int>(std::move(helper), {"x0", "x1", "x2", "x3"});
myTHnT->Print();
}
double
typedef void(GLAPIENTRYP _GLUfuncptr)(void)
#define d(i)
Definition RSha256.hxx:102
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
char name[80]
Definition TGX11.cxx:110
float xmax
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTree,...
Templated implementation of the abstract base THn.
Definition THn.h:219
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition TRandom.cxx:672
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition TTreeReader.h:44
CPYCPPYY_EXTERN bool Exec(const std::string &cmd)
Definition API.cxx:333
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
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:527
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:558
UInt_t GetThreadPoolSize()
Returns the size of ROOT's thread pool.
Definition TROOT.cxx:565
TSeq< unsigned int > TSeqU
Definition TSeq.hxx:202
void Initialize(Bool_t useTMVAStyle=kTRUE)
Definition tmvaglob.cxx:176
THnT<float> (*0x55ff9fd5c290): "myThN" "A THn with 4 dimensions"
4 dimensions, 128 entries in 1440 filled bins
Date
April 2018
Authors
Enrico Guiraud, Danilo Piparo (CERN)

Definition in file df018_customActions.C.