Logo ROOT   6.16/01
Reference Guide
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 <ROOT/RDF/RAction.hxx>
15#include <ROOT/RDF/ActionHelpers.hxx> // for BuildAction
18#include <ROOT/RDF/RFilter.hxx>
19#include <ROOT/RDF/Utils.hxx>
25#include <ROOT/RMakeUnique.hxx>
26#include <ROOT/RStringView.hxx>
27#include <ROOT/TypeTraits.hxx>
28#include <TError.h> // gErrorIgnoreLevel
29#include <TH1.h>
30
31#include <deque>
32#include <functional>
33#include <map>
34#include <memory>
35#include <string>
36#include <type_traits>
37#include <typeinfo>
38#include <vector>
39
40class TObjArray;
41class TTree;
42namespace ROOT {
43namespace Detail {
44namespace RDF {
45class RNodeBase;
46}
47}
48namespace RDF {
49template <typename T>
50class RResultPtr;
51template<typename T, typename V>
52class RInterface;
54class RDataSource;
55} // namespace RDF
56
57} // namespace ROOT
58
59/// \cond HIDDEN_SYMBOLS
60
61namespace ROOT {
62namespace Internal {
63namespace RDF {
64using namespace ROOT::Detail::RDF;
65using namespace ROOT::RDF;
66namespace TTraits = ROOT::TypeTraits;
68
70HeadNode_t CreateSnaphotRDF(const ColumnNames_t &validCols,
71 const std::string &fullTreeName,
72 const std::string &fileName,
73 bool isLazy,
74 RLoopManager &loopManager,
75 std::unique_ptr<RDFInternal::RActionBase> actionPtr);
76
77std::string DemangleTypeIdName(const std::type_info &typeInfo);
78
80 std::string_view columnNameRegexp,
81 std::string_view callerName);
82
83/// An helper object that sets and resets gErrorIgnoreLevel via RAII.
84class RIgnoreErrorLevelRAII {
85private:
86 int fCurIgnoreErrorLevel = gErrorIgnoreLevel;
87
88public:
89 RIgnoreErrorLevelRAII(int errorIgnoreLevel) { gErrorIgnoreLevel = errorIgnoreLevel; }
90 RIgnoreErrorLevelRAII() { gErrorIgnoreLevel = fCurIgnoreErrorLevel; }
91};
92
93/****** BuildAction overloads *******/
94
95// clang-format off
96/// This namespace defines types to be used for tag dispatching in RInterface.
97namespace ActionTags {
98struct Histo1D{};
99struct Histo2D{};
100struct Histo3D{};
101struct Graph{};
102struct Profile1D{};
103struct Profile2D{};
104struct Min{};
105struct Max{};
106struct Sum{};
107struct Mean{};
108struct Fill{};
109struct StdDev{};
110struct Display{};
111}
112// clang-format on
113
114template <int D, typename P, template <int, typename, template <typename> class> class... S>
115class THist;
116
117/// Check whether a histogram type is a classic or v7 histogram.
118template <typename T>
119struct IsV7Hist : public std::false_type {
120 static_assert(std::is_base_of<TH1, T>::value, "not implemented for this type");
121};
122
123template <int D, typename P, template <int, typename, template <typename> class> class... S>
124struct IsV7Hist<THist<D, P, S...>> : public std::true_type {
125};
126
127template <typename T, bool ISV7HISTO = IsV7Hist<T>::value>
128struct HistoUtils {
129 static void SetCanExtendAllAxes(T &h) { h.SetCanExtend(::TH1::kAllAxes); }
130 static bool HasAxisLimits(T &h)
131 {
132 auto xaxis = h.GetXaxis();
133 return !(xaxis->GetXmin() == 0. && xaxis->GetXmax() == 0.);
134 }
135};
136
137template <typename T>
138struct HistoUtils<T, true> {
139 static void SetCanExtendAllAxes(T &) {}
140 static bool HasAxisLimits(T &) { return true; }
141};
142
143// Generic filling (covers Histo2D, Histo3D, Profile1D and Profile2D actions, with and without weights)
144template <typename... BranchTypes, typename ActionTag, typename ActionResultType, typename PrevNodeType>
145std::unique_ptr<RActionBase>
146BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &h, const unsigned int nSlots,
147 std::shared_ptr<PrevNodeType> prevNode, ActionTag, RDFInternal::RBookedCustomColumns customColumns)
148{
149 using Helper_t = FillParHelper<ActionResultType>;
150 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
151 return std::make_unique<Action_t>(Helper_t(h, nSlots), bl, std::move(prevNode), customColumns);
152}
153
154// Histo1D filling (must handle the special case of distinguishing FillParHelper and FillHelper
155template <typename... BranchTypes, typename PrevNodeType>
156std::unique_ptr<RActionBase> BuildAction(const ColumnNames_t &bl, const std::shared_ptr<::TH1D> &h,
157 const unsigned int nSlots, std::shared_ptr<PrevNodeType> prevNode,
158 ActionTags::Histo1D, RDFInternal::RBookedCustomColumns customColumns)
159{
160 auto hasAxisLimits = HistoUtils<::TH1D>::HasAxisLimits(*h);
161
162 if (hasAxisLimits) {
163 using Helper_t = FillParHelper<::TH1D>;
164 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
165 return std::make_unique<Action_t>(Helper_t(h, nSlots), bl, std::move(prevNode), customColumns);
166 } else {
167 using Helper_t = FillHelper;
168 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
169 return std::make_unique<Action_t>(Helper_t(h, nSlots), bl, std::move(prevNode), customColumns);
170 }
171}
172
173template <typename... BranchTypes, typename PrevNodeType>
174std::unique_ptr<RActionBase>
175BuildAction(const ColumnNames_t &bl, const std::shared_ptr<TGraph> &g, const unsigned int nSlots,
176 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Graph, RDFInternal::RBookedCustomColumns customColumns)
177{
178 using Helper_t = FillTGraphHelper;
179 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
180 return std::make_unique<Action_t>(Helper_t(g, nSlots), bl, std::move(prevNode), customColumns);
181}
182
183// Min action
184template <typename BranchType, typename PrevNodeType, typename ActionResultType>
185std::unique_ptr<RActionBase>
186BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &minV, const unsigned int nSlots,
187 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Min, RDFInternal::RBookedCustomColumns customColumns)
188{
189 using Helper_t = MinHelper<ActionResultType>;
190 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
191 return std::make_unique<Action_t>(Helper_t(minV, nSlots), bl, std::move(prevNode), customColumns);
192}
193
194// Max action
195template <typename BranchType, typename PrevNodeType, typename ActionResultType>
196std::unique_ptr<RActionBase>
197BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &maxV, const unsigned int nSlots,
198 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Max, RDFInternal::RBookedCustomColumns customColumns)
199{
200 using Helper_t = MaxHelper<ActionResultType>;
201 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
202 return std::make_unique<Action_t>(Helper_t(maxV, nSlots), bl, std::move(prevNode), customColumns);
203}
204
205// Sum action
206template <typename BranchType, typename PrevNodeType, typename ActionResultType>
207std::unique_ptr<RActionBase>
208BuildAction(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &sumV, const unsigned int nSlots,
209 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Sum, RDFInternal::RBookedCustomColumns customColumns)
210{
211 using Helper_t = SumHelper<ActionResultType>;
212 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
213 return std::make_unique<Action_t>(Helper_t(sumV, nSlots), bl, std::move(prevNode), customColumns);
214}
215
216// Mean action
217template <typename BranchType, typename PrevNodeType>
218std::unique_ptr<RActionBase>
219BuildAction(const ColumnNames_t &bl, const std::shared_ptr<double> &meanV, const unsigned int nSlots,
220 std::shared_ptr<PrevNodeType> prevNode, ActionTags::Mean, RDFInternal::RBookedCustomColumns customColumns)
221{
222 using Helper_t = MeanHelper;
223 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
224 return std::make_unique<Action_t>(Helper_t(meanV, nSlots), bl, std::move(prevNode), customColumns);
225}
226
227// Standard Deviation action
228template <typename BranchType, typename PrevNodeType>
229std::unique_ptr<RActionBase>
230BuildAction(const ColumnNames_t &bl, const std::shared_ptr<double> &stdDeviationV, const unsigned int nSlots,
231 std::shared_ptr<PrevNodeType> prevNode, ActionTags::StdDev, RDFInternal::RBookedCustomColumns customColumns)
232{
233 using Helper_t = StdDevHelper;
234 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
235 return std::make_unique<Action_t>(Helper_t(stdDeviationV, nSlots), bl, prevNode, customColumns);
236}
237
238// Display action
239template <typename... BranchTypes, typename PrevNodeType>
240std::unique_ptr<RActionBase> BuildAction(const ColumnNames_t &bl, const std::shared_ptr<RDisplay> &d,
241 const unsigned int, std::shared_ptr<PrevNodeType> prevNode,
242 ActionTags::Display, RDFInternal::RBookedCustomColumns customColumns)
243{
244 using Helper_t = DisplayHelper<PrevNodeType>;
245 using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
246 return std::make_unique<Action_t>(Helper_t(d, prevNode), bl, prevNode, customColumns);
247}
248
249/****** end BuildAndBook ******/
250
251template <typename Filter>
252void CheckFilter(Filter &)
253{
254 using FilterRet_t = typename RDF::CallableTraits<Filter>::ret_type;
255 static_assert(std::is_same<FilterRet_t, bool>::value, "filter functions must return a bool");
256}
257
258void CheckCustomColumn(std::string_view definedCol, TTree *treePtr, const ColumnNames_t &customCols,
259 const std::map<std::string, std::string> &aliasMap, const ColumnNames_t &dataSourceColumns);
260
261std::string PrettyPrintAddr(const void *const addr);
262
263void BookFilterJit(RJittedFilter *jittedFilter, void *prevNodeOnHeap, std::string_view name,
264 std::string_view expression, const std::map<std::string, std::string> &aliasMap,
265 const ColumnNames_t &branches, const RDFInternal::RBookedCustomColumns &customCols, TTree *tree,
266 RDataSource *ds, unsigned int namespaceID);
267
269 const std::shared_ptr<RJittedCustomColumn> &jittedCustomColumn,
271
272std::string JitBuildAction(const ColumnNames_t &bl, void *prevNode, const std::type_info &art, const std::type_info &at,
273 void *r, TTree *tree, const unsigned int nSlots,
274 const RDFInternal::RBookedCustomColumns &customColumns, RDataSource *ds,
275 std::shared_ptr<RJittedAction> *jittedActionOnHeap, unsigned int namespaceID);
276
277// allocate a shared_ptr on the heap, return a reference to it. the user is responsible of deleting the shared_ptr*.
278// this function is meant to only be used by RInterface's action methods, and should be deprecated as soon as we find
279// a better way to make jitting work: the problem it solves is that we need to pass the same shared_ptr to the Helper
280// object of each action and to the RResultPtr returned by the action. While the former is only instantiated when
281// the event loop is about to start, the latter has to be returned to the user as soon as the action is booked.
282// a heap allocated shared_ptr will stay alive long enough that at jitting time its address is still valid.
283template <typename T>
284std::shared_ptr<T> *MakeSharedOnHeap(const std::shared_ptr<T> &shPtr)
285{
286 return new std::shared_ptr<T>(shPtr);
287}
288
289bool AtLeastOneEmptyString(const std::vector<std::string_view> strings);
290
291/// Take a shared_ptr<AnyNodeType> and return a shared_ptr<RNodeBase>.
292/// This works for RLoopManager nodes as well as filters and ranges.
293std::shared_ptr<RNodeBase> UpcastNode(std::shared_ptr<RNodeBase> ptr);
294
295ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns,
296 const ColumnNames_t &validCustomColumns, RDataSource *ds);
297
298std::vector<bool> FindUndefinedDSColumns(const ColumnNames_t &requestedCols, const ColumnNames_t &definedDSCols);
299
301
302template <typename T>
303void AddDSColumnsHelper(RLoopManager &lm, std::string_view name, RDFInternal::RBookedCustomColumns &currentCols,
304 RDataSource &ds, unsigned int nSlots)
305{
306 auto readers = ds.GetColumnReaders<T>(name);
307 auto getValue = [readers](unsigned int slot) { return *readers[slot]; };
308 using NewCol_t = RCustomColumn<decltype(getValue), CustomColExtraArgs::Slot>;
309
310 auto newCol = std::make_shared<NewCol_t>(&lm, name, std::move(getValue), ColumnNames_t{}, nSlots, currentCols,
311 /*isDSColumn=*/true);
312
313 lm.RegisterCustomColumn(newCol.get());
314 currentCols.AddName(name);
315 currentCols.AddColumn(newCol, name);
316}
317
318/// Take list of column names that must be defined, current map of custom columns, current list of defined column names,
319/// and return a new map of custom columns (with the new datasource columns added to it)
320template <typename... ColumnTypes, std::size_t... S>
322AddDSColumns(RLoopManager &lm, const std::vector<std::string> &requiredCols,
323 const RDFInternal::RBookedCustomColumns &currentCols, RDataSource &ds, unsigned int nSlots,
325{
326
327 const auto mustBeDefined = FindUndefinedDSColumns(requiredCols, currentCols.GetNames());
328 if (std::none_of(mustBeDefined.begin(), mustBeDefined.end(), [](bool b) { return b; })) {
329 // no need to define any column
330 return currentCols;
331 } else {
332 auto newColumns(currentCols);
333
334 // hack to expand a template parameter pack without c++17 fold expressions.
335 int expander[] = {(mustBeDefined[S] ? AddDSColumnsHelper<ColumnTypes>(lm, requiredCols[S], newColumns, ds, nSlots)
336 : /*no-op*/ ((void)0),
337 0)...,
338 0};
339 (void)expander; // avoid unused variable warnings
340 (void)nSlots; // avoid unused variable warnings
341 return newColumns;
342 }
343}
344
345// this function is meant to be called by the jitted code generated by BookFilterJit
346template <typename F, typename PrevNode>
347void JitFilterHelper(F &&f, const ColumnNames_t &cols, std::string_view name, RJittedFilter *jittedFilter,
348 std::shared_ptr<PrevNode> *prevNodeOnHeap, RDFInternal::RBookedCustomColumns *customColumns)
349{
350 // mock Filter logic -- validity checks and Define-ition of RDataSource columns
351 using F_t = RFilter<F, PrevNode>;
352 using ColTypes_t = typename TTraits::CallableTraits<F>::arg_types;
353 constexpr auto nColumns = ColTypes_t::list_size;
354 RDFInternal::CheckFilter(f);
355
356 auto &lm = *jittedFilter->GetLoopManagerUnchecked(); // RLoopManager must exist at this time
357 auto ds = lm.GetDataSource();
358
359 auto newColumns = ds ? RDFInternal::AddDSColumns(lm, cols, *customColumns, *ds, lm.GetNSlots(),
361 : *customColumns;
362
363 // customColumns points to the columns structure in the heap, created before the jitted call so that the jitter can
364 // share data after it has lazily compiled the code. Here the data has been used and the memory can be freed.
365 delete customColumns;
366
367 jittedFilter->SetFilter(std::make_unique<F_t>(std::move(f), cols, *prevNodeOnHeap, newColumns, name));
368 delete prevNodeOnHeap;
369}
370
371template <typename F>
372void JitDefineHelper(F &&f, const ColumnNames_t &cols, std::string_view name, RLoopManager *lm,
373 RJittedCustomColumn &jittedCustomCol, RDFInternal::RBookedCustomColumns *customColumns)
374{
376 using ColTypes_t = typename TTraits::CallableTraits<F>::arg_types;
377 constexpr auto nColumns = ColTypes_t::list_size;
378
379 auto ds = lm->GetDataSource();
380 auto newColumns = ds ? RDFInternal::AddDSColumns(*lm, cols, *customColumns, *ds, lm->GetNSlots(),
382 : *customColumns;
383
384 // customColumns points to the columns structure in the heap, created before the jitted call so that the jitter can
385 // share data after it has lazily compiled the code. Here the data has been used and the memory can be freed.
386 delete customColumns;
387
388 jittedCustomCol.SetCustomColumn(
389 std::make_unique<NewCol_t>(lm, name, std::move(f), cols, lm->GetNSlots(), newColumns));
390}
391
392/// Convenience function invoked by jitted code to build action nodes at runtime
393template <typename ActionTag, typename... BranchTypes, typename PrevNodeType, typename ActionResultType>
394void CallBuildAction(std::shared_ptr<PrevNodeType> *prevNodeOnHeap, const ColumnNames_t &bl, const unsigned int nSlots,
395 const std::shared_ptr<ActionResultType> *rOnHeap,
396 std::shared_ptr<RJittedAction> *jittedActionOnHeap,
398{
399 // if we are here it means we are jitting, if we are jitting the loop manager must be alive
400 auto &prevNodePtr = *prevNodeOnHeap;
401 auto &loopManager = *prevNodePtr->GetLoopManagerUnchecked();
402 using ColTypes_t = TypeList<BranchTypes...>;
403 constexpr auto nColumns = ColTypes_t::list_size;
404 auto ds = loopManager.GetDataSource();
405 auto newColumns = ds ? RDFInternal::AddDSColumns(loopManager, bl, *customColumns, *ds, loopManager.GetNSlots(),
407 : *customColumns;
408
409 auto actionPtr = BuildAction<BranchTypes...>(bl, *rOnHeap, nSlots, std::move(prevNodePtr), ActionTag{}, newColumns);
410 (*jittedActionOnHeap)->SetAction(std::move(actionPtr));
411
412 // customColumns points to the columns structure in the heap, created before the jitted call so that the jitter can
413 // share data after it has lazily compiled the code. Here the data has been used and the memory can be freed.
414 delete customColumns;
415
416 delete rOnHeap;
417 delete prevNodeOnHeap;
418 delete jittedActionOnHeap;
419}
420
421/// The contained `type` alias is `double` if `T == RInferredType`, `U` if `T == std::container<U>`, `T` otherwise.
422template <typename T, bool Container = TTraits::IsContainer<T>::value && !std::is_same<T, std::string>::value>
423struct TMinReturnType {
424 using type = T;
425};
426
427template <>
428struct TMinReturnType<RInferredType, false> {
429 using type = double;
430};
431
432template <typename T>
433struct TMinReturnType<T, true> {
434 using type = TTraits::TakeFirstParameter_t<T>;
435};
436
437// return wrapper around f that prepends an `unsigned int slot` parameter
438template <typename R, typename F, typename... Args>
439std::function<R(unsigned int, Args...)> AddSlotParameter(F &f, TypeList<Args...>)
440{
441 return [f](unsigned int, Args... a) -> R { return f(a...); };
442}
443
444template <typename BranchType, typename... Rest>
445struct TNeedJitting {
446 static constexpr bool value = TNeedJitting<Rest...>::value;
447};
448
449template <typename... Rest>
450struct TNeedJitting<RInferredType, Rest...> {
451 static constexpr bool value = true;
452};
453
454template <typename T>
455struct TNeedJitting<T> {
456 static constexpr bool value = false;
457};
458
459template <>
460struct TNeedJitting<RInferredType> {
461 static constexpr bool value = true;
462};
463
465
466///////////////////////////////////////////////////////////////////////////////
467/// Check preconditions for RInterface::Aggregate:
468/// - the aggregator callable must have signature `U(U,T)` or `void(U&,T)`.
469/// - the merge callable must have signature `U(U,U)` or `void(std::vector<U>&)`
471 typename mergeArgsNoDecay_t = typename CallableTraits<Merge>::arg_types_nodecay,
472 typename mergeArgs_t = typename CallableTraits<Merge>::arg_types,
473 typename mergeRet_t = typename CallableTraits<Merge>::ret_type>
474void CheckAggregate(TypeList<U, T>)
475{
476 constexpr bool isAggregatorOk =
477 (std::is_same<R, decayedU>::value) || (std::is_same<R, void>::value && std::is_lvalue_reference<U>::value);
478 static_assert(isAggregatorOk, "aggregator function must have signature `U(U,T)` or `void(U&,T)`");
479 constexpr bool isMergeOk =
480 (std::is_same<TypeList<decayedU, decayedU>, mergeArgs_t>::value && std::is_same<decayedU, mergeRet_t>::value) ||
481 (std::is_same<TypeList<std::vector<decayedU> &>, mergeArgsNoDecay_t>::value &&
482 std::is_same<void, mergeRet_t>::value);
483 static_assert(isMergeOk, "merge function must have signature `U(U,U)` or `void(std::vector<U>&)`");
484}
485
486///////////////////////////////////////////////////////////////////////////////
487/// This overload of CheckAggregate is called when the aggregator takes more than two arguments
488template <typename R, typename T>
489void CheckAggregate(T)
490{
491 static_assert(sizeof(T) == 0, "aggregator function must take exactly two arguments");
492}
493
494///////////////////////////////////////////////////////////////////////////////
495/// Check as many template parameters were passed as the number of column names, throw if this is not the case.
496void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames);
497
498/// Return local BranchNames or default BranchNames according to which one should be used
499const ColumnNames_t SelectColumns(unsigned int nArgs, const ColumnNames_t &bl, const ColumnNames_t &defBl);
500
501/// Check whether column names refer to a valid branch of a TTree or have been `Define`d. Return invalid column names.
502ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns,
503 const ColumnNames_t &definedCols, const ColumnNames_t &dataSourceColumns);
504
506
507/// Returns the list of Filters defined in the whole graph
508std::vector<std::string> GetFilterNames(const std::shared_ptr<RLoopManager> &loopManager);
509
510/// Returns the list of Filters defined in the branch
511template <typename NodeType>
512std::vector<std::string> GetFilterNames(const std::shared_ptr<NodeType> &node)
513{
514 std::vector<std::string> filterNames;
515 node->AddFilterName(filterNames);
516 return filterNames;
517}
518
519// Check if a condition is true for all types
520template <bool...>
521struct TBoolPack;
522
523template <bool... bs>
524using IsTrueForAllImpl_t = typename std::is_same<TBoolPack<bs..., true>, TBoolPack<true, bs...>>;
525
526template <bool... Conditions>
527struct TEvalAnd {
528 static constexpr bool value = IsTrueForAllImpl_t<Conditions...>::value;
529};
530
531// Check if a class is a specialisation of stl containers templates
532// clang-format off
533
534template <typename>
535struct IsList_t : std::false_type {};
536
537template <typename T>
538struct IsList_t<std::list<T>> : std::true_type {};
539
540template <typename>
541struct IsDeque_t : std::false_type {};
542
543template <typename T>
544struct IsDeque_t<std::deque<T>> : std::true_type {};
545// clang-format on
546
547} // namespace RDF
548} // namespace Internal
549
550namespace Detail {
551namespace RDF {
552
553/// The aliased type is `double` if `T == RInferredType`, `U` if `T == container<U>`, `T` otherwise.
554template <typename T>
555using MinReturnType_t = typename RDFInternal::TMinReturnType<T>::type;
556
557template <typename T>
558using MaxReturnType_t = MinReturnType_t<T>;
559
560template <typename T>
561using SumReturnType_t = MinReturnType_t<T>;
562
563} // namespace RDF
564} // namespace Detail
565} // namespace ROOT
566
567/// \endcond
568
569#endif
ROOT::R::TRInterface & r
Definition: Object.C:4
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
#define g(i)
Definition: RSha256.hxx:105
#define h(i)
Definition: RSha256.hxx:106
#define R(a, b, c, d, e, f, g, h, i)
Definition: RSha256.hxx:110
R__EXTERN Int_t gErrorIgnoreLevel
Definition: TError.h:105
int type
Definition: TGX11.cxx:120
typedef void((*Func_t)())
A wrapper around a concrete RCustomColumn, which forwards all calls to it RJittedCustomColumn is a pl...
void SetCustomColumn(std::unique_ptr< RCustomColumnBase > c)
A wrapper around a concrete RFilter, which forwards all calls to it RJittedFilter is the type of the ...
void SetFilter(std::unique_ptr< RFilterBase > f)
The head node of a RDF computation graph.
void RegisterCustomColumn(RCustomColumnBase *column)
RDataSource * GetDataSource() const
unsigned int GetNSlots() const
virtual RLoopManager * GetLoopManagerUnchecked()
Definition: RNodeBase.hxx:64
Encapsulates the columns defined by the user.
ColumnNames_t GetNames() const
Returns the list of the names of the defined columns.
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
std::vector< T ** > GetColumnReaders(std::string_view columnName)
Called at most once per column by RDF.
The public interface to the RDataFrame federation of classes.
Definition: RInterface.hxx:87
Smart pointer for the return type of actions.
Definition: RResultPtr.hxx:72
@ kAllAxes
Definition: TH1.h:73
An array of TObjects.
Definition: TObjArray.h:37
A TTree object has a header with a name and a title.
Definition: TTree.h:71
#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.
ColumnNames_t ConvertRegexToColumns(const ROOT::RDF::RNode &node, std::string_view columnNameRegexp, std::string_view callerName)
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)
ColumnNames_t GetTopLevelBranchNames(TTree &t)
Get all the top-level branches names, including the ones of the friend trees.
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)
ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns, const ColumnNames_t &definedCols, const ColumnNames_t &dataSourceColumns)
std::string JitBuildAction(const ColumnNames_t &bl, void *prevNode, const std::type_info &art, const std::type_info &at, void *rOnHeap, TTree *tree, const unsigned int nSlots, const RDFInternal::RBookedCustomColumns &customCols, RDataSource *ds, std::shared_ptr< RJittedAction > *jittedActionOnHeap, unsigned int namespaceID)
HeadNode_t CreateSnaphotRDF(const ColumnNames_t &validCols, const std::string &fullTreeName, const std::string &fileName, bool isLazy, RLoopManager &loopManager, std::unique_ptr< RDFInternal::RActionBase > actionPtr)
bool IsInternalColumn(std::string_view colName)
ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns, const ColumnNames_t &validCustomColumns, RDataSource *ds)
Given the desired number of columns and the user-provided list of columns:
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 ...
void BookDefineJit(std::string_view name, std::string_view expression, RLoopManager &lm, RDataSource *ds, const std::shared_ptr< RJittedCustomColumn > &jittedCustomColumn, const RDFInternal::RBookedCustomColumns &customCols, const ColumnNames_t &branches)
void BookFilterJit(RJittedFilter *jittedFilter, void *prevNodeOnHeap, std::string_view name, std::string_view expression, const std::map< std::string, std::string > &aliasMap, const ColumnNames_t &branches, const RDFInternal::RBookedCustomColumns &customCols, TTree *tree, RDataSource *ds, unsigned int namespaceID)
void CheckCustomColumn(std::string_view definedCol, TTree *treePtr, const ColumnNames_t &customCols, const std::map< std::string, std::string > &aliasMap, const ColumnNames_t &dataSourceColumns)
static double P[]
double T(double x)
Definition: ChebyshevPol.h:34
RInterface<::ROOT::Detail::RDF::RNodeBase, void > RNode
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
ROOT type_traits extensions.
Definition: TypeTraits.hxx:23
T Sum(const RVec< T > &v)
Sum elements of an RVec.
Definition: RVec.hxx:702
RVec< T > Filter(const RVec< T > &v, F &&f)
Create a new collection with the elements passing the filter expressed by the predicate.
Definition: RVec.hxx:795
Namespace for new ROOT classes and functions.
Definition: StringConv.hxx:21
ROOT::Detail::RDF::ColumnNames_t ColumnNames_t
Definition: RDataFrame.cxx:790
RooArgSet S(const RooAbsArg &v1)
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Return the weighted mean of an array a with length n.
Definition: TMath.h:1061
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
Double_t StdDev(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:516
STL namespace.
basic_string_view< char > string_view
Definition: RStringView.hxx:35
make_integer_sequence< size_t, _Np > make_index_sequence
Definition: tree.py:1
Lightweight storage for a collection of types.
Definition: TypeTraits.hxx:27
auto * a
Definition: textangle.C:12