Logo ROOT  
Reference Guide
RLoopManager.cxx
Go to the documentation of this file.
1#include "RConfigure.h" // R__USE_IMT
9#include "RtypesCore.h" // Long64_t
10#include "TBranchElement.h"
11#include "TBranchObject.h"
12#include "TEntryList.h"
13#include "TError.h"
14#include "TInterpreter.h"
15#include "TROOT.h" // IsImplicitMTEnabled
16#include "TTreeReader.h"
17
18#ifdef R__USE_IMT
20#endif
21
22#include <atomic>
23#include <functional>
24#include <memory>
25#include <stdexcept>
26#include <string>
27#include <vector>
28
29using namespace ROOT::Detail::RDF;
30using namespace ROOT::Internal::RDF;
31
32bool ContainsLeaf(const std::set<TLeaf *> &leaves, TLeaf *leaf)
33{
34 return (leaves.find(leaf) != leaves.end());
35}
36
37///////////////////////////////////////////////////////////////////////////////
38/// This overload does not perform any check on the duplicates.
39/// It is used for TBranch objects.
40void UpdateList(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
41 const std::string &friendName)
42{
43
44 if (!friendName.empty()) {
45 // In case of a friend tree, users might prepend its name/alias to the branch names
46 const auto friendBName = friendName + "." + branchName;
47 if (bNamesReg.insert(friendBName).second)
48 bNames.push_back(friendBName);
49 }
50
51 if (bNamesReg.insert(branchName).second)
52 bNames.push_back(branchName);
53}
54
55///////////////////////////////////////////////////////////////////////////////
56/// This overloads makes sure that the TLeaf has not been already inserted.
57void UpdateList(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
58 const std::string &friendName, std::set<TLeaf *> &foundLeaves, TLeaf *leaf, bool allowDuplicates)
59{
60 const bool canAdd = allowDuplicates ? true : !ContainsLeaf(foundLeaves, leaf);
61 if (!canAdd) {
62 return;
63 }
64
65 UpdateList(bNamesReg, bNames, branchName, friendName);
66
67 foundLeaves.insert(leaf);
68}
69
70void ExploreBranch(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames, TBranch *b, std::string prefix,
71 std::string &friendName)
72{
73 for (auto sb : *b->GetListOfBranches()) {
74 TBranch *subBranch = static_cast<TBranch *>(sb);
75 auto subBranchName = std::string(subBranch->GetName());
76 auto fullName = prefix + subBranchName;
77
78 std::string newPrefix;
79 if (!prefix.empty())
80 newPrefix = fullName + ".";
81
82 ExploreBranch(t, bNamesReg, bNames, subBranch, newPrefix, friendName);
83
84 if (t.GetBranch(fullName.c_str()) || t.FindBranch(fullName.c_str()))
85 UpdateList(bNamesReg, bNames, fullName, friendName);
86
87 if (t.GetBranch(subBranchName.c_str()))
88 UpdateList(bNamesReg, bNames, subBranchName, friendName);
89 }
90}
91
92void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
93 std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
94{
95 std::set<TLeaf *> foundLeaves;
96 if (!analysedTrees.insert(&t).second) {
97 return;
98 }
99
100 const auto branches = t.GetListOfBranches();
101 // Getting the branches here triggered the read of the first file of the chain if t is a chain.
102 // We check if a tree has been successfully read, otherwise we throw (see ROOT-9984) to avoid further
103 // operations
104 if (!t.GetTree()) {
105 std::string err("GetBranchNames: error in opening the tree ");
106 err += t.GetName();
107 throw std::runtime_error(err);
108 }
109 if (branches) {
110 for (auto b : *branches) {
111 TBranch *branch = static_cast<TBranch *>(b);
112 const auto branchName = std::string(branch->GetName());
113 if (branch->IsA() == TBranch::Class()) {
114 // Leaf list
115 auto listOfLeaves = branch->GetListOfLeaves();
116 if (listOfLeaves->GetEntries() == 1) {
117 auto leaf = static_cast<TLeaf *>(listOfLeaves->At(0));
118 const auto leafName = std::string(leaf->GetName());
119 if (leafName == branchName) {
120 UpdateList(bNamesReg, bNames, branchName, friendName, foundLeaves, leaf, allowDuplicates);
121 }
122 }
123
124 for (auto leaf : *listOfLeaves) {
125 auto castLeaf = static_cast<TLeaf *>(leaf);
126 const auto leafName = std::string(leaf->GetName());
127 const auto fullName = branchName + "." + leafName;
128 UpdateList(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
129 }
130 } else if (branch->IsA() == TBranchObject::Class()) {
131 // TBranchObject
132 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName);
133 UpdateList(bNamesReg, bNames, branchName, friendName);
134 } else {
135 // TBranchElement
136 // Check if there is explicit or implicit dot in the name
137
138 bool dotIsImplied = false;
139 auto be = dynamic_cast<TBranchElement *>(b);
140 if (!be)
141 throw std::runtime_error("GetBranchNames: unsupported branch type");
142 // TClonesArray (3) and STL collection (4)
143 if (be->GetType() == 3 || be->GetType() == 4)
144 dotIsImplied = true;
145
146 if (dotIsImplied || branchName.back() == '.')
147 ExploreBranch(t, bNamesReg, bNames, branch, "", friendName);
148 else
149 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName);
150
151 UpdateList(bNamesReg, bNames, branchName, friendName);
152 }
153 }
154 }
155
156 auto friendTrees = t.GetListOfFriends();
157
158 if (!friendTrees)
159 return;
160
161 for (auto friendTreeObj : *friendTrees) {
162 auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
163
164 std::string frName;
165 auto alias = t.GetFriendAlias(friendTree);
166 if (alias != nullptr)
167 frName = std::string(alias);
168 else
169 frName = std::string(friendTree->GetName());
170
171 GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
172 }
173}
174
175///////////////////////////////////////////////////////////////////////////////
176/// Get all the branches names, including the ones of the friend trees
178{
179 std::set<std::string> bNamesSet;
180 ColumnNames_t bNames;
181 std::set<TTree *> analysedTrees;
182 std::string emptyFrName = "";
183 GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
184 return bNames;
185}
186
187
188RLoopManager::RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
189 : fTree(std::shared_ptr<TTree>(tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
190 fNSlots(RDFInternal::GetNSlots()),
191 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles)
192{
193}
194
196 : fNEmptyEntries(nEmptyEntries), fNSlots(RDFInternal::GetNSlots()),
197 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles)
198{
199}
200
201RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
202 : fDefaultColumns(defaultBranches), fNSlots(RDFInternal::GetNSlots()),
203 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
204 fDataSource(std::move(ds))
205{
206 fDataSource->SetNSlots(fNSlots);
207}
208
209// ROOT-9559: we cannot handle indexed friends
211{
212 auto friends = fTree->GetListOfFriends();
213 if (!friends)
214 return;
215 for (auto friendElObj : *friends) {
216 auto friendEl = static_cast<TFriendElement *>(friendElObj);
217 auto friendTree = friendEl->GetTree();
218 if (friendTree && friendTree->GetTreeIndex()) {
219 std::string err = fTree->GetName();
220 err += " has a friend, \"";
221 err += friendTree->GetName();
222 err += "\", which has an index. This is not supported.";
223 throw std::runtime_error(err);
224 }
225 }
226}
227
228/// Run event loop with no source files, in parallel.
230{
231#ifdef R__USE_IMT
232 RSlotStack slotStack(fNSlots);
233 // Working with an empty tree.
234 // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
235 const auto nEntriesPerSlot = fNEmptyEntries / (fNSlots * 2);
236 auto remainder = fNEmptyEntries % (fNSlots * 2);
237 std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
238 ULong64_t start = 0;
239 while (start < fNEmptyEntries) {
240 ULong64_t end = start + nEntriesPerSlot;
241 if (remainder > 0) {
242 ++end;
243 --remainder;
244 }
245 entryRanges.emplace_back(start, end);
246 start = end;
247 }
248
249 // Each task will generate a subrange of entries
250 auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
251 auto slot = slotStack.GetSlot();
252 InitNodeSlots(nullptr, slot);
253 for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
254 RunAndCheckFilters(slot, currEntry);
255 }
256 CleanUpTask(slot);
257 slotStack.ReturnSlot(slot);
258 };
259
261 pool.Foreach(genFunction, entryRanges);
262
263#endif // not implemented otherwise
264}
265
266/// Run event loop with no source files, in sequence.
268{
269 InitNodeSlots(nullptr, 0);
270 for (ULong64_t currEntry = 0; currEntry < fNEmptyEntries && fNStopsReceived < fNChildren; ++currEntry) {
271 RunAndCheckFilters(0, currEntry);
272 }
273 CleanUpTask(0u);
274}
275
276/// Run event loop over one or multiple ROOT files, in parallel.
278{
279#ifdef R__USE_IMT
281 RSlotStack slotStack(fNSlots);
282 const auto &entryList = fTree->GetEntryList() ? *fTree->GetEntryList() : TEntryList();
283 auto tp = std::make_unique<ROOT::TTreeProcessorMT>(*fTree, entryList);
284
285 std::atomic<ULong64_t> entryCount(0ull);
286
287 tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
288 auto slot = slotStack.GetSlot();
289 InitNodeSlots(&r, slot);
290 const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
291 const auto nEntries = entryRange.second - entryRange.first;
292 auto count = entryCount.fetch_add(nEntries);
293 // recursive call to check filters and conditionally execute actions
294 while (r.Next()) {
295 RunAndCheckFilters(slot, count++);
296 }
297 CleanUpTask(slot);
298 slotStack.ReturnSlot(slot);
299 });
300#endif // no-op otherwise (will not be called)
301}
302
303/// Run event loop over one or multiple ROOT files, in sequence.
305{
307 TTreeReader r(fTree.get(), fTree->GetEntryList());
308 if (0 == fTree->GetEntriesFast())
309 return;
310 InitNodeSlots(&r, 0);
311
312 // recursive call to check filters and conditionally execute actions
313 // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
314 while (r.Next() && fNStopsReceived < fNChildren) {
315 RunAndCheckFilters(0, r.GetCurrentEntry());
316 }
317 CleanUpTask(0u);
318}
319
320/// Run event loop over data accessed through a DataSource, in sequence.
322{
323 R__ASSERT(fDataSource != nullptr);
324 fDataSource->Initialise();
325 auto ranges = fDataSource->GetEntryRanges();
326 while (!ranges.empty()) {
327 InitNodeSlots(nullptr, 0u);
328 fDataSource->InitSlot(0u, 0ull);
329 for (const auto &range : ranges) {
330 auto end = range.second;
331 for (auto entry = range.first; entry < end; ++entry) {
332 if (fDataSource->SetEntry(0u, entry)) {
333 RunAndCheckFilters(0u, entry);
334 }
335 }
336 }
337 CleanUpTask(0u);
338 fDataSource->FinaliseSlot(0u);
339 ranges = fDataSource->GetEntryRanges();
340 }
341 fDataSource->Finalise();
342}
343
344/// Run event loop over data accessed through a DataSource, in parallel.
346{
347#ifdef R__USE_IMT
348 R__ASSERT(fDataSource != nullptr);
349 RSlotStack slotStack(fNSlots);
351
352 // Each task works on a subrange of entries
353 auto runOnRange = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
354 const auto slot = slotStack.GetSlot();
355 InitNodeSlots(nullptr, slot);
356 fDataSource->InitSlot(slot, range.first);
357 const auto end = range.second;
358 for (auto entry = range.first; entry < end; ++entry) {
359 if (fDataSource->SetEntry(slot, entry)) {
360 RunAndCheckFilters(slot, entry);
361 }
362 }
363 CleanUpTask(slot);
364 fDataSource->FinaliseSlot(slot);
365 slotStack.ReturnSlot(slot);
366 };
367
368 fDataSource->Initialise();
369 auto ranges = fDataSource->GetEntryRanges();
370 while (!ranges.empty()) {
371 pool.Foreach(runOnRange, ranges);
372 ranges = fDataSource->GetEntryRanges();
373 }
374 fDataSource->Finalise();
375#endif // not implemented otherwise (never called)
376}
377
378/// Execute actions and make sure named filters are called for each event.
379/// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
380void RLoopManager::RunAndCheckFilters(unsigned int slot, Long64_t entry)
381{
382 for (auto &actionPtr : fBookedActions)
383 actionPtr->Run(slot, entry);
384 for (auto &namedFilterPtr : fBookedNamedFilters)
385 namedFilterPtr->CheckFilters(slot, entry);
386 for (auto &callback : fCallbacks)
387 callback(slot);
388}
389
390/// Build TTreeReaderValues for all nodes
391/// This method loops over all filters, actions and other booked objects and
392/// calls their `InitRDFValues` methods. It is called once per node per slot, before
393/// running the event loop. It also informs each node of the TTreeReader that
394/// a particular slot will be using.
396{
397 for (auto &ptr : fBookedActions)
398 ptr->InitSlot(r, slot);
399 for (auto &ptr : fBookedFilters)
400 ptr->InitSlot(r, slot);
401 for (auto &callback : fCallbacksOnce)
402 callback(slot);
403}
404
405/// Initialize all nodes of the functional graph before running the event loop.
406/// This method is called once per event-loop and performs generic initialization
407/// operations that do not depend on the specific processing slot (i.e. operations
408/// that are common for all threads).
410{
412 for (auto column : fCustomColumns)
413 column->InitNode();
414 for (auto &filter : fBookedFilters)
415 filter->InitNode();
416 for (auto &range : fBookedRanges)
417 range->InitNode();
418 for (auto &ptr : fBookedActions)
419 ptr->Initialize();
420}
421
422/// Perform clean-up operations. To be called at the end of each event loop.
424{
425 fMustRunNamedFilters = false;
426
427 // forget RActions and detach TResultProxies
428 for (auto &ptr : fBookedActions)
429 ptr->Finalize();
430
431 fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
432 fBookedActions.clear();
433
434 // reset children counts
435 fNChildren = 0;
436 fNStopsReceived = 0;
437 for (auto &ptr : fBookedFilters)
438 ptr->ResetChildrenCount();
439 for (auto &ptr : fBookedRanges)
440 ptr->ResetChildrenCount();
441
442 fCallbacks.clear();
443 fCallbacksOnce.clear();
444}
445
446/// Perform clean-up operations. To be called at the end of each task execution.
447void RLoopManager::CleanUpTask(unsigned int slot)
448{
449 for (auto &ptr : fBookedActions)
450 ptr->FinalizeSlot(slot);
451 for (auto &ptr : fBookedFilters)
452 ptr->ClearTask(slot);
453}
454
455/// Declare to the interpreter type aliases and other entities required by RDF jitted nodes.
456/// This method clears the `fToJitDeclare` member variable.
458{
459 if (fToJitDeclare.empty())
460 return;
461
463 fToJitDeclare.clear();
464}
465
466/// Add RDF nodes that require just-in-time compilation to the computation graph.
467/// This method also invokes JitDeclarations() if needed, and clears the `fToJitExec` member variable.
469{
470 if (fToJitExec.empty())
471 return;
472
474 RDFInternal::InterpreterCalc(fToJitExec, "RLoopManager::Run");
475 fToJitExec.clear();
476}
477
478/// Trigger counting of number of children nodes for each node of the functional graph.
479/// This is done once before starting the event loop. Each action sends an `increase children count` signal
480/// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
481/// children counter. Each node only propagates the signal once, even if it receives it multiple times.
482/// Named filters also send an `increase children count` signal, just like actions, as they always execute during
483/// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
485{
486 for (auto &actionPtr : fBookedActions)
487 actionPtr->TriggerChildrenCount();
488 for (auto &namedFilterPtr : fBookedNamedFilters)
489 namedFilterPtr->TriggerChildrenCount();
490}
491
493{
494 static unsigned int id = 0;
495 ++id;
496 return id;
497}
498
499/// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
500/// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
502{
503 Jit();
504
505 InitNodes();
506
507 switch (fLoopType) {
511 case ELoopType::kNoFiles: RunEmptySource(); break;
514 }
515
516 CleanUpNodes();
517}
518
519/// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
521{
522 return fDefaultColumns;
523}
524
526{
527 return fTree.get();
528}
529
531{
532 fBookedActions.emplace_back(actionPtr);
533}
534
536{
537 RDFInternal::Erase(actionPtr, fRunActions);
538 RDFInternal::Erase(actionPtr, fBookedActions);
539}
540
542{
543 fBookedFilters.emplace_back(filterPtr);
544 if (filterPtr->HasName()) {
545 fBookedNamedFilters.emplace_back(filterPtr);
547 }
548}
549
551{
552 RDFInternal::Erase(filterPtr, fBookedFilters);
553 RDFInternal::Erase(filterPtr, fBookedNamedFilters);
554}
555
557{
558 fBookedRanges.emplace_back(rangePtr);
559}
560
562{
563 RDFInternal::Erase(rangePtr, fBookedRanges);
564}
565
566// dummy call, end of recursive chain of calls
568{
569 return true;
570}
571
572/// Call `FillReport` on all booked filters
574{
575 for (const auto &fPtr : fBookedNamedFilters)
576 fPtr->FillReport(rep);
577}
578
579void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
580{
581 if (everyNEvents == 0ull)
582 fCallbacksOnce.emplace_back(std::move(f), fNSlots);
583 else
584 fCallbacks.emplace_back(everyNEvents, std::move(f), fNSlots);
585}
586
587std::vector<std::string> RLoopManager::GetFiltersNames()
588{
589 std::vector<std::string> filters;
590 for (auto &filter : fBookedFilters) {
591 auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
592 filters.push_back(name);
593 }
594 return filters;
595}
596
597std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions()
598{
599 std::vector<RDFInternal::RActionBase *> actions;
600 actions.insert(actions.begin(), fBookedActions.begin(), fBookedActions.end());
601 actions.insert(actions.begin(), fRunActions.begin(), fRunActions.end());
602 return actions;
603}
604
605std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph()
606{
607 std::string name;
608 if (fDataSource) {
609 name = fDataSource->GetLabel();
610 } else if (fTree) {
611 name = fTree->GetName();
612 } else {
613 name = std::to_string(fNEmptyEntries);
614 }
615
616 auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(name);
617 thisNode->SetRoot();
618 thisNode->SetCounter(0);
619 return thisNode;
620}
621
622////////////////////////////////////////////////////////////////////////////
623/// Return all valid TTree::Branch names (caching results for subsequent calls).
624/// Never use fBranchNames directy, always request it through this method.
626{
627 if (fValidBranchNames.empty() && fTree) {
628 fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
629 }
630 return fValidBranchNames;
631}
void Class()
Definition: Class.C:29
ROOT::R::TRInterface & r
Definition: Object.C:4
void GetBranchNamesImpl(TTree &t, std::set< std::string > &bNamesReg, ColumnNames_t &bNames, std::set< TTree * > &analysedTrees, std::string &friendName, bool allowDuplicates)
void UpdateList(std::set< std::string > &bNamesReg, ColumnNames_t &bNames, const std::string &branchName, const std::string &friendName)
This overload does not perform any check on the duplicates.
void ExploreBranch(TTree &t, std::set< std::string > &bNamesReg, ColumnNames_t &bNames, TBranch *b, std::string prefix, std::string &friendName)
bool ContainsLeaf(const std::set< TLeaf * > &leaves, TLeaf *leaf)
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
long long Long64_t
Definition: RtypesCore.h:69
unsigned long long ULong64_t
Definition: RtypesCore.h:70
#define R__ASSERT(e)
Definition: TError.h:96
const char * filters[]
XFontStruct * id
Definition: TGX11.cxx:108
char name[80]
Definition: TGX11.cxx:109
RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
bool CheckFilters(unsigned int, Long64_t) final
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.
void Report(ROOT::RDF::RCutFlowReport &rep) const final
Call FillReport on all booked filters.
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.
const ColumnNames_t & GetBranchNames()
Return all valid TTree::Branch names (caching results for subsequent calls).
std::shared_ptr< TTree > fTree
Shared pointer to the input TTree.
std::string fToJitExec
Code that should be just-in-time executed right before the event loop.
void RunTreeReader()
Run event loop over one or multiple ROOT files, in sequence.
std::vector< RDFInternal::RActionBase * > GetAllActions()
For all the actions, either booked or run.
void CleanUpTask(unsigned int slot)
Perform clean-up operations. To be called at the end of each task execution.
std::vector< RDFInternal::RActionBase * > fRunActions
Non-owning pointers to actions already run.
void Run()
Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
std::vector< RRangeBase * > fBookedRanges
std::vector< TCallback > fCallbacks
Registered callbacks.
void RunAndCheckFilters(unsigned int slot, Long64_t entry)
Execute actions and make sure named filters are called for each event.
std::vector< RFilterBase * > fBookedFilters
std::vector< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
void JitDeclarations()
Declare to the interpreter type aliases and other entities required by RDF jitted nodes.
std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > GetGraph()
const ELoopType fLoopType
The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
ColumnNames_t fValidBranchNames
Cache of the tree/chain branch names. Never access directy, always use GetBranchNames().
std::vector< RCustomColumnBase * > fCustomColumns
Non-owning container of all custom columns created so far.
const ColumnNames_t & GetDefaultColumnNames() const
Return the list of default columns – empty if none was provided when constructing the RDataFrame.
std::vector< TOneTimeCallback > fCallbacksOnce
Registered callbacks to invoke just once before running the loop.
void RunDataSourceMT()
Run event loop over data accessed through a DataSource, in parallel.
std::string fToJitDeclare
Code that should be just-in-time declared right before the event loop.
std::vector< std::string > GetFiltersNames()
For each booked filter, returns either the name or "Unnamed Filter".
const std::unique_ptr< RDataSource > fDataSource
Owning pointer to a data-source object. Null if no data-source.
static unsigned int GetNextID()
const ColumnNames_t fDefaultColumns
void Book(RDFInternal::RActionBase *actionPtr)
void InitNodeSlots(TTreeReader *r, unsigned int slot)
Build TTreeReaderValues for all nodes This method loops over all filters, actions and other booked ob...
void RegisterCallback(ULong64_t everyNEvents, std::function< void(unsigned int)> &&f)
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.
void RunTreeProcessorMT()
Run event loop over one or multiple ROOT files, in parallel.
void Deregister(RDFInternal::RActionBase *actionPtr)
void InitNodes()
Initialize all nodes of the functional graph before running the event loop.
unsigned int fNStopsReceived
Number of times that a children node signaled to stop processing entries.
Definition: RNodeBase.hxx:45
unsigned int fNChildren
Number of nodes of the functional graph hanging from this object.
Definition: RNodeBase.hxx:44
This is an helper class to allow to pick a slot resorting to a map indexed by thread ids.
Definition: RSlotStack.hxx:26
void ReturnSlot(unsigned int slotNumber)
Definition: RSlotStack.cxx:23
This class provides a simple interface to execute the same task multiple times in parallel,...
void Foreach(F func, unsigned nTimes, unsigned nChunks=0)
Execute func (with no arguments) nTimes in parallel.
A Branch for the case of an object.
A TTree is a list of TBranches.
Definition: TBranch.h:91
TObjArray * GetListOfLeaves()
Definition: TBranch.h:245
A List of entry numbers in a TTree or TChain.
Definition: TEntryList.h:26
A TFriendElement TF describes a TTree object TF in a file.
virtual TTree * GetTree()
Return pointer to friend TTree.
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:49
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition: TTreeReader.h:43
A TTree represents a columnar dataset.
Definition: TTree.h:72
virtual TBranch * FindBranch(const char *name)
Return the branch that correspond to the path 'branchname', which can include the name of the tree or...
Definition: TTree.cxx:4723
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:5170
virtual TObjArray * GetListOfBranches()
Definition: TTree.h:475
virtual TTree * GetTree() const
Definition: TTree.h:504
virtual TList * GetListOfFriends() const
Definition: TTree.h:477
virtual const char * GetFriendAlias(TTree *) const
If the 'tree' is a friend, this method returns its alias name.
Definition: TTree.cxx:5892
ColumnNames_t GetBranchNames(TTree &t, bool allowDuplicates=true)
Get all the branches names, including the ones of the friend trees.
unsigned int GetNSlots()
Definition: RDFUtils.cxx:270
Long64_t InterpreterCalc(const std::string &code, const std::string &context)
Definition: RDFUtils.cxx:310
void InterpreterDeclare(const std::string &code)
Definition: RDFUtils.cxx:302
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
VSD Structures.
Definition: StringConv.hxx:21
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:611
ROOT::Detail::RDF::ColumnNames_t ColumnNames_t
Definition: RDataFrame.cxx:788
Definition: tree.py:1