Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
InterfaceUtils.hxx
Go to the documentation of this file.
1// Author: Enrico Guiraud, Danilo Piparo CERN 02/2018
2
3/*************************************************************************
4 * Copyright (C) 1995-2018, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11#ifndef ROOT_RDF_TINTERFACE_UTILS
12#define ROOT_RDF_TINTERFACE_UTILS
13
14#include "RColumnRegister.hxx"
15#include <ROOT/RDF/RAction.hxx>
16#include <ROOT/RDF/ActionHelpers.hxx> // for BuildAction
18#include <ROOT/RDF/RDefine.hxx>
20#include <ROOT/RDF/RFilter.hxx>
21#include <ROOT/RDF/Utils.hxx>
27#include <string_view>
29#include <ROOT/TypeTraits.hxx>
30#include <TError.h> // gErrorIgnoreLevel
31#include <TH1.h>
32#include <TROOT.h> // IsImplicitMTEnabled
33
34#include <deque>
35#include <functional>
36#include <map>
37#include <memory>
38#include <string>
39#include <type_traits>
40#include <typeinfo>
41#include <vector>
42#include <unordered_map>
43
44class TObjArray;
45class TTree;
46namespace ROOT {
47namespace Detail {
48namespace RDF {
49class RNodeBase;
50}
51}
52namespace RDF {
53template <typename T>
54class RResultPtr;
55template<typename T, typename V>
56class RInterface;
58class RDataSource;
59} // namespace RDF
60
61} // namespace ROOT
62
63/// \cond HIDDEN_SYMBOLS
64
65namespace ROOT {
66namespace Internal {
67namespace RDF {
68using namespace ROOT::Detail::RDF;
69using namespace ROOT::RDF;
70namespace TTraits = ROOT::TypeTraits;
71
72std::string DemangleTypeIdName(const std::type_info &typeInfo);
73
74ColumnNames_t
75ConvertRegexToColumns(const ColumnNames_t &colNames, std::string_view columnNameRegexp, std::string_view callerName);
76
77/// An helper object that sets and resets gErrorIgnoreLevel via RAII.
79private:
81
82public:
85};
86
87/****** BuildAction overloads *******/
88
89// clang-format off
90/// This namespace defines types to be used for tag dispatching in RInterface.
91namespace ActionTags {
92struct Histo1D{};
93struct Histo2D{};
94struct Histo3D{};
95struct HistoND{};
96struct Graph{};
97struct GraphAsymmErrors{};
98struct Profile1D{};
99struct Profile2D{};
100struct Min{};
101struct Max{};
102struct Sum{};
103struct Mean{};
104struct Fill{};
105struct StdDev{};
106struct Display{};
107struct Snapshot{};
108struct Book{};
109}
110// clang-format on
111
112template <typename T, bool ISV6HISTO = std::is_base_of<TH1, std::decay_t<T>>::value>
113struct HistoUtils {
114 static void SetCanExtendAllAxes(T &h) { h.SetCanExtend(::TH1::kAllAxes); }
115 static bool HasAxisLimits(T &h)
116 {
117 auto xaxis = h.GetXaxis();
118 return !(xaxis->GetXmin() == 0. && xaxis->GetXmax() == 0.);
119 }
120};
121
122template <typename T>
123struct HistoUtils<T, false> {
124 static void SetCanExtendAllAxes(T &) {}
125 static bool HasAxisLimits(T &) { return true; }
126};
127
128// Generic filling (covers Histo2D, Histo3D, HistoND, Profile1D and Profile2D actions, with and without weights)
129template <typename... ColTypes, typename ActionTag, typename ActionResultType, typename PrevNodeType>
130std::unique_ptr<RActionBase>
131BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &h, const unsigned int nSlots,
132 std::shared_ptr<PrevNodeType> prevNode, ActionTag, const RColumnRegister &colRegister)
133{
135 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
136 return std::make_unique<Action_t>(Helper_t(h, nSlots), bl, std::move(prevNode), colRegister);
137}
138
139// Histo1D filling (must handle the special case of distinguishing FillHelper and BufferedFillHelper
140template <typename... ColTypes, typename PrevNodeType>
141std::unique_ptr<RActionBase>
142BuildAction(const ColumnNames_t &bl, const std::shared_ptr<::TH1D> &h, const unsigned int nSlots,
143 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Histo1D, const RColumnRegister &colRegister)
144{
146
149 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
150 return std::make_unique<Action_t>(Helper_t(h, nSlots), bl, std::move(prevNode), colRegister);
151 } else {
153 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
154 return std::make_unique<Action_t>(Helper_t(h, nSlots), bl, std::move(prevNode), colRegister);
155 }
156}
157
158template <typename... ColTypes, typename PrevNodeType>
159std::unique_ptr<RActionBase>
160BuildAction(const ColumnNames_t &bl, const std::shared_ptr<TGraph> &g, const unsigned int nSlots,
161 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Graph, const RColumnRegister &colRegister)
162{
164 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
165 return std::make_unique<Action_t>(Helper_t(g, nSlots), bl, std::move(prevNode), colRegister);
166}
167
168template <typename... ColTypes, typename PrevNodeType>
169std::unique_ptr<RActionBase>
170BuildAction(const ColumnNames_t &bl, const std::shared_ptr<TGraphAsymmErrors> &g, const unsigned int nSlots,
171 std::shared_ptr<PrevNodeType> prevNode, ActionTags::GraphAsymmErrors, const RColumnRegister &colRegister)
172{
174 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
175 return std::make_unique<Action_t>(Helper_t(g, nSlots), bl, std::move(prevNode), colRegister);
176}
177
178// Min action
179template <typename ColType, typename PrevNodeType, typename ActionResultType>
180std::unique_ptr<RActionBase>
181BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &minV, const unsigned int nSlots,
182 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Min, const RColumnRegister &colRegister)
183{
186 return std::make_unique<Action_t>(Helper_t(minV, nSlots), bl, std::move(prevNode), colRegister);
187}
188
189// Max action
190template <typename ColType, typename PrevNodeType, typename ActionResultType>
191std::unique_ptr<RActionBase>
192BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &maxV, const unsigned int nSlots,
193 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Max, const RColumnRegister &colRegister)
194{
197 return std::make_unique<Action_t>(Helper_t(maxV, nSlots), bl, std::move(prevNode), colRegister);
198}
199
200// Sum action
201template <typename ColType, typename PrevNodeType, typename ActionResultType>
202std::unique_ptr<RActionBase>
203BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &sumV, const unsigned int nSlots,
204 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Sum, const RColumnRegister &colRegister)
205{
208 return std::make_unique<Action_t>(Helper_t(sumV, nSlots), bl, std::move(prevNode), colRegister);
209}
210
211// Mean action
212template <typename ColType, typename PrevNodeType>
213std::unique_ptr<RActionBase>
214BuildAction(const ColumnNames_t &bl, const std::shared_ptr<double> &meanV, const unsigned int nSlots,
215 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Mean, const RColumnRegister &colRegister)
216{
217 using Helper_t = MeanHelper;
219 return std::make_unique<Action_t>(Helper_t(meanV, nSlots), bl, std::move(prevNode), colRegister);
220}
221
222// Standard Deviation action
223template <typename ColType, typename PrevNodeType>
224std::unique_ptr<RActionBase>
225BuildAction(const ColumnNames_t &bl, const std::shared_ptr<double> &stdDeviationV, const unsigned int nSlots,
226 std::shared_ptr<PrevNodeType> prevNode, ActionTags::StdDev, const RColumnRegister &colRegister)
227{
228 using Helper_t = StdDevHelper;
230 return std::make_unique<Action_t>(Helper_t(stdDeviationV, nSlots), bl, prevNode, colRegister);
231}
232
233using displayHelperArgs_t = std::pair<size_t, std::shared_ptr<ROOT::RDF::RDisplay>>;
234
235// Display action
236template <typename... ColTypes, typename PrevNodeType>
237std::unique_ptr<RActionBase>
238BuildAction(const ColumnNames_t &bl, const std::shared_ptr<displayHelperArgs_t> &helperArgs, const unsigned int,
239 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Display, const RColumnRegister &colRegister)
240{
242 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
243 return std::make_unique<Action_t>(Helper_t(helperArgs->first, helperArgs->second, prevNode), bl, prevNode,
245}
246
247struct SnapshotHelperArgs {
248 std::string fFileName;
249 std::string fDirName;
250 std::string fTreeName;
251 std::vector<std::string> fOutputColNames;
253 RDFDetail::RLoopManager *fLoopManager;
254 bool fToNTuple;
255};
256
257// SnapshotTTree action
258template <typename... ColTypes, typename PrevNodeType>
259std::unique_ptr<RActionBase>
260BuildAction(const ColumnNames_t &colNames, const std::shared_ptr<SnapshotHelperArgs> &snapHelperArgs,
261 const unsigned int nSlots, std::shared_ptr<PrevNodeType> prevNode, ActionTags::Snapshot,
262 const RColumnRegister &colRegister)
263{
264 const auto &filename = snapHelperArgs->fFileName;
265 const auto &dirname = snapHelperArgs->fDirName;
266 const auto &treename = snapHelperArgs->fTreeName;
267 const auto &outputColNames = snapHelperArgs->fOutputColNames;
268 const auto &options = snapHelperArgs->fOptions;
269
270 auto sz = sizeof...(ColTypes);
271 std::vector<bool> isDefine(sz);
272 for (auto i = 0u; i < sz; ++i)
273 isDefine[i] = colRegister.IsDefineOrAlias(colNames[i]);
274
275 std::unique_ptr<RActionBase> actionPtr;
276 if (snapHelperArgs->fToNTuple) {
277#ifdef R__HAS_ROOT7
279 // single-thread snapshot
282
283 auto loopManager = snapHelperArgs->fLoopManager;
284
285 actionPtr.reset(new Action_t(
287 colNames, prevNode, colRegister));
288 } else {
289 // multi-thread snapshot to RNTuple is not yet supported
290 // TODO(fdegeus) Add MT snapshotting
291 throw std::runtime_error("Snapshot: Snapshotting to RNTuple with IMT enabled is not supported yet.");
292 }
293
294 return actionPtr;
295#else
296 throw std::runtime_error(
297 "RDataFrame: Cannot snapshot to RNTuple - this installation of ROOT has not been build with ROOT7 "
298 "components enabled.");
299#endif
300 } else {
302 // single-thread snapshot
305 actionPtr.reset(
307 colNames, prevNode, colRegister));
308 } else {
309 // multi-thread snapshot
312 actionPtr.reset(new Action_t(
314 colNames, prevNode, colRegister));
315 }
316 }
317 return actionPtr;
318}
319
320// Book with custom helper type
321template <typename... ColTypes, typename PrevNodeType, typename Helper_t>
322std::unique_ptr<RActionBase>
323BuildAction(const ColumnNames_t &bl, const std::shared_ptr<Helper_t> &h, const unsigned int /*nSlots*/,
324 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Book, const RColumnRegister &colRegister)
325{
326 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<ColTypes...>>;
327 return std::make_unique<Action_t>(Helper_t(std::move(*h)), bl, std::move(prevNode), colRegister);
328}
329
330/****** end BuildAndBook ******/
331
332template <typename Filter>
333void CheckFilter(Filter &)
334{
335 using FilterRet_t = typename RDF::CallableTraits<Filter>::ret_type;
336 static_assert(std::is_convertible<FilterRet_t, bool>::value,
337 "filter expression returns a type that is not convertible to bool");
338}
339
340ColumnNames_t FilterArraySizeColNames(const ColumnNames_t &columnNames, const std::string &action);
341
342void CheckValidCppVarName(std::string_view var, const std::string &where);
343
344void CheckForRedefinition(const std::string &where, std::string_view definedCol, const RColumnRegister &colRegister,
345 const ColumnNames_t &treeColumns, const ColumnNames_t &dataSourceColumns);
346
347void CheckForDefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister,
348 const ColumnNames_t &treeColumns, const ColumnNames_t &dataSourceColumns);
349
350void CheckForNoVariations(const std::string &where, std::string_view definedColView,
351 const RColumnRegister &colRegister);
352
353std::string PrettyPrintAddr(const void *const addr);
354
355std::shared_ptr<RJittedFilter> BookFilterJit(std::shared_ptr<RNodeBase> *prevNodeOnHeap, std::string_view name,
356 std::string_view expression, const ColumnNames_t &branches,
357 const RColumnRegister &colRegister, TTree *tree, RDataSource *ds);
358
359std::shared_ptr<RJittedDefine> BookDefineJit(std::string_view name, std::string_view expression, RLoopManager &lm,
360 RDataSource *ds, const RColumnRegister &colRegister,
361 const ColumnNames_t &branches, std::shared_ptr<RNodeBase> *prevNodeOnHeap);
362
363std::shared_ptr<RJittedDefine> BookDefinePerSampleJit(std::string_view name, std::string_view expression,
364 RLoopManager &lm, const RColumnRegister &colRegister,
365 std::shared_ptr<RNodeBase> *upcastNodeOnHeap);
366
367std::shared_ptr<RJittedVariation>
368BookVariationJit(const std::vector<std::string> &colNames, std::string_view variationName,
369 const std::vector<std::string> &variationTags, std::string_view expression, RLoopManager &lm,
370 RDataSource *ds, const RColumnRegister &colRegister, const ColumnNames_t &branches,
371 std::shared_ptr<RNodeBase> *upcastNodeOnHeap, bool isSingleColumn);
372
373std::string JitBuildAction(const ColumnNames_t &bl, std::shared_ptr<RDFDetail::RNodeBase> *prevNode,
374 const std::type_info &art, const std::type_info &at, void *rOnHeap, TTree *tree,
375 const unsigned int nSlots, const RColumnRegister &colRegister, RDataSource *ds,
376 std::weak_ptr<RJittedAction> *jittedActionOnHeap, const bool vector2RVec = true);
377
378// Allocate a weak_ptr on the heap, return a pointer to it. The user is responsible for deleting this weak_ptr.
379// This function is meant to be used by RInterface's methods that book code for jitting.
380// The problem it solves is that we generate code to be lazily jitted with the addresses of certain objects in them,
381// and we need to check those objects are still alive when the generated code is finally jitted and executed.
382// So we pass addresses to weak_ptrs allocated on the heap to the jitted code, which is then responsible for
383// the deletion of the weak_ptr object.
384template <typename T>
385std::weak_ptr<T> *MakeWeakOnHeap(const std::shared_ptr<T> &shPtr)
386{
387 return new std::weak_ptr<T>(shPtr);
388}
389
390// Same as MakeWeakOnHeap, but create a shared_ptr that makes sure the object is definitely kept alive.
391template <typename T>
392std::shared_ptr<T> *MakeSharedOnHeap(const std::shared_ptr<T> &shPtr)
393{
394 return new std::shared_ptr<T>(shPtr);
395}
396
397bool AtLeastOneEmptyString(const std::vector<std::string_view> strings);
398
399/// Take a shared_ptr<AnyNodeType> and return a shared_ptr<RNodeBase>.
400/// This works for RLoopManager nodes as well as filters and ranges.
401std::shared_ptr<RNodeBase> UpcastNode(std::shared_ptr<RNodeBase> ptr);
402
403ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns,
404 const RColumnRegister &validDefines, RDataSource *ds);
405
406std::vector<std::string> GetValidatedArgTypes(const ColumnNames_t &colNames, const RColumnRegister &colRegister,
407 TTree *tree, RDataSource *ds, const std::string &context,
408 bool vector2RVec);
409
410std::vector<bool> FindUndefinedDSColumns(const ColumnNames_t &requestedCols, const ColumnNames_t &definedDSCols);
411
412template <typename T>
413void AddDSColumnsHelper(const std::string &colName, RLoopManager &lm, RDataSource &ds, RColumnRegister &colRegister)
414{
415 if (colRegister.IsDefineOrAlias(colName) || !ds.HasColumn(colName) ||
416 lm.HasDataSourceColumnReaders(colName, typeid(T)))
417 return;
418
419 const auto nSlots = lm.GetNSlots();
420 std::vector<std::unique_ptr<RColumnReaderBase>> colReaders;
421 colReaders.reserve(nSlots);
422
423 const auto valuePtrs = ds.GetColumnReaders<T>(colName);
424 if (!valuePtrs.empty()) { // we are using the old GetColumnReaders mechanism in this RDataSource
425 for (auto *ptr : valuePtrs)
426 colReaders.emplace_back(new RDSColumnReader<T>(ptr));
427
428 } else { // using the new GetColumnReaders mechanism
429 // TODO consider changing the interface so we return all of these for all slots in one go
430 for (auto slot = 0u; slot < lm.GetNSlots(); ++slot)
431 colReaders.emplace_back(ds.GetColumnReaders(slot, colName, typeid(T)));
432 }
433
434 lm.AddDataSourceColumnReaders(colName, std::move(colReaders), typeid(T));
435}
436
437/// Take list of column names that must be defined, current map of custom columns, current list of defined column names,
438/// and return a new map of custom columns (with the new datasource columns added to it)
439template <typename... ColumnTypes>
440void AddDSColumns(const std::vector<std::string> &requiredCols, RLoopManager &lm, RDataSource &ds,
442{
443 // hack to expand a template parameter pack without c++17 fold expressions.
444 using expander = int[];
445 int i = 0;
447}
448
449// this function is meant to be called by the jitted code generated by BookFilterJit
450template <typename F, typename PrevNode>
451void JitFilterHelper(F &&f, const char **colsPtr, std::size_t colsSize, std::string_view name,
452 std::weak_ptr<RJittedFilter> *wkJittedFilter, std::shared_ptr<PrevNode> *prevNodeOnHeap,
453 RColumnRegister *colRegister) noexcept
454{
455 if (wkJittedFilter->expired()) {
456 // The branch of the computation graph that needed this jitted code went out of scope between the type
457 // jitting was booked and the time jitting actually happened. Nothing to do other than cleaning up.
458 delete wkJittedFilter;
459 delete colRegister;
460 delete prevNodeOnHeap;
461 return;
462 }
463
464 const ColumnNames_t cols(colsPtr, colsPtr + colsSize);
465 delete[] colsPtr;
466
467 const auto jittedFilter = wkJittedFilter->lock();
468
469 // mock Filter logic -- validity checks and Define-ition of RDataSource columns
470 using Callable_t = std::decay_t<F>;
472 using ColTypes_t = typename TTraits::CallableTraits<Callable_t>::arg_types;
473 constexpr auto nColumns = ColTypes_t::list_size;
474 CheckFilter(f);
475
476 auto &lm = *jittedFilter->GetLoopManagerUnchecked(); // RLoopManager must exist at this time
477 auto ds = lm.GetDataSource();
478
479 if (ds != nullptr)
481
482 jittedFilter->SetFilter(
483 std::unique_ptr<RFilterBase>(new F_t(std::forward<F>(f), cols, *prevNodeOnHeap, *colRegister, name)));
484 // colRegister points to the columns structure in the heap, created before the jitted call so that the jitter can
485 // share data after it has lazily compiled the code. Here the data has been used and the memory can be freed.
486 delete colRegister;
487 delete prevNodeOnHeap;
488 delete wkJittedFilter;
489}
490
491namespace DefineTypes {
492struct RDefineTag {};
493struct RDefinePerSampleTag {};
494}
495
496template <typename F>
497auto MakeDefineNode(DefineTypes::RDefineTag, std::string_view name, std::string_view dummyType, F &&f,
498 const ColumnNames_t &cols, RColumnRegister &colRegister, RLoopManager &lm)
499{
500 return std::unique_ptr<RDefineBase>(new RDefine<std::decay_t<F>, ExtraArgsForDefine::None>(
501 name, dummyType, std::forward<F>(f), cols, colRegister, lm));
502}
503
504template <typename F>
505auto MakeDefineNode(DefineTypes::RDefinePerSampleTag, std::string_view name, std::string_view dummyType, F &&f,
506 const ColumnNames_t &, RColumnRegister &, RLoopManager &lm)
507{
508 return std::unique_ptr<RDefineBase>(
509 new RDefinePerSample<std::decay_t<F>>(name, dummyType, std::forward<F>(f), lm));
510}
511
512// Build a RDefine or a RDefinePerSample object and attach it to an existing RJittedDefine
513// This function is meant to be called by jitted code right before starting the event loop.
514// If colsPtr is null, build a RDefinePerSample (it has no input columns), otherwise a RDefine.
515template <typename RDefineTypeTag, typename F>
516void JitDefineHelper(F &&f, const char **colsPtr, std::size_t colsSize, std::string_view name, RLoopManager *lm,
517 std::weak_ptr<RJittedDefine> *wkJittedDefine, RColumnRegister *colRegister,
518 std::shared_ptr<RNodeBase> *prevNodeOnHeap) noexcept
519{
520 // a helper to delete objects allocated before jitting, so that the jitter can share data with lazily jitted code
521 auto doDeletes = [&] {
522 delete wkJittedDefine;
523 delete colRegister;
524 delete prevNodeOnHeap;
525 delete[] colsPtr;
526 };
527
528 if (wkJittedDefine->expired()) {
529 // The branch of the computation graph that needed this jitted code went out of scope between the type
530 // jitting was booked and the time jitting actually happened. Nothing to do other than cleaning up.
531 doDeletes();
532 return;
533 }
534
535 const ColumnNames_t cols(colsPtr, colsPtr + colsSize);
536
537 auto jittedDefine = wkJittedDefine->lock();
538
539 using Callable_t = std::decay_t<F>;
540 using ColTypes_t = typename TTraits::CallableTraits<Callable_t>::arg_types;
541
542 auto ds = lm->GetDataSource();
543 if (ds != nullptr)
545
546 // will never actually be used (trumped by jittedDefine->GetTypeName()), but we set it to something meaningful
547 // to help devs debugging
548 const auto dummyType = "jittedCol_t";
549 // use unique_ptr<RDefineBase> instead of make_unique<NewCol_t> to reduce jit/compile-times
550 std::unique_ptr<RDefineBase> newCol{
551 MakeDefineNode(RDefineTypeTag{}, name, dummyType, std::forward<F>(f), cols, *colRegister, *lm)};
552 jittedDefine->SetDefine(std::move(newCol));
553
554 doDeletes();
555}
556
557template <bool IsSingleColumn, typename F>
558void JitVariationHelper(F &&f, const char **colsPtr, std::size_t colsSize, const char **variedCols,
559 std::size_t variedColsSize, const char **variationTags, std::size_t variationTagsSize,
560 std::string_view variationName, RLoopManager *lm,
561 std::weak_ptr<RJittedVariation> *wkJittedVariation, RColumnRegister *colRegister,
562 std::shared_ptr<RNodeBase> *prevNodeOnHeap) noexcept
563{
564 // a helper to delete objects allocated before jitting, so that the jitter can share data with lazily jitted code
565 auto doDeletes = [&] {
566 delete[] colsPtr;
567 delete[] variedCols;
568 delete[] variationTags;
569
570 delete wkJittedVariation;
571 delete colRegister;
572 delete prevNodeOnHeap;
573 };
574
575 if (wkJittedVariation->expired()) {
576 // The branch of the computation graph that needed this jitted variation went out of scope between the type
577 // jitting was booked and the time jitting actually happened. Nothing to do other than cleaning up.
578 doDeletes();
579 return;
580 }
581
582 const ColumnNames_t inputColNames(colsPtr, colsPtr + colsSize);
583 std::vector<std::string> variedColNames(variedCols, variedCols + variedColsSize);
584 std::vector<std::string> tags(variationTags, variationTags + variationTagsSize);
585
586 auto jittedVariation = wkJittedVariation->lock();
587
588 using Callable_t = std::decay_t<F>;
589 using ColTypes_t = typename TTraits::CallableTraits<Callable_t>::arg_types;
590
591 auto ds = lm->GetDataSource();
592 if (ds != nullptr)
594
595 // use unique_ptr<RDefineBase> instead of make_unique<NewCol_t> to reduce jit/compile-times
596 std::unique_ptr<RVariationBase> newVariation{new RVariation<std::decay_t<F>, IsSingleColumn>(
597 std::move(variedColNames), variationName, std::forward<F>(f), std::move(tags), jittedVariation->GetTypeName(),
599 jittedVariation->SetVariation(std::move(newVariation));
600
601 doDeletes();
602}
603
604/// Convenience function invoked by jitted code to build action nodes at runtime
605template <typename ActionTag, typename... ColTypes, typename PrevNodeType, typename HelperArgType>
606void CallBuildAction(std::shared_ptr<PrevNodeType> *prevNodeOnHeap, const char **colsPtr, std::size_t colsSize,
607 const unsigned int nSlots, std::shared_ptr<HelperArgType> *helperArgOnHeap,
608 std::weak_ptr<RJittedAction> *wkJittedActionOnHeap, RColumnRegister *colRegister) noexcept
609{
610 // a helper to delete objects allocated before jitting, so that the jitter can share data with lazily jitted code
611 auto doDeletes = [&] {
612 delete[] colsPtr;
613 delete helperArgOnHeap;
615 // colRegister must be deleted before prevNodeOnHeap because their dtor needs the RLoopManager to be alive
616 // and prevNodeOnHeap is what keeps it alive if the rest of the computation graph is already out of scope
617 delete colRegister;
618 delete prevNodeOnHeap;
619 };
620
621 if (wkJittedActionOnHeap->expired()) {
622 // The branch of the computation graph that needed this jitted variation went out of scope between the type
623 // jitting was booked and the time jitting actually happened. Nothing to do other than cleaning up.
624 doDeletes();
625 return;
626 }
627
628 const ColumnNames_t cols(colsPtr, colsPtr + colsSize);
629
631
632 // if we are here it means we are jitting, if we are jitting the loop manager must be alive
634 auto &loopManager = *prevNodePtr->GetLoopManagerUnchecked();
635 using ColTypes_t = TypeList<ColTypes...>;
636 constexpr auto nColumns = ColTypes_t::list_size;
637 auto ds = loopManager.GetDataSource();
638 if (ds != nullptr)
640
641 auto actionPtr = BuildAction<ColTypes...>(cols, std::move(*helperArgOnHeap), nSlots, std::move(prevNodePtr),
643 jittedActionOnHeap->SetAction(std::move(actionPtr));
644
645 doDeletes();
646}
647
648/// The contained `type` alias is `double` if `T == RInferredType`, `U` if `T == std::container<U>`, `T` otherwise.
649template <typename T, bool Container = IsDataContainer<T>::value && !std::is_same<T, std::string>::value>
650struct RMinReturnType {
651 using type = T;
652};
653
654template <>
656 using type = double;
657};
658
659template <typename T>
660struct RMinReturnType<T, true> {
661 using type = TTraits::TakeFirstParameter_t<T>;
662};
663
664// return wrapper around f that prepends an `unsigned int slot` parameter
665template <typename R, typename F, typename... Args>
666std::function<R(unsigned int, Args...)> AddSlotParameter(F &f, TypeList<Args...>)
667{
668 return [f](unsigned int, Args... a) mutable -> R { return f(a...); };
669}
670
671template <typename ColType, typename... Rest>
672struct RNeedJittingHelper {
673 static constexpr bool value = RNeedJittingHelper<Rest...>::value;
674};
675
676template <typename... Rest>
678 static constexpr bool value = true;
679};
680
681template <typename T>
682struct RNeedJittingHelper<T> {
683 static constexpr bool value = false;
684};
685
686template <>
688 static constexpr bool value = true;
689};
690
691template <typename ...ColTypes>
692struct RNeedJitting {
693 static constexpr bool value = RNeedJittingHelper<ColTypes...>::value;
694};
695
696template <>
697struct RNeedJitting<> {
698 static constexpr bool value = false;
699};
700
701///////////////////////////////////////////////////////////////////////////////
702/// Check preconditions for RInterface::Aggregate:
703/// - the aggregator callable must have signature `U(U,T)` or `void(U&,T)`.
704/// - the merge callable must have signature `U(U,U)` or `void(std::vector<U>&)`
705template <typename R, typename Merge, typename U, typename T, typename decayedU = std::decay_t<U>,
706 typename mergeArgsNoDecay_t = typename CallableTraits<Merge>::arg_types_nodecay,
707 typename mergeArgs_t = typename CallableTraits<Merge>::arg_types,
708 typename mergeRet_t = typename CallableTraits<Merge>::ret_type>
710{
711 constexpr bool isAggregatorOk =
712 (std::is_same<R, decayedU>::value) || (std::is_same<R, void>::value && std::is_lvalue_reference<U>::value);
713 static_assert(isAggregatorOk, "aggregator function must have signature `U(U,T)` or `void(U&,T)`");
714 constexpr bool isMergeOk =
715 (std::is_same<TypeList<decayedU, decayedU>, mergeArgs_t>::value && std::is_same<decayedU, mergeRet_t>::value) ||
716 (std::is_same<TypeList<std::vector<decayedU> &>, mergeArgsNoDecay_t>::value &&
717 std::is_same<void, mergeRet_t>::value);
718 static_assert(isMergeOk, "merge function must have signature `U(U,U)` or `void(std::vector<U>&)`");
719}
720
721///////////////////////////////////////////////////////////////////////////////
722/// This overload of CheckAggregate is called when the aggregator takes more than two arguments
723template <typename R, typename T>
724void CheckAggregate(T)
725{
726 static_assert(sizeof(T) == 0, "aggregator function must take exactly two arguments");
727}
728
729///////////////////////////////////////////////////////////////////////////////
730/// Check as many template parameters were passed as the number of column names, throw if this is not the case.
731void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames);
732
733/// Return local BranchNames or default BranchNames according to which one should be used
734const ColumnNames_t SelectColumns(unsigned int nArgs, const ColumnNames_t &bl, const ColumnNames_t &defBl);
735
736/// Check whether column names refer to a valid branch of a TTree or have been `Define`d. Return invalid column names.
737ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns,
738 const RColumnRegister &definedCols, const ColumnNames_t &dataSourceColumns);
739
740/// Returns the list of Filters defined in the whole graph
741std::vector<std::string> GetFilterNames(const std::shared_ptr<RLoopManager> &loopManager);
742
743/// Returns the list of Filters defined in the branch
744template <typename NodeType>
745std::vector<std::string> GetFilterNames(const std::shared_ptr<NodeType> &node)
746{
747 std::vector<std::string> filterNames;
748 node->AddFilterName(filterNames);
749 return filterNames;
750}
751
752struct ParsedTreePath {
753 std::string fTreeName;
754 std::string fDirName;
755};
756
758
759// Check if a condition is true for all types
760template <bool...>
761struct TBoolPack;
762
763template <bool... bs>
764using IsTrueForAllImpl_t = typename std::is_same<TBoolPack<bs..., true>, TBoolPack<true, bs...>>;
765
766template <bool... Conditions>
767struct TEvalAnd {
768 static constexpr bool value = IsTrueForAllImpl_t<Conditions...>::value;
769};
770
771// Check if a class is a specialisation of stl containers templates
772// clang-format off
773
774template <typename>
775struct IsList_t : std::false_type {};
776
777template <typename T>
778struct IsList_t<std::list<T>> : std::true_type {};
779
780template <typename>
781struct IsDeque_t : std::false_type {};
782
783template <typename T>
784struct IsDeque_t<std::deque<T>> : std::true_type {};
785// clang-format on
786
787void CheckForDuplicateSnapshotColumns(const ColumnNames_t &cols);
788
789template <typename T>
790struct InnerValueType {
791 using type = T; // fallback for when T is not a nested RVec
792};
793
794template <typename Elem>
795struct InnerValueType<ROOT::VecOps::RVec<ROOT::VecOps::RVec<Elem>>> {
796 using type = Elem;
797};
798
799template <typename T>
801
802std::pair<std::vector<std::string>, std::vector<std::string>>
803AddSizeBranches(const std::vector<std::string> &branches, TTree *tree, std::vector<std::string> &&colsWithoutAliases,
804 std::vector<std::string> &&colsWithAliases);
805
806void RemoveDuplicates(ColumnNames_t &columnNames);
807
808#ifdef R__HAS_ROOT7
809void RemoveRNTupleSubFields(ColumnNames_t &columnNames);
810#endif
811
812} // namespace RDF
813} // namespace Internal
814
815namespace Detail {
816namespace RDF {
817
818/// The aliased type is `double` if `T == RInferredType`, `U` if `T == container<U>`, `T` otherwise.
819template <typename T>
820using MinReturnType_t = typename RDFInternal::RMinReturnType<T>::type;
821
822template <typename T>
824
825template <typename T>
827
828} // namespace RDF
829} // namespace Detail
830} // namespace ROOT
831
832/// \endcond
833
834#endif
#define f(i)
Definition RSha256.hxx:104
#define g(i)
Definition RSha256.hxx:105
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Int_t gErrorIgnoreLevel
Error handling routines.
Definition TError.cxx:31
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
The head node of a RDF computation graph.
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
The public interface to the RDataFrame federation of classes.
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition RVec.hxx:1529
@ kAllAxes
Definition TH1.h:76
An array of TObjects.
Definition TObjArray.h:31
A TTree represents a columnar dataset.
Definition TTree.h:79
T Sum(const RVec< T > &v, const T zero=T(0))
Sum elements of an RVec.
Definition RVec.hxx:1954
#define F(x, y, z)
const ColumnNames_t SelectColumns(unsigned int nRequiredNames, const ColumnNames_t &names, const ColumnNames_t &defaultNames)
Choose between local column names or default column names, throw in case of errors.
void CheckForNoVariations(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister)
Throw if the column has systematic variations attached.
ParsedTreePath ParseTreePath(std::string_view fullTreeName)
void CheckForRedefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister, const ColumnNames_t &treeColumns, const ColumnNames_t &dataSourceColumns)
Throw if column definedColView is already there.
void CheckForDefinition(const std::string &where, std::string_view definedColView, const RColumnRegister &colRegister, const ColumnNames_t &treeColumns, const ColumnNames_t &dataSourceColumns)
Throw if column definedColView is not already there.
std::shared_ptr< RJittedDefine > BookDefineJit(std::string_view name, std::string_view expression, RLoopManager &lm, RDataSource *ds, const RColumnRegister &colRegister, const ColumnNames_t &branches, std::shared_ptr< RNodeBase > *upcastNodeOnHeap)
Book the jitting of a Define call.
void CheckValidCppVarName(std::string_view var, const std::string &where)
void RemoveDuplicates(ColumnNames_t &columnNames)
ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns, const RColumnRegister &colRegister, RDataSource *ds)
Given the desired number of columns and the user-provided list of columns:
std::shared_ptr< RNodeBase > UpcastNode(std::shared_ptr< RNodeBase > ptr)
std::vector< std::string > GetFilterNames(const std::shared_ptr< RLoopManager > &loopManager)
std::string PrettyPrintAddr(const void *const addr)
void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames)
std::string DemangleTypeIdName(const std::type_info &typeInfo)
bool AtLeastOneEmptyString(const std::vector< std::string_view > strings)
std::shared_ptr< RDFDetail::RJittedFilter > BookFilterJit(std::shared_ptr< RDFDetail::RNodeBase > *prevNodeOnHeap, std::string_view name, std::string_view expression, const ColumnNames_t &branches, const RColumnRegister &colRegister, TTree *tree, RDataSource *ds)
Book the jitting of a Filter call.
ColumnNames_t FilterArraySizeColNames(const ColumnNames_t &columnNames, const std::string &action)
Take a list of column names, return that list with entries starting by '#' filtered out.
std::shared_ptr< RJittedVariation > BookVariationJit(const std::vector< std::string > &colNames, std::string_view variationName, const std::vector< std::string > &variationTags, std::string_view expression, RLoopManager &lm, RDataSource *ds, const RColumnRegister &colRegister, const ColumnNames_t &branches, std::shared_ptr< RNodeBase > *upcastNodeOnHeap, bool isSingleColumn)
Book the jitting of a Vary call.
std::vector< std::string > GetValidatedArgTypes(const ColumnNames_t &colNames, const RColumnRegister &colRegister, TTree *tree, RDataSource *ds, const std::string &context, bool vector2RVec)
void CheckForDuplicateSnapshotColumns(const ColumnNames_t &cols)
ColumnNames_t ConvertRegexToColumns(const ColumnNames_t &colNames, std::string_view columnNameRegexp, std::string_view callerName)
std::pair< std::vector< std::string >, std::vector< std::string > > AddSizeBranches(const std::vector< std::string > &branches, TTree *tree, std::vector< std::string > &&colsWithoutAliases, std::vector< std::string > &&colsWithAliases)
Return copies of colsWithoutAliases and colsWithAliases with size branches for variable-sized array b...
std::vector< bool > FindUndefinedDSColumns(const ColumnNames_t &requestedCols, const ColumnNames_t &definedCols)
Return a bitset each element of which indicates whether the corresponding element in selectedColumns ...
std::shared_ptr< RJittedDefine > BookDefinePerSampleJit(std::string_view name, std::string_view expression, RLoopManager &lm, const RColumnRegister &colRegister, std::shared_ptr< RNodeBase > *upcastNodeOnHeap)
Book the jitting of a DefinePerSample call.
std::string JitBuildAction(const ColumnNames_t &cols, std::shared_ptr< RDFDetail::RNodeBase > *prevNode, const std::type_info &helperArgType, const std::type_info &at, void *helperArgOnHeap, TTree *tree, const unsigned int nSlots, const RColumnRegister &colRegister, RDataSource *ds, std::weak_ptr< RJittedAction > *jittedActionOnHeap, const bool vector2RVec)
ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns, const RColumnRegister &definedCols, const ColumnNames_t &dataSourceColumns)
ROOT type_traits extensions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:570
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Mean(Long64_t n, const T *a, const Double_t *w=nullptr)
Returns the weighted mean of an array a with length n.
Definition TMath.h:1169
Double_t StdDev(Long64_t n, const T *a, const Double_t *w=nullptr)
Definition TMath.h:531
A collection of options to steer the creation of the dataset on file.
Lightweight storage for a collection of types.