Logo ROOT   6.16/01
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, std::string &branchName,
41 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 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, std::string &branchName,
58 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())) {
85 UpdateList(bNamesReg, bNames, fullName, friendName);
86
87 } else if (t.GetBranch(subBranchName.c_str())) {
88 UpdateList(bNamesReg, bNames, subBranchName, friendName);
89 }
90 }
91}
92
93void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
94 std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
95{
96 std::set<TLeaf *> foundLeaves;
97 if (!analysedTrees.insert(&t).second) {
98 return;
99 }
100
101 const auto branches = t.GetListOfBranches();
102 if (branches) {
103 std::string prefix = "";
104 for (auto b : *branches) {
105 TBranch *branch = static_cast<TBranch *>(b);
106 auto branchName = std::string(branch->GetName());
107 if (branch->IsA() == TBranch::Class()) {
108 // Leaf list
109 auto listOfLeaves = branch->GetListOfLeaves();
110 if (listOfLeaves->GetEntries() == 1) {
111 auto leaf = static_cast<TLeaf *>(listOfLeaves->At(0));
112 auto leafName = std::string(leaf->GetName());
113 if (leafName == branchName) {
114 UpdateList(bNamesReg, bNames, branchName, friendName, foundLeaves, leaf, allowDuplicates);
115 }
116 }
117
118 for (auto leaf : *listOfLeaves) {
119 auto castLeaf = static_cast<TLeaf *>(leaf);
120 auto leafName = std::string(leaf->GetName());
121 auto fullName = branchName + "." + leafName;
122 UpdateList(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
123 }
124 } else if (branch->IsA() == TBranchObject::Class()) {
125 // TBranchObject
126 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName);
127 UpdateList(bNamesReg, bNames, branchName, friendName);
128 } else {
129 // TBranchElement
130 // Check if there is explicit or implicit dot in the name
131
132 bool dotIsImplied = false;
133 auto be = dynamic_cast<TBranchElement *>(b);
134 if (!be)
135 throw std::runtime_error("GetBranchNames: unsupported branch type");
136 // TClonesArray (3) and STL collection (4)
137 if (be->GetType() == 3 || be->GetType() == 4)
138 dotIsImplied = true;
139
140 if (dotIsImplied || branchName.back() == '.')
141 ExploreBranch(t, bNamesReg, bNames, branch, "", friendName);
142 else
143 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName);
144
145 UpdateList(bNamesReg, bNames, branchName, friendName);
146 }
147 }
148 }
149
150 auto friendTrees = t.GetListOfFriends();
151
152 if (!friendTrees)
153 return;
154
155 for (auto friendTreeObj : *friendTrees) {
156 auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
157
158 std::string frName;
159 auto alias = t.GetFriendAlias(friendTree);
160 if (alias != nullptr)
161 frName = std::string(alias);
162 else
163 frName = std::string(friendTree->GetName());
164
165 GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
166 }
167}
168
169///////////////////////////////////////////////////////////////////////////////
170/// Get all the branches names, including the ones of the friend trees
172{
173 std::set<std::string> bNamesSet;
174 ColumnNames_t bNames;
175 std::set<TTree *> analysedTrees;
176 std::string emptyFrName = "";
177 GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
178 return bNames;
179}
180
181
182RLoopManager::RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
183 : fTree(std::shared_ptr<TTree>(tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
184 fNSlots(RDFInternal::GetNSlots()),
185 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles)
186{
187}
188
190 : fNEmptyEntries(nEmptyEntries), fNSlots(RDFInternal::GetNSlots()),
191 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles)
192{
193}
194
195RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
196 : fDefaultColumns(defaultBranches), fNSlots(RDFInternal::GetNSlots()),
197 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
198 fDataSource(std::move(ds))
199{
200 fDataSource->SetNSlots(fNSlots);
201}
202
203/// Run event loop with no source files, in parallel.
205{
206#ifdef R__USE_IMT
207 RSlotStack slotStack(fNSlots);
208 // Working with an empty tree.
209 // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
210 const auto nEntriesPerSlot = fNEmptyEntries / (fNSlots * 2);
211 auto remainder = fNEmptyEntries % (fNSlots * 2);
212 std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
213 ULong64_t start = 0;
214 while (start < fNEmptyEntries) {
215 ULong64_t end = start + nEntriesPerSlot;
216 if (remainder > 0) {
217 ++end;
218 --remainder;
219 }
220 entryRanges.emplace_back(start, end);
221 start = end;
222 }
223
224 // Each task will generate a subrange of entries
225 auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
226 auto slot = slotStack.GetSlot();
227 InitNodeSlots(nullptr, slot);
228 for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
229 RunAndCheckFilters(slot, currEntry);
230 }
231 CleanUpTask(slot);
232 slotStack.ReturnSlot(slot);
233 };
234
236 pool.Foreach(genFunction, entryRanges);
237
238#endif // not implemented otherwise
239}
240
241/// Run event loop with no source files, in sequence.
243{
244 InitNodeSlots(nullptr, 0);
245 for (ULong64_t currEntry = 0; currEntry < fNEmptyEntries && fNStopsReceived < fNChildren; ++currEntry) {
246 RunAndCheckFilters(0, currEntry);
247 }
248}
249
250/// Run event loop over one or multiple ROOT files, in parallel.
252{
253#ifdef R__USE_IMT
254 RSlotStack slotStack(fNSlots);
255 auto tp = std::make_unique<ROOT::TTreeProcessorMT>(*fTree);
256
257 std::atomic<ULong64_t> entryCount(0ull);
258
259 tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
260 auto slot = slotStack.GetSlot();
261 InitNodeSlots(&r, slot);
262 const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
263 const auto nEntries = entryRange.second - entryRange.first;
264 auto count = entryCount.fetch_add(nEntries);
265 // recursive call to check filters and conditionally execute actions
266 while (r.Next()) {
267 RunAndCheckFilters(slot, count++);
268 }
269 CleanUpTask(slot);
270 slotStack.ReturnSlot(slot);
271 });
272#endif // no-op otherwise (will not be called)
273}
274
275/// Run event loop over one or multiple ROOT files, in sequence.
277{
278 TTreeReader r(fTree.get());
279 if (0 == fTree->GetEntriesFast())
280 return;
281 InitNodeSlots(&r, 0);
282
283 // recursive call to check filters and conditionally execute actions
284 // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
285 while (r.Next() && fNStopsReceived < fNChildren) {
286 RunAndCheckFilters(0, r.GetCurrentEntry());
287 }
288}
289
290/// Run event loop over data accessed through a DataSource, in sequence.
292{
293 R__ASSERT(fDataSource != nullptr);
294 fDataSource->Initialise();
295 auto ranges = fDataSource->GetEntryRanges();
296 while (!ranges.empty()) {
297 InitNodeSlots(nullptr, 0u);
298 fDataSource->InitSlot(0u, 0ull);
299 for (const auto &range : ranges) {
300 auto end = range.second;
301 for (auto entry = range.first; entry < end; ++entry) {
302 if (fDataSource->SetEntry(0u, entry)) {
303 RunAndCheckFilters(0u, entry);
304 }
305 }
306 }
307 fDataSource->FinaliseSlot(0u);
308 ranges = fDataSource->GetEntryRanges();
309 }
310 fDataSource->Finalise();
311}
312
313/// Run event loop over data accessed through a DataSource, in parallel.
315{
316#ifdef R__USE_IMT
317 R__ASSERT(fDataSource != nullptr);
318 RSlotStack slotStack(fNSlots);
320
321 // Each task works on a subrange of entries
322 auto runOnRange = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
323 const auto slot = slotStack.GetSlot();
324 InitNodeSlots(nullptr, slot);
325 fDataSource->InitSlot(slot, range.first);
326 const auto end = range.second;
327 for (auto entry = range.first; entry < end; ++entry) {
328 if (fDataSource->SetEntry(slot, entry)) {
329 RunAndCheckFilters(slot, entry);
330 }
331 }
332 CleanUpTask(slot);
333 fDataSource->FinaliseSlot(slot);
334 slotStack.ReturnSlot(slot);
335 };
336
337 fDataSource->Initialise();
338 auto ranges = fDataSource->GetEntryRanges();
339 while (!ranges.empty()) {
340 pool.Foreach(runOnRange, ranges);
341 ranges = fDataSource->GetEntryRanges();
342 }
343 fDataSource->Finalise();
344#endif // not implemented otherwise (never called)
345}
346
347/// Execute actions and make sure named filters are called for each event.
348/// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
349void RLoopManager::RunAndCheckFilters(unsigned int slot, Long64_t entry)
350{
351 for (auto &actionPtr : fBookedActions)
352 actionPtr->Run(slot, entry);
353 for (auto &namedFilterPtr : fBookedNamedFilters)
354 namedFilterPtr->CheckFilters(slot, entry);
355 for (auto &callback : fCallbacks)
356 callback(slot);
357}
358
359/// Build TTreeReaderValues for all nodes
360/// This method loops over all filters, actions and other booked objects and
361/// calls their `InitRDFValues` methods. It is called once per node per slot, before
362/// running the event loop. It also informs each node of the TTreeReader that
363/// a particular slot will be using.
365{
366 for (auto &ptr : fBookedActions)
367 ptr->InitSlot(r, slot);
368 for (auto &ptr : fBookedFilters)
369 ptr->InitSlot(r, slot);
370 for (auto &callback : fCallbacksOnce)
371 callback(slot);
372}
373
374/// Initialize all nodes of the functional graph before running the event loop.
375/// This method is called once per event-loop and performs generic initialization
376/// operations that do not depend on the specific processing slot (i.e. operations
377/// that are common for all threads).
379{
381 for (auto column : fCustomColumns)
382 column->InitNode();
383 for (auto &filter : fBookedFilters)
384 filter->InitNode();
385 for (auto &range : fBookedRanges)
386 range->InitNode();
387 for (auto &ptr : fBookedActions)
388 ptr->Initialize();
389}
390
391/// Perform clean-up operations. To be called at the end of each event loop.
393{
394 fMustRunNamedFilters = false;
395
396 // forget RActions and detach TResultProxies
397 for (auto &ptr : fBookedActions)
398 ptr->Finalize();
399
400 fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
401 fBookedActions.clear();
402
403 // reset children counts
404 fNChildren = 0;
405 fNStopsReceived = 0;
406 for (auto &ptr : fBookedFilters)
407 ptr->ResetChildrenCount();
408 for (auto &ptr : fBookedRanges)
409 ptr->ResetChildrenCount();
410
411 fCallbacks.clear();
412 fCallbacksOnce.clear();
413}
414
415/// Perform clean-up operations. To be called at the end of each task execution.
416void RLoopManager::CleanUpTask(unsigned int slot)
417{
418 for (auto &ptr : fBookedActions)
419 ptr->FinalizeSlot(slot);
420 for (auto &ptr : fBookedFilters)
421 ptr->ClearTask(slot);
422}
423
424/// Jit all actions that required runtime column type inference, and clean the `fToJit` member variable.
426{
427 auto error = TInterpreter::EErrorCode::kNoError;
428 gInterpreter->Calc(fToJit.c_str(), &error);
429 if (TInterpreter::EErrorCode::kNoError != error) {
430 std::string exceptionText =
431 "An error occurred while jitting. The lines above might indicate the cause of the crash\n";
432 throw std::runtime_error(exceptionText.c_str());
433 }
434 fToJit.clear();
435}
436
437/// Trigger counting of number of children nodes for each node of the functional graph.
438/// This is done once before starting the event loop. Each action sends an `increase children count` signal
439/// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
440/// children counter. Each node only propagates the signal once, even if it receives it multiple times.
441/// Named filters also send an `increase children count` signal, just like actions, as they always execute during
442/// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
444{
445 for (auto &actionPtr : fBookedActions)
446 actionPtr->TriggerChildrenCount();
447 for (auto &namedFilterPtr : fBookedNamedFilters)
448 namedFilterPtr->TriggerChildrenCount();
449}
450
452{
453 static unsigned int id = 0;
454 ++id;
455 return id;
456}
457
458/// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
459/// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
461{
462 if (!fToJit.empty())
464
465 InitNodes();
466
467 switch (fLoopType) {
471 case ELoopType::kNoFiles: RunEmptySource(); break;
474 }
475
476 CleanUpNodes();
477}
478
479/// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
481{
482 return fDefaultColumns;
483}
484
486{
487 return fTree.get();
488}
489
491{
492 fBookedActions.emplace_back(actionPtr);
493}
494
496{
497 RDFInternal::Erase(actionPtr, fRunActions);
498 RDFInternal::Erase(actionPtr, fBookedActions);
499}
500
502{
503 fBookedFilters.emplace_back(filterPtr);
504 if (filterPtr->HasName()) {
505 fBookedNamedFilters.emplace_back(filterPtr);
507 }
508}
509
511{
512 RDFInternal::Erase(filterPtr, fBookedFilters);
513 RDFInternal::Erase(filterPtr, fBookedNamedFilters);
514}
515
517{
518 fBookedRanges.emplace_back(rangePtr);
519}
520
522{
523 RDFInternal::Erase(rangePtr, fBookedRanges);
524}
525
526// dummy call, end of recursive chain of calls
528{
529 return true;
530}
531
532/// Call `FillReport` on all booked filters
534{
535 for (const auto &fPtr : fBookedNamedFilters)
536 fPtr->FillReport(rep);
537}
538
539void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
540{
541 if (everyNEvents == 0ull)
542 fCallbacksOnce.emplace_back(std::move(f), fNSlots);
543 else
544 fCallbacks.emplace_back(everyNEvents, std::move(f), fNSlots);
545}
546
547std::vector<std::string> RLoopManager::GetFiltersNames()
548{
549 std::vector<std::string> filters;
550 for (auto &filter : fBookedFilters) {
551 auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
552 filters.push_back(name);
553 }
554 return filters;
555}
556
557std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions()
558{
559 std::vector<RDFInternal::RActionBase *> actions;
560 actions.insert(actions.begin(), fBookedActions.begin(), fBookedActions.end());
561 actions.insert(actions.begin(), fRunActions.begin(), fRunActions.end());
562 return actions;
563}
564
565std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph()
566{
567 std::string name;
568 if (fDataSource) {
569 name = fDataSource->GetLabel();
570 } else if (fTree) {
571 name = fTree->GetName();
572 } else {
573 name = std::to_string(fNEmptyEntries);
574 }
575
576 auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(name);
577 thisNode->SetRoot();
578 thisNode->SetCounter(0);
579 return thisNode;
580}
581
582////////////////////////////////////////////////////////////////////////////
583/// Return all valid TTree::Branch names (caching results for subsequent calls).
584/// Never use fBranchNames directy, always request it through this method.
586{
587 if (fValidBranchNames.empty() && fTree) {
588 fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
589 }
590 return fValidBranchNames;
591}
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, std::string &branchName, 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[]
#define gInterpreter
Definition: TInterpreter.h:538
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.
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.
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::string fToJit
code that should be jitted and executed right before the event loop
void BuildJittedNodes()
Jit all actions that required runtime column type inference, and clean the fToJit member variable.
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::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 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:20
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:64
TObjArray * GetListOfLeaves()
Definition: TBranch.h:204
A TFriendElement TF describes a TTree object TF in a file.
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:32
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
A simple, robust and fast interface to read values from ROOT colmnar datasets such as TTree,...
Definition: TTreeReader.h:44
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:246
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
Namespace for new ROOT classes and functions.
Definition: StringConv.hxx:21
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:607
ROOT::Detail::RDF::ColumnNames_t ColumnNames_t
Definition: RDataFrame.cxx:790
STL namespace.
Definition: tree.py:1