Logo ROOT   6.14/05
Reference Guide
RDFInterfaceUtils.hxx
Go to the documentation of this file.
1 // Author: Enrico Guiraud, Danilo Piparo CERN 02/2018
2 
3 /*************************************************************************
4  * Copyright (C) 1995-2016, 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 
15 #include <ROOT/RMakeUnique.hxx>
16 #include <ROOT/RStringView.hxx>
17 #include <ROOT/RDFActionHelpers.hxx> // for BuildAndBook
18 #include <ROOT/RDFNodes.hxx>
19 #include <ROOT/RDFUtils.hxx>
20 #include <ROOT/TypeTraits.hxx>
21 #include <ROOT/TSeq.hxx>
22 #include <algorithm>
23 #include <deque>
24 #include <functional>
25 #include <initializer_list>
26 #include <list>
27 #include <map>
28 #include <memory>
29 #include <string>
30 #include <type_traits>
31 #include <typeinfo>
32 #include <vector>
33 
34 #include "RtypesCore.h"
35 #include "TError.h"
36 #include "TH1.h"
37 
38 class TObjArray;
39 class TTree;
40 namespace ROOT {
41 
42 namespace RDF {
43 class RDataSource;
44 } // namespace RDF
45 
46 } // namespace ROOT
47 
48 /// \cond HIDDEN_SYMBOLS
49 
50 namespace ROOT {
51 namespace Internal {
52 namespace RDF {
53 using namespace ROOT::Detail::RDF;
54 using namespace ROOT::RDF;
55 namespace TTraits = ROOT::TypeTraits;
57 
58 /// An helper object that sets and resets gErrorIgnoreLevel via RAII.
59 class RIgnoreErrorLevelRAII {
60 private:
61  int fCurIgnoreErrorLevel = gErrorIgnoreLevel;
62 
63 public:
64  RIgnoreErrorLevelRAII(int errorIgnoreLevel) { gErrorIgnoreLevel = errorIgnoreLevel; }
65  RIgnoreErrorLevelRAII() { gErrorIgnoreLevel = fCurIgnoreErrorLevel; }
66 };
67 
68 /****** BuildAndBook overloads *******/
69 // BuildAndBook builds a RAction with the right operation and books it with the RLoopManager
70 
71 // clang-format off
72 /// This namespace defines types to be used for tag dispatching in RInterface.
73 namespace ActionTypes {
74 // they cannot just be forward declared: we need concrete types for jitting and to use them with TClass::GetClass
75 struct Histo1D {};
76 struct Histo2D {};
77 struct Histo3D {};
78 struct Profile1D {};
79 struct Profile2D {};
80 struct Min {};
81 struct Max {};
82 struct Sum {};
83 struct Mean {};
84 struct Fill {};
85 }
86 // clang-format on
87 
88 template <int D, typename P, template <int, typename, template <typename> class> class... S>
89 class THist;
90 
91 /// Check whether a histogram type is a classic or v7 histogram.
92 template <typename T>
93 struct IsV7Hist : public std::false_type {
94  static_assert(std::is_base_of<TH1, T>::value, "not implemented for this type");
95 };
96 
97 template <int D, typename P, template <int, typename, template <typename> class> class... S>
98 struct IsV7Hist<THist<D, P, S...>> : public std::true_type {
99 };
100 
101 template <typename T, bool ISV7HISTO = IsV7Hist<T>::value>
102 struct HistoUtils {
103  static void SetCanExtendAllAxes(T &h) { h.SetCanExtend(::TH1::kAllAxes); }
104  static bool HasAxisLimits(T &h)
105  {
106  auto xaxis = h.GetXaxis();
107  return !(xaxis->GetXmin() == 0. && xaxis->GetXmax() == 0.);
108  }
109 };
110 
111 template <typename T>
112 struct HistoUtils<T, true> {
113  static void SetCanExtendAllAxes(T &) {}
114  static bool HasAxisLimits(T &) { return true; }
115 };
116 
117 // Generic filling (covers Histo2D, Histo3D, Profile1D and Profile2D actions, with and without weights)
118 template <typename... BranchTypes, typename ActionType, typename ActionResultType, typename PrevNodeType>
119 RActionBase *BuildAndBook(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &h,
120  const unsigned int nSlots, RLoopManager &loopManager, PrevNodeType &prevNode, ActionType *)
121 {
122  using Helper_t = FillTOHelper<ActionResultType>;
123  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
124  auto action = std::make_shared<Action_t>(Helper_t(h, nSlots), bl, prevNode);
125  loopManager.Book(action);
126  return action.get();
127 }
128 
129 // Histo1D filling (must handle the special case of distinguishing FillTOHelper and FillHelper
130 template <typename... BranchTypes, typename PrevNodeType>
131 RActionBase *BuildAndBook(const ColumnNames_t &bl, const std::shared_ptr<::TH1D> &h, const unsigned int nSlots,
132  RLoopManager &loopManager, PrevNodeType &prevNode, ActionTypes::Histo1D *)
133 {
134  auto hasAxisLimits = HistoUtils<::TH1D>::HasAxisLimits(*h);
135 
136  RActionBase *actionBase;
137  if (hasAxisLimits) {
138  using Helper_t = FillTOHelper<::TH1D>;
139  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
140  auto action = std::make_shared<Action_t>(Helper_t(h, nSlots), bl, prevNode);
141  loopManager.Book(action);
142  actionBase = action.get();
143  } else {
144  using Helper_t = FillHelper;
145  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchTypes...>>;
146  auto action = std::make_shared<Action_t>(Helper_t(h, nSlots), bl, prevNode);
147  loopManager.Book(action);
148  actionBase = action.get();
149  }
150 
151  return actionBase;
152 }
153 
154 // Min action
155 template <typename BranchType, typename PrevNodeType, typename ActionResultType>
156 RActionBase *
157 BuildAndBook(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &minV, const unsigned int nSlots,
158  RLoopManager &loopManager, PrevNodeType &prevNode, ActionTypes::Min *)
159 {
160  using Helper_t = MinHelper<ActionResultType>;
161  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
162  auto action = std::make_shared<Action_t>(Helper_t(minV, nSlots), bl, prevNode);
163  loopManager.Book(action);
164  return action.get();
165 }
166 
167 // Max action
168 template <typename BranchType, typename PrevNodeType, typename ActionResultType>
169 RActionBase *
170 BuildAndBook(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &maxV, const unsigned int nSlots,
171  RLoopManager &loopManager, PrevNodeType &prevNode, ActionTypes::Max *)
172 {
173  using Helper_t = MaxHelper<ActionResultType>;
174  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
175  auto action = std::make_shared<Action_t>(Helper_t(maxV, nSlots), bl, prevNode);
176  loopManager.Book(action);
177  return action.get();
178 }
179 
180 // Sum action
181 template <typename BranchType, typename PrevNodeType, typename ActionResultType>
182 RActionBase *
183 BuildAndBook(const ColumnNames_t &bl, const std::shared_ptr<ActionResultType> &sumV, const unsigned int nSlots,
184  RLoopManager &loopManager, PrevNodeType &prevNode, ActionTypes::Sum *)
185 {
186  using Helper_t = SumHelper<ActionResultType>;
187  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
188  auto action = std::make_shared<Action_t>(Helper_t(sumV, nSlots), bl, prevNode);
189  loopManager.Book(action);
190  return action.get();
191 }
192 
193 // Mean action
194 template <typename BranchType, typename PrevNodeType>
195 RActionBase *BuildAndBook(const ColumnNames_t &bl, const std::shared_ptr<double> &meanV, const unsigned int nSlots,
196  RLoopManager &loopManager, PrevNodeType &prevNode, ActionTypes::Mean *)
197 {
198  using Helper_t = MeanHelper;
199  using Action_t = RAction<Helper_t, PrevNodeType, TTraits::TypeList<BranchType>>;
200  auto action = std::make_shared<Action_t>(Helper_t(meanV, nSlots), bl, prevNode);
201  loopManager.Book(action);
202  return action.get();
203 }
204 /****** end BuildAndBook ******/
205 
206 template <typename Filter>
207 void CheckFilter(Filter &)
208 {
209  using FilterRet_t = typename RDF::CallableTraits<Filter>::ret_type;
210  static_assert(std::is_same<FilterRet_t, bool>::value, "filter functions must return a bool");
211 }
212 
213 void CheckCustomColumn(std::string_view definedCol, TTree *treePtr, const ColumnNames_t &customCols,
214  const ColumnNames_t &dataSourceColumns);
215 
216 using TmpBranchBasePtr_t = std::shared_ptr<RCustomColumnBase>;
217 
218 void BookFilterJit(RJittedFilter *jittedFilter, void *prevNode, std::string_view prevNodeTypeName,
220  const std::map<std::string, std::string> &aliasMap, const ColumnNames_t &branches,
221  const ColumnNames_t &customColumns, TTree *tree, RDataSource *ds, unsigned int namespaceID);
222 
224 
225 std::string JitBuildAndBook(const ColumnNames_t &bl, const std::string &prevNodeTypename, void *prevNode,
226  const std::type_info &art, const std::type_info &at, const void *r, TTree *tree,
227  const unsigned int nSlots, const ColumnNames_t &customColumns, RDataSource *ds,
228  const std::shared_ptr<RActionBase *> *const actionPtrPtr, unsigned int namespaceID);
229 
230 // allocate a shared_ptr on the heap, return a reference to it. the user is responsible of deleting the shared_ptr*.
231 // this function is meant to only be used by RInterface's action methods, and should be deprecated as soon as we find
232 // a better way to make jitting work: the problem it solves is that we need to pass the same shared_ptr to the Helper
233 // object of each action and to the RResultPtr returned by the action. While the former is only instantiated when
234 // the event loop is about to start, the latter has to be returned to the user as soon as the action is booked.
235 // a heap allocated shared_ptr will stay alive long enough that at jitting time its address is still valid.
236 template <typename T>
237 std::shared_ptr<T> *MakeSharedOnHeap(const std::shared_ptr<T> &shPtr)
238 {
239  return new std::shared_ptr<T>(shPtr);
240 }
241 
242 bool AtLeastOneEmptyString(const std::vector<std::string_view> strings);
243 
244 /* The following functions upcast shared ptrs to RFilter, RCustomColumn, RRange to their parent class (***Base).
245  * Shared ptrs to RLoopManager are just copied, as well as shared ptrs to ***Base classes. */
246 std::shared_ptr<RFilterBase> UpcastNode(const std::shared_ptr<RFilterBase> ptr);
247 std::shared_ptr<RCustomColumnBase> UpcastNode(const std::shared_ptr<RCustomColumnBase> ptr);
248 std::shared_ptr<RRangeBase> UpcastNode(const std::shared_ptr<RRangeBase> ptr);
249 std::shared_ptr<RLoopManager> UpcastNode(const std::shared_ptr<RLoopManager> ptr);
250 std::shared_ptr<RJittedFilter> UpcastNode(const std::shared_ptr<RJittedFilter> ptr);
251 
252 ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns,
253  const ColumnNames_t &datasetColumns, const ColumnNames_t &validCustomColumns,
254  RDataSource *ds);
255 
256 std::vector<bool> FindUndefinedDSColumns(const ColumnNames_t &requestedCols, const ColumnNames_t &definedDSCols);
257 
258 /// Helper function to be used by `DefineDataSourceColumns`
259 template <typename T>
260 void DefineDSColumnHelper(std::string_view name, RLoopManager &lm, RDataSource &ds)
261 {
262  auto readers = ds.GetColumnReaders<T>(name);
263  auto getValue = [readers](unsigned int slot) { return *readers[slot]; };
265  lm.Book(std::make_shared<NewCol_t>(name, std::move(getValue), ColumnNames_t{}, &lm, /*isDSColumn=*/true));
266  lm.AddCustomColumnName(name);
267  lm.AddDataSourceColumn(name);
268 }
269 
270 /// Take a list of data-source column names and define the ones that haven't been defined yet.
271 template <typename... ColumnTypes, std::size_t... S>
272 void DefineDataSourceColumns(const std::vector<std::string> &columns, RLoopManager &lm, RDataSource &ds,
273  std::index_sequence<S...>, TTraits::TypeList<ColumnTypes...>)
274 {
275  const auto mustBeDefined = FindUndefinedDSColumns(columns, lm.GetCustomColumnNames());
276  if (std::none_of(mustBeDefined.begin(), mustBeDefined.end(), [](bool b) { return b; })) {
277  // no need to define any column
278  return;
279  } else {
280  // hack to expand a template parameter pack without c++17 fold expressions.
281  std::initializer_list<int> expander{
282  (mustBeDefined[S] ? DefineDSColumnHelper<ColumnTypes>(columns[S], lm, ds) : /*no-op*/ ((void)0), 0)...};
283  (void)expander; // avoid unused variable warnings
284  }
285 }
286 
287 // this function is meant to be called by the jitted code generated by BookFilterJit
288 template <typename F, typename PrevNode>
289 void JitFilterHelper(F &&f, const ColumnNames_t &cols, std::string_view name, RJittedFilter *jittedFilter,
290  PrevNode *prevNode)
291 {
292  // mock Filter logic -- validity checks and Define-ition of RDataSource columns
293  using F_t = RFilter<F, PrevNode>;
294  using ColTypes_t = typename TTraits::CallableTraits<F>::arg_types;
295  constexpr auto nColumns = ColTypes_t::list_size;
296  RDFInternal::CheckFilter(f);
297  auto &lm = *jittedFilter->GetLoopManagerUnchecked(); // RLoopManager must exist at this time
298  auto ds = lm.GetDataSource();
299  if (ds)
300  RDFInternal::DefineDataSourceColumns(cols, lm, *ds, std::make_index_sequence<nColumns>(), ColTypes_t());
301 
302  jittedFilter->SetFilter(std::make_unique<F_t>(std::move(f), cols, *prevNode, name));
303 }
304 
305 template <typename F>
306 void JitDefineHelper(F &&f, const ColumnNames_t &cols, std::string_view name, RLoopManager *lm)
307 {
309  using ColTypes_t = typename TTraits::CallableTraits<F>::arg_types;
310  constexpr auto nColumns = ColTypes_t::list_size;
311 
312  auto ds = lm->GetDataSource();
313  if (ds)
314  RDFInternal::DefineDataSourceColumns(cols, *lm, *ds, std::make_index_sequence<nColumns>(), ColTypes_t());
315 
316  lm->Book(std::make_shared<NewCol_t>(name, std::move(f), cols, lm));
317 }
318 
319 /// Convenience function invoked by jitted code to build action nodes at runtime
320 template <typename ActionType, typename... BranchTypes, typename PrevNodeType, typename ActionResultType>
321 void CallBuildAndBook(PrevNodeType &prevNode, const ColumnNames_t &bl, const unsigned int nSlots,
322  const std::shared_ptr<ActionResultType> *rOnHeap,
323  const std::shared_ptr<RActionBase *> *actionPtrPtrOnHeap)
324 {
325  // if we are here it means we are jitting, if we are jitting the loop manager must be alive
326  auto &loopManager = *prevNode.GetLoopManagerUnchecked();
327  using ColTypes_t = TypeList<BranchTypes...>;
328  constexpr auto nColumns = ColTypes_t::list_size;
329  auto ds = loopManager.GetDataSource();
330  if (ds)
331  DefineDataSourceColumns(bl, loopManager, *ds, std::make_index_sequence<nColumns>(), ColTypes_t());
332  RActionBase *actionPtr =
333  BuildAndBook<BranchTypes...>(bl, *rOnHeap, nSlots, loopManager, prevNode, (ActionType *)nullptr);
334  **actionPtrPtrOnHeap = actionPtr;
335  delete rOnHeap;
336  delete actionPtrPtrOnHeap;
337 }
338 
339 /// The contained `type` alias is `double` if `T == TInferType`, `U` if `T == std::container<U>`, `T` otherwise.
340 template <typename T, bool Container = TTraits::IsContainer<T>::value && !std::is_same<T, std::string>::value>
341 struct TMinReturnType {
342  using type = T;
343 };
344 
345 template <>
346 struct TMinReturnType<TInferType, false> {
347  using type = double;
348 };
349 
350 template <typename T>
351 struct TMinReturnType<T, true> {
353 };
354 
355 // return wrapper around f that prepends an `unsigned int slot` parameter
356 template <typename R, typename F, typename... Args>
357 std::function<R(unsigned int, Args...)> AddSlotParameter(F &f, TypeList<Args...>)
358 {
359  return [f](unsigned int, Args... a) -> R { return f(a...); };
360 }
361 
362 template <typename BranchType, typename... Rest>
363 struct TNeedJitting {
364  static constexpr bool value = TNeedJitting<Rest...>::value;
365 };
366 
367 template <typename... Rest>
368 struct TNeedJitting<TInferType, Rest...> {
369  static constexpr bool value = true;
370 };
371 
372 template <typename T>
373 struct TNeedJitting<T> {
374  static constexpr bool value = false;
375 };
376 
377 template <>
378 struct TNeedJitting<TInferType> {
379  static constexpr bool value = true;
380 };
381 
382 ColumnNames_t GetBranchNames(TTree &t);
383 
384 ColumnNames_t GetTopLevelBranchNames(TTree &t);
385 
386 ///////////////////////////////////////////////////////////////////////////////
387 /// Check preconditions for RInterface::Aggregate:
388 /// - the aggregator callable must have signature `U(U,T)` or `void(U&,T)`.
389 /// - the merge callable must have signature `U(U,U)` or `void(std::vector<U>&)`
391  typename mergeArgsNoDecay_t = typename CallableTraits<Merge>::arg_types_nodecay,
392  typename mergeArgs_t = typename CallableTraits<Merge>::arg_types,
393  typename mergeRet_t = typename CallableTraits<Merge>::ret_type>
394 void CheckAggregate(TypeList<U, T>)
395 {
396  constexpr bool isAggregatorOk =
397  (std::is_same<R, decayedU>::value) || (std::is_same<R, void>::value && std::is_lvalue_reference<U>::value);
398  static_assert(isAggregatorOk, "aggregator function must have signature `U(U,T)` or `void(U&,T)`");
399  constexpr bool isMergeOk =
400  (std::is_same<TypeList<decayedU, decayedU>, mergeArgs_t>::value && std::is_same<decayedU, mergeRet_t>::value) ||
401  (std::is_same<TypeList<std::vector<decayedU> &>, mergeArgsNoDecay_t>::value &&
402  std::is_same<void, mergeRet_t>::value);
403  static_assert(isMergeOk, "merge function must have signature `U(U,U)` or `void(std::vector<U>&)`");
404 }
405 
406 ///////////////////////////////////////////////////////////////////////////////
407 /// This overload of CheckAggregate is called when the aggregator takes more than two arguments
408 template <typename R, typename T>
409 void CheckAggregate(T)
410 {
411  static_assert(sizeof(T) == 0, "aggregator function must take exactly two arguments");
412 }
413 
414 ///////////////////////////////////////////////////////////////////////////////
415 /// Check as many template parameters were passed as the number of column names, throw if this is not the case.
416 void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames);
417 
418 /// Return local BranchNames or default BranchNames according to which one should be used
419 const ColumnNames_t SelectColumns(unsigned int nArgs, const ColumnNames_t &bl, const ColumnNames_t &defBl);
420 
421 /// Check whether column names refer to a valid branch of a TTree or have been `Define`d. Return invalid column names.
422 ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns,
423  const ColumnNames_t &definedCols, const ColumnNames_t &dataSourceColumns);
424 
425 bool IsInternalColumn(std::string_view colName);
426 
427 // Check if a condition is true for all types
428 template <bool...>
429 struct TBoolPack;
430 
431 template <bool... bs>
432 using IsTrueForAllImpl_t = typename std::is_same<TBoolPack<bs..., true>, TBoolPack<true, bs...>>;
433 
434 template <bool... Conditions>
435 struct TEvalAnd {
436  static constexpr bool value = IsTrueForAllImpl_t<Conditions...>::value;
437 };
438 
439 // Check if a class is a specialisation of stl containers templates
440 // clang-format off
441 
442 template <typename>
443 struct IsList_t : std::false_type {};
444 
445 template <typename T>
446 struct IsList_t<std::list<T>> : std::true_type {};
447 
448 template <typename>
449 struct IsDeque_t : std::false_type {};
450 
451 template <typename T>
452 struct IsDeque_t<std::deque<T>> : std::true_type {};
453 // clang-format on
454 
455 } // namespace RDF
456 } // namespace Internal
457 
458 namespace Detail {
459 namespace RDF {
460 
461 /// The aliased type is `double` if `T == TInferType`, `U` if `T == container<U>`, `T` otherwise.
462 template <typename T>
463 using MinReturnType_t = typename RDFInternal::TMinReturnType<T>::type;
464 
465 template <typename T>
466 using MaxReturnType_t = MinReturnType_t<T>;
467 
468 template <typename T>
469 using SumReturnType_t = MinReturnType_t<T>;
470 
471 } // namespace RDF
472 } // namespace Detail
473 } // namespace ROOT
474 
475 /// \endcond
476 
477 #endif
std::string JitBuildAndBook(const ColumnNames_t &bl, const std::string &prevNodeTypename, void *prevNode, const std::type_info &art, const std::type_info &at, const void *rOnHeap, TTree *tree, const unsigned int nSlots, const ColumnNames_t &customColumns, RDataSource *ds, const std::shared_ptr< RActionBase *> *const actionPtrPtr, unsigned int namespaceID)
void BookDefineJit(std::string_view name, std::string_view expression, RLoopManager &lm, RDataSource *ds)
An array of TObjects.
Definition: TObjArray.h:37
void AddCustomColumnName(std::string_view name)
Definition: RDFNodes.hxx:209
ColumnNames_t GetTopLevelBranchNames(TTree &t)
Get all the top-level branches names, including the ones of the friend trees.
bool AtLeastOneEmptyString(const std::vector< std::string_view > strings)
R__EXTERN Int_t gErrorIgnoreLevel
Definition: TError.h:105
Namespace for new ROOT classes and functions.
Definition: StringConv.hxx:21
bool IsInternalColumn(std::string_view colName)
typename TakeFirstParameter< T >::type TakeFirstParameter_t
Definition: TypeTraits.hxx:191
double T(double x)
Definition: ChebyshevPol.h:34
#define f(i)
Definition: RSha256.hxx:104
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:168
STL namespace.
#define R(a, b, c, d, e, f, g, h, i)
Definition: RSha256.hxx:110
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 AddDataSourceColumn(std::string_view name)
Definition: RDFNodes.hxx:207
RLoopManager * GetLoopManagerUnchecked()
Definition: RDFNodes.cxx:523
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:146
A wrapper around a concrete RFilter, which forwards all calls to it RJittedFilter is the type of the ...
Definition: RDFNodes.hxx:611
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:636
ColumnNames_t FindUnknownColumns(const ColumnNames_t &requiredCols, const ColumnNames_t &datasetColumns, const ColumnNames_t &definedCols, const ColumnNames_t &dataSourceColumns)
RLoopManager * GetLoopManagerUnchecked() const
Definition: RDFNodes.cxx:88
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 ...
RooArgSet S(const RooAbsArg &v1)
#define F(x, y, z)
const ColumnNames_t & GetCustomColumnNames() const
Definition: RDFNodes.hxx:185
ROOT::R::TRInterface & r
Definition: Object.C:4
std::shared_ptr< RFilterBase > UpcastNode(const std::shared_ptr< RFilterBase > ptr)
auto * a
Definition: textangle.C:12
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:1098
void BookFilterJit(RJittedFilter *jittedFilter, void *prevNode, std::string_view prevNodeTypeName, std::string_view name, std::string_view expression, const std::map< std::string, std::string > &aliasMap, const ColumnNames_t &branches, const ColumnNames_t &customCols, TTree *tree, RDataSource *ds, unsigned int namespaceID)
static double P[]
ColumnNames_t GetValidatedColumnNames(RLoopManager &lm, const unsigned int nColumns, const ColumnNames_t &columns, const ColumnNames_t &datasetColumns, const ColumnNames_t &validCustomColumns, RDataSource *ds)
Given the desired number of columns and the user-provided list of columns:
#define h(i)
Definition: RSha256.hxx:106
Lightweight storage for a collection of types.
Definition: TypeTraits.hxx:27
int type
Definition: TGX11.cxx:120
basic_string_view< char > string_view
Definition: RStringView.hxx:35
ROOT type_traits extensions.
Definition: TypeTraits.hxx:23
void CheckCustomColumn(std::string_view definedCol, TTree *treePtr, const ColumnNames_t &customCols, const ColumnNames_t &dataSourceColumns)
Extract types from the signature of a callable object. See CallableTraits.
Definition: TypeTraits.hxx:38
T Sum(const RVec< T > &v)
Sum elements.
Definition: RVec.hxx:591
typedef void((*Func_t)())
RDataSource * GetDataSource() const
Definition: RDFNodes.hxx:190
void CheckTypesAndPars(unsigned int nTemplateParams, unsigned int nColumnNames)
void SetFilter(std::unique_ptr< RFilterBase > f)
Definition: RDFNodes.cxx:114
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
std::vector< T ** > GetColumnReaders(std::string_view columnName)
Called at most once per column by RDF.
Definition: tree.py:1
A TTree object has a header with a name and a title.
Definition: TTree.h:70
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
Definition: RDataSource.hxx:91
make_integer_sequence< size_t, _Np > make_index_sequence
char name[80]
Definition: TGX11.cxx:109
ColumnNames_t GetBranchNames(TTree &t)
Get all the branches names, including the ones of the friend trees.
void Book(const ActionBasePtr_t &actionPtr)
Definition: RDFNodes.cxx:544