Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RLoopManager.hxx
Go to the documentation of this file.
1// Author: Enrico Guiraud, Danilo Piparo CERN 03/2017
2
3/*************************************************************************
4 * Copyright (C) 1995-2022, 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_RLOOPMANAGER
12#define ROOT_RLOOPMANAGER
13
14#include "ROOT/RDataSource.hxx"
20#include "ROOT/RDF/Utils.hxx"
21
22#include <functional>
23#include <limits>
24#include <map>
25#include <memory>
26#include <set>
27#include <string>
28#include <string_view>
29#include <unordered_map>
30#include <unordered_set>
31#include <vector>
32#include <any>
33
34// forward declarations
35class TTree;
36class TTreeReader;
37class TDirectory;
38
39namespace ROOT {
40namespace RDF {
41class RCutFlowReport;
42class RDataSource;
43} // ns RDF
44
45namespace Internal {
46class RSlotStack;
47namespace RDF {
48class GraphNode;
49class RActionBase;
50class RVariationBase;
51class RDefinesWithReaders;
52class RVariationsWithReaders;
53class RColumnRegister;
54
55namespace GraphDrawing {
57} // ns GraphDrawing
58
59using Callback_t = std::function<void(unsigned int)>;
60
61class RCallback {
64 std::vector<ULong64_t> fCounters;
65
66public:
69 {
70 }
71
72 void operator()(unsigned int slot)
73 {
74 auto &c = fCounters[slot];
75 ++c;
76 if (c == fEveryN) {
77 c = 0ull;
78 fFun(slot);
79 }
80 }
81};
82
85 std::vector<int> fHasBeenCalled; // std::vector<bool> is thread-unsafe for our purposes (and generally evil)
86
87public:
89
90 void operator()(unsigned int slot)
91 {
92 if (fHasBeenCalled[slot] == 1)
93 return;
94 fFun(slot);
96 }
97};
98
99struct RDSRangeRAII;
100
101} // namespace RDF
102} // namespace Internal
103} // namespace ROOT
104
105namespace ROOT {
106namespace Detail {
107namespace RDF {
109
110class RFilterBase;
111class RRangeBase;
112class RDefineBase;
114
115/// The head node of a RDF computation graph.
116/// This class is responsible of running the event loop.
117class RLoopManager : public RNodeBase {
118 using ColumnNames_t = std::vector<std::string>;
119 enum class ELoopType {
120 kInvalid,
121 kNoFiles,
125 };
126
127 friend struct RCallCleanUpTask;
129
130 /*
131 A wrapped reference to a TTree dataset that can be shared by many computation graphs. Ensures lifetime management.
132 It needs to be destroyed *after* the fDataSource data member, in case it is holding the only Python reference to the
133 input dataset and the data source needs to access it before it is destroyed.
134 */
135 std::any fTTreeLifeline{};
136
137 std::vector<RDFInternal::RActionBase *> fBookedActions; ///< Non-owning pointers to actions to be run
138 std::vector<RDFInternal::RActionBase *> fRunActions; ///< Non-owning pointers to actions already run
139 std::vector<RFilterBase *> fBookedFilters;
140 std::vector<RFilterBase *> fBookedNamedFilters; ///< Contains a subset of fBookedFilters, i.e. only the named filters
141 std::vector<RRangeBase *> fBookedRanges;
142 std::vector<RDefineBase *> fBookedDefines;
143 std::vector<RDFInternal::RVariationBase *> fBookedVariations;
144
146 Long64_t fEndEntry{std::numeric_limits<Long64_t>::max()};
147
148 /// Keys are `fname + "/" + treename` as RSampleInfo::fID; Values are pointers to the corresponding sample
149 std::unordered_map<std::string, ROOT::RDF::Experimental::RSample *> fSampleMap;
150 /// Samples need to survive throughout the whole event loop, hence stored as an attribute
151 std::vector<ROOT::RDF::Experimental::RSample> fSamples;
152
154 /// Range of entries created when no data source is specified.
155 std::pair<ULong64_t, ULong64_t> fEmptyEntryRange{};
156 unsigned int fNSlots{1};
158 /// The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
160 std::unique_ptr<RDataSource> fDataSource{}; ///< Owning pointer to a data-source object. Null if no data-source
161 /// Registered callbacks to be executed every N events.
162 /// The registration happens via the RegisterCallback method.
163 std::vector<RDFInternal::RCallback> fCallbacksEveryNEvents;
164 /// Registered callbacks to invoke just once before running the loop.
165 /// The registration happens via the RegisterCallback method.
166 std::vector<RDFInternal::ROneTimeCallback> fCallbacksOnce;
167 /// Registered callbacks to call at the beginning of each "data block".
168 /// The key is the pointer of the corresponding node in the computation graph (a RDefinePerSample or a RAction).
169 std::unordered_map<void *, ROOT::RDF::SampleCallback_t> fSampleCallbacks;
171 std::vector<ROOT::RDF::RSampleInfo> fSampleInfos;
172 unsigned int fNRuns{0}; ///< Number of event loops run
173
174 /// Readers for TTree/RDataSource columns (one per slot), shared by all nodes in the computation graph.
175 std::vector<std::unordered_map<std::string, std::unique_ptr<RColumnReaderBase>>> fDatasetColumnReaders;
176
177 /// Pointer to a shared slot stack in case this instance runs concurrently with others:
178 std::weak_ptr<ROOT::Internal::RSlotStack> fSlotStack;
179
180 void RunEmptySourceMT();
181 void RunEmptySource();
182 void RunDataSourceMT();
183 void RunDataSource();
184 void RunAndCheckFilters(unsigned int slot, Long64_t entry);
185 void InitNodeSlots(TTreeReader *r, unsigned int slot);
186 void InitNodes();
187 void CleanUpNodes();
188 void CleanUpTask(TTreeReader *r, unsigned int slot);
189 void EvalChildrenCounts();
190 void SetupSampleCallbacks(TTreeReader *r, unsigned int slot);
191 void UpdateSampleInfo(unsigned int slot, const std::pair<ULong64_t, ULong64_t> &range);
192 void UpdateSampleInfo(unsigned int slot, TTreeReader &r);
193 std::shared_ptr<ROOT::Internal::RSlotStack> SlotStack() const;
194
195 // List of branches for which we want to suppress the printed error about
196 // missing branch when switching to a new tree. This is modified by readers,
197 // so must be declared before them in this class.
198 std::set<std::string> fSuppressErrorsForMissingBranches{};
200 std::set<std::pair<std::string_view, std::unique_ptr<ROOT::Internal::RDF::RDefinesWithReaders>>>
202 std::set<std::pair<std::string_view, std::unique_ptr<ROOT::Internal::RDF::RVariationsWithReaders>>>
204
205 // deferred function calls to Jitted functions
207 std::size_t fFunctionId{};
208 std::unique_ptr<ROOT::Internal::RDF::RColumnRegister> fColRegister;
209 std::vector<std::string> fColNames;
210 std::shared_ptr<void> fJittedNode;
211 // Extra arguments to be passed to the jitted function, in one value. Each function will need to unpack it
212 // accordingly.
213 std::shared_ptr<void> fExtraArgs;
214 DeferredJitCall(std::size_t id, std::unique_ptr<ROOT::Internal::RDF::RColumnRegister> cols,
215 const std::vector<std::string> &colNamesArg, std::shared_ptr<void> jittedNode,
216 std::shared_ptr<void> arg);
217
223 };
225 std::hash<std::string> fStringHasher{};
226
227public:
229 RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches);
231 RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches);
233
234 // Rule of five
235
236 RLoopManager(const RLoopManager &) = delete;
240 ~RLoopManager() override;
241
242 void Jit();
243 void RunDeferredCalls();
245 void Run(bool jit = true);
246 const ColumnNames_t &GetDefaultColumnNames() const;
248 RDataSource *GetDataSource() const { return fDataSource.get(); }
259 bool CheckFilters(unsigned int, Long64_t) final;
260 unsigned int GetNSlots() const { return fNSlots; }
261 void Report(ROOT::RDF::RCutFlowReport &rep) const final;
262 /// End of recursive chain of calls, does nothing
266 void ToJitExec(const std::string &) const;
267 void RegisterJitHelperCall(const std::string &funcBody,
268 std::unique_ptr<ROOT::Internal::RDF::RColumnRegister> colRegister,
269 const std::vector<std::string> &colnames, std::shared_ptr<void> jittedNode,
270 std::shared_ptr<void> argument = nullptr);
271 void RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f);
272 unsigned int GetNRuns() const { return fNRuns; }
273 bool HasDataSourceColumnReaders(std::string_view col, const std::type_info &ti) const;
274 void AddDataSourceColumnReaders(std::string_view col, std::vector<std::unique_ptr<RColumnReaderBase>> &&readers,
275 const std::type_info &ti);
276 RColumnReaderBase *GetDatasetColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti) const;
277 RColumnReaderBase *AddDataSourceColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti,
279
280 /// End of recursive chain of calls, does nothing
281 void AddFilterName(std::vector<std::string> &) final {}
282 /// For each booked filter, returns either the name or "Unnamed Filter"
283 std::vector<std::string> GetFiltersNames();
284
285 /// Return all graph edges known to RLoopManager
286 /// This includes Filters and Ranges but not Defines.
287 std::vector<RNodeBase *> GetGraphEdges() const;
288
289 /// Return all actions, either booked or already run
290 std::vector<RDFInternal::RActionBase *> GetAllActions() const;
291
292 std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode>
293 GetGraph(std::unordered_map<void *, std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode>> &visitedMap) final;
294
296
297 void SetEmptyEntryRange(std::pair<ULong64_t, ULong64_t> &&newRange);
300
302 std::set<std::pair<std::string_view, std::unique_ptr<ROOT::Internal::RDF::RDefinesWithReaders>>> &
307 std::set<std::pair<std::string_view, std::unique_ptr<ROOT::Internal::RDF::RVariationsWithReaders>>> &
312
313 void SetTTreeLifeline(std::any lifeline);
314 /// Register a slot stack to be used by this RLoopManager. This allows for sharing RDataFrame helpers safely in the
315 /// context of RunGraphs(). Note that the loop manager only stores a weak_ptr, in between runs.
316 void SetSlotStack(const std::shared_ptr<ROOT::Internal::RSlotStack> &slotStack) { fSlotStack = slotStack; }
317
318 void SetDataSource(std::unique_ptr<ROOT::RDF::RDataSource> dataSource);
319
320 void InsertSuppressErrorsForMissingBranch(const std::string &branchName)
321 {
322 fSuppressErrorsForMissingBranches.insert(branchName);
323 }
324 void EraseSuppressErrorsForMissingBranch(const std::string &branchName)
325 {
326 fSuppressErrorsForMissingBranches.erase(branchName);
327 }
328 const std::set<std::string> &GetSuppressErrorsForMissingBranches() const
329 {
331 }
332
333 /// The task run by every thread on the input entry range, for the generic RDataSource.
334 void DataSourceThreadTask(const std::pair<ULong64_t, ULong64_t> &entryRange, ROOT::Internal::RSlotStack &slotStack,
335 std::atomic<ULong64_t> &entryCount);
336 /// The task run by every thread on an entry range (known by the input TTreeReader), for the TTree data source.
337 void
339};
340
341/// \brief Create an RLoopManager that reads a TChain.
342/// \param[in] datasetName Name of the TChain
343/// \param[in] fileNameGlob File name (or glob) in which the TChain is stored.
344/// \param[in] defaultColumns List of default columns, see
345/// \param[in] checkFile file validator boolean
346/// \ref https://root.cern/doc/master/classROOT_1_1RDataFrame.html#default-branches "Default column lists"
347/// \return the RLoopManager instance.
348std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
349CreateLMFromTTree(std::string_view datasetName, std::string_view fileNameGlob,
350 const std::vector<std::string> &defaultColumns, bool checkFile = true);
351
352/// \brief Create an RLoopManager that reads a TChain.
353/// \param[in] datasetName Name of the TChain
354/// \param[in] fileNameGlobs List of file names (potentially globs).
355/// \param[in] defaultColumns List of default columns, see
356/// \param[in] checkFile file validator boolean
357/// \ref https://root.cern/doc/master/classROOT_1_1RDataFrame.html#default-branches "Default column lists"
358/// \return the RLoopManager instance.
359std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
360CreateLMFromTTree(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
361 const std::vector<std::string> &defaultColumns, bool checkFile = true);
362
363/// \brief Create an RLoopManager that reads an RNTuple.
364/// \param[in] datasetName Name of the RNTuple
365/// \param[in] fileNameGlob File name (or glob) in which the RNTuple is stored.
366/// \param[in] defaultColumns List of default columns, see
367/// \ref https://root.cern/doc/master/classROOT_1_1RDataFrame.html#default-branches "Default column lists"
368/// \return the RLoopManager instance.
369std::shared_ptr<ROOT::Detail::RDF::RLoopManager> CreateLMFromRNTuple(std::string_view datasetName,
370 std::string_view fileNameGlob,
371 const std::vector<std::string> &defaultColumns);
372
373/// \brief Create an RLoopManager that reads multiple RNTuples chained vertically.
374/// \param[in] datasetName Name of the RNTuple
375/// \param[in] fileNameGlobs List of file names (potentially globs).
376/// \param[in] defaultColumns List of default columns, see
377/// \ref https://root.cern/doc/master/classROOT_1_1RDataFrame.html#default-branches "Default column lists"
378/// \return the RLoopManager instance.
379std::shared_ptr<ROOT::Detail::RDF::RLoopManager> CreateLMFromRNTuple(std::string_view datasetName,
380 const std::vector<std::string> &fileNameGlobs,
381 const std::vector<std::string> &defaultColumns);
382
383/// \brief Create an RLoopManager opening a file and checking the data format of the dataset.
384/// \param[in] datasetName Name of the dataset in the file.
385/// \param[in] fileNameGlob File name (or glob) in which the dataset is stored.
386/// \param[in] defaultColumns List of default columns, see
387/// \ref https://root.cern/doc/master/classROOT_1_1RDataFrame.html#default-branches "Default column lists"
388/// \throws std::invalid_argument if the file could not be opened.
389/// \return an RLoopManager of the appropriate data source.
390std::shared_ptr<ROOT::Detail::RDF::RLoopManager> CreateLMFromFile(std::string_view datasetName,
391 std::string_view fileNameGlob,
392 const std::vector<std::string> &defaultColumns);
393
394/// \brief Create an RLoopManager that reads many files. The first is opened to infer the data source type.
395/// \param[in] datasetName Name of the dataset.
396/// \param[in] fileNameGlobs List of file names (potentially globs).
397/// \param[in] defaultColumns List of default columns, see
398/// \ref https://root.cern/doc/master/classROOT_1_1RDataFrame.html#default-branches "Default column lists"
399/// \throws std::invalid_argument if the file could not be opened.
400/// \return an RLoopManager of the appropriate data source.
401std::shared_ptr<ROOT::Detail::RDF::RLoopManager> CreateLMFromFile(std::string_view datasetName,
402 const std::vector<std::string> &fileNameGlobs,
403 const std::vector<std::string> &defaultColumns);
404
405} // namespace RDF
406} // namespace Detail
407} // namespace ROOT
408
409#endif
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
unsigned long long ULong64_t
Portable unsigned long integer 8 bytes.
Definition RtypesCore.h:84
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 r
The head node of a RDF computation graph.
RColumnReaderBase * AddDataSourceColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti, TTreeReader *treeReader)
void UpdateSampleInfo(unsigned int slot, const std::pair< ULong64_t, ULong64_t > &range)
unsigned int fNRuns
Number of event loops run.
bool CheckFilters(unsigned int, Long64_t) final
void RegisterJitHelperCall(const std::string &funcBody, std::unique_ptr< ROOT::Internal::RDF::RColumnRegister > colRegister, const std::vector< std::string > &colnames, std::shared_ptr< void > jittedNode, std::shared_ptr< void > argument=nullptr)
void EvalChildrenCounts()
Trigger counting of number of children nodes for each node of the functional graph.
void CleanUpNodes()
Perform clean-up operations. To be called at the end of each event loop.
void RunEmptySource()
Run event loop with no source files, in sequence.
const std::set< std::string > & GetSuppressErrorsForMissingBranches() const
void SetEmptyEntryRange(std::pair< ULong64_t, ULong64_t > &&newRange)
void Report(ROOT::RDF::RCutFlowReport &rep) const final
Call FillReport on all booked filters.
void AddSampleCallback(void *nodePtr, ROOT::RDF::SampleCallback_t &&callback)
std::vector< RFilterBase * > fBookedNamedFilters
Contains a subset of fBookedFilters, i.e. only the named filters.
void RunEmptySourceMT()
Run event loop with no source files, in parallel.
RLoopManager & operator=(RLoopManager &&)=delete
std::hash< std::string > fStringHasher
std::unordered_map< std::string, ROOT::RDF::Experimental::RSample * > fSampleMap
Keys are fname + "/" + treename as RSampleInfo::fID; Values are pointers to the corresponding sample.
void AddDataSourceColumnReaders(std::string_view col, std::vector< std::unique_ptr< RColumnReaderBase > > &&readers, const std::type_info &ti)
std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > GetGraph(std::unordered_map< void *, std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > > &visitedMap) final
RLoopManager & operator=(const RLoopManager &)=delete
void ToJitExec(const std::string &) const
std::vector< RDFInternal::RActionBase * > GetAllActions() const
Return all actions, either booked or already run.
std::vector< ROOT::RDF::RSampleInfo > fSampleInfos
std::set< std::string > fSuppressErrorsForMissingBranches
void ChangeSpec(ROOT::RDF::Experimental::RDatasetSpec &&spec)
Changes the internal TTree held by the RLoopManager.
std::weak_ptr< ROOT::Internal::RSlotStack > fSlotStack
Pointer to a shared slot stack in case this instance runs concurrently with others:
std::set< std::pair< std::string_view, std::unique_ptr< ROOT::Internal::RDF::RVariationsWithReaders > > > & GetUniqueVariationsWithReaders()
std::vector< RDefineBase * > fBookedDefines
void TTreeThreadTask(TTreeReader &treeReader, ROOT::Internal::RSlotStack &slotStack, std::atomic< ULong64_t > &entryCount)
The task run by every thread on an entry range (known by the input TTreeReader), for the TTree data s...
void InsertSuppressErrorsForMissingBranch(const std::string &branchName)
std::vector< RDFInternal::RActionBase * > fRunActions
Non-owning pointers to actions already run.
RLoopManager(const ColumnNames_t &defaultColumns={})
void AddFilterName(std::vector< std::string > &) final
End of recursive chain of calls, does nothing.
std::vector< RRangeBase * > fBookedRanges
std::vector< ROOT::RDF::Experimental::RSample > fSamples
Samples need to survive throughout the whole event loop, hence stored as an attribute.
std::vector< std::string > ColumnNames_t
void RunAndCheckFilters(unsigned int slot, Long64_t entry)
Execute actions and make sure named filters are called for each event.
void ChangeBeginAndEndEntries(Long64_t begin, Long64_t end)
std::vector< RFilterBase * > fBookedFilters
void Run(bool jit=true)
Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
std::set< std::pair< std::string_view, std::unique_ptr< ROOT::Internal::RDF::RVariationsWithReaders > > > fUniqueVariationsWithReaders
std::unordered_map< void *, ROOT::RDF::SampleCallback_t > fSampleCallbacks
Registered callbacks to call at the beginning of each "data block".
std::vector< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
RLoopManager(RLoopManager &&)=delete
void SetupSampleCallbacks(TTreeReader *r, unsigned int slot)
void CleanUpTask(TTreeReader *r, unsigned int slot)
Perform clean-up operations. To be called at the end of each task execution.
std::vector< RDFInternal::RCallback > fCallbacksEveryNEvents
Registered callbacks to be executed every N events.
std::vector< std::unordered_map< std::string, std::unique_ptr< RColumnReaderBase > > > fDatasetColumnReaders
Readers for TTree/RDataSource columns (one per slot), shared by all nodes in the computation graph.
void Register(RDFInternal::RActionBase *actionPtr)
std::vector< DeferredJitCall > fJitHelperCalls
const ColumnNames_t & GetDefaultColumnNames() const
Return the list of default columns – empty if none was provided when constructing the RDataFrame.
std::vector< RDFInternal::RVariationBase * > fBookedVariations
std::vector< RNodeBase * > GetGraphEdges() const
Return all graph edges known to RLoopManager This includes Filters and Ranges but not Defines.
std::set< std::pair< std::string_view, std::unique_ptr< ROOT::Internal::RDF::RDefinesWithReaders > > > & GetUniqueDefinesWithReaders()
void EraseSuppressErrorsForMissingBranch(const std::string &branchName)
RDataSource * GetDataSource() const
void RunDataSourceMT()
Run event loop over data accessed through a DataSource, in parallel.
void PartialReport(ROOT::RDF::RCutFlowReport &) const final
End of recursive chain of calls, does nothing.
void SetSlotStack(const std::shared_ptr< ROOT::Internal::RSlotStack > &slotStack)
Register a slot stack to be used by this RLoopManager.
std::vector< std::string > GetFiltersNames()
For each booked filter, returns either the name or "Unnamed Filter".
RLoopManager(const RLoopManager &)=delete
RDFInternal::RNewSampleNotifier fNewSampleNotifier
std::pair< ULong64_t, ULong64_t > fEmptyEntryRange
Range of entries created when no data source is specified.
std::set< std::pair< std::string_view, std::unique_ptr< ROOT::Internal::RDF::RDefinesWithReaders > > > fUniqueDefinesWithReaders
ROOT::Internal::RDF::RStringCache fCachedColNames
ROOT::Internal::RDF::RStringCache & GetColumnNamesCache()
std::unique_ptr< RDataSource > fDataSource
Owning pointer to a data-source object.
void DataSourceThreadTask(const std::pair< ULong64_t, ULong64_t > &entryRange, ROOT::Internal::RSlotStack &slotStack, std::atomic< ULong64_t > &entryCount)
The task run by every thread on the input entry range, for the generic RDataSource.
void InitNodeSlots(TTreeReader *r, unsigned int slot)
Build TTreeReaderValues for all nodes This method loops over all filters, actions and other booked ob...
std::vector< RDFInternal::ROneTimeCallback > fCallbacksOnce
Registered callbacks to invoke just once before running the loop.
void SetDataSource(std::unique_ptr< ROOT::RDF::RDataSource > dataSource)
void RegisterCallback(ULong64_t everyNEvents, std::function< void(unsigned int)> &&f)
void SetTTreeLifeline(std::any lifeline)
void RunDataSource()
Run event loop over data accessed through a DataSource, in sequence.
void Jit()
Add RDF nodes that require just-in-time compilation to the computation graph.
RColumnReaderBase * GetDatasetColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti) const
std::shared_ptr< ROOT::Internal::RSlotStack > SlotStack() const
Create a slot stack with the desired number of slots or reuse a shared instance.
void Deregister(RDFInternal::RActionBase *actionPtr)
ELoopType fLoopType
The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
void InitNodes()
Initialize all nodes of the functional graph before running the event loop.
bool HasDataSourceColumnReaders(std::string_view col, const std::type_info &ti) const
Return true if AddDataSourceColumnReaders was called for column name col.
RLoopManager * GetLoopManagerUnchecked() final
Base class for non-leaf nodes of the computational graph.
Definition RNodeBase.hxx:43
unsigned int fNStopsReceived
Number of times that a children node signaled to stop processing entries.
Definition RNodeBase.hxx:47
unsigned int fNChildren
Number of nodes of the functional graph hanging from this object.
Definition RNodeBase.hxx:46
Helper class that provides the operation graph nodes.
RCallback(ULong64_t everyN, Callback_t &&f, unsigned int nSlots)
std::vector< ULong64_t > fCounters
void operator()(unsigned int slot)
ROneTimeCallback(Callback_t &&f, unsigned int nSlots)
A Thread-safe cache for strings.
Definition Utils.hxx:290
This type includes all parts of RVariation that do not depend on the callable signature.
A thread-safe list of N indexes (0 to size - 1).
The dataset specification for RDataFrame.
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
Describe directory structure in memory.
Definition TDirectory.h:45
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition TTreeReader.h:46
A TTree represents a columnar dataset.
Definition TTree.h:89
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > CreateLMFromTTree(std::string_view datasetName, std::string_view fileNameGlob, const std::vector< std::string > &defaultColumns, bool checkFile=true)
Create an RLoopManager that reads a TChain.
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > CreateLMFromFile(std::string_view datasetName, std::string_view fileNameGlob, const std::vector< std::string > &defaultColumns)
Create an RLoopManager opening a file and checking the data format of the dataset.
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > CreateLMFromRNTuple(std::string_view datasetName, std::string_view fileNameGlob, const std::vector< std::string > &defaultColumns)
Create an RLoopManager that reads an RNTuple.
std::function< void(unsigned int)> Callback_t
std::function< void(unsigned int, const ROOT::RDF::RSampleInfo &)> SampleCallback_t
The type of a data-block callback, registered with an RDataFrame computation graph via e....
A RAII object that calls RLoopManager::CleanUpTask at destruction.
std::unique_ptr< ROOT::Internal::RDF::RColumnRegister > fColRegister
DeferredJitCall(const DeferredJitCall &)=delete
DeferredJitCall(DeferredJitCall &&) noexcept
DeferredJitCall(std::size_t id, std::unique_ptr< ROOT::Internal::RDF::RColumnRegister > cols, const std::vector< std::string > &colNamesArg, std::shared_ptr< void > jittedNode, std::shared_ptr< void > arg)
DeferredJitCall & operator=(const DeferredJitCall &)=delete