Logo ROOT  
Reference Guide
RLoopManager.cxx
Go to the documentation of this file.
1/*************************************************************************
2 * Copyright (C) 1995-2021, Rene Brun and Fons Rademakers. *
3 * All rights reserved. *
4 * *
5 * For the licensing terms see $ROOTSYS/LICENSE. *
6 * For the list of contributors see $ROOTSYS/README/CREDITS. *
7 *************************************************************************/
8
9#include "RConfigure.h" // R__USE_IMT
10#include "ROOT/RDataSource.hxx"
12#include "ROOT/InternalTreeUtils.hxx" // GetTreeFullPaths
20#include "ROOT/RLogger.hxx"
21#include "RtypesCore.h" // Long64_t
22#include "TStopwatch.h"
23#include "TBranchElement.h"
24#include "TBranchObject.h"
25#include "TChain.h"
26#include "TEntryList.h"
27#include "TFile.h"
28#include "TFriendElement.h"
29#include "TInterpreter.h"
30#include "TROOT.h" // IsImplicitMTEnabled
31#include "TTreeReader.h"
32#include "TTree.h" // For MaxTreeSizeRAII. Revert when #6640 will be solved.
33
34#ifdef R__USE_IMT
37#endif
38
39#include <algorithm>
40#include <atomic>
41#include <cassert>
42#include <exception>
43#include <functional>
44#include <iostream>
45#include <memory>
46#include <stdexcept>
47#include <string>
48#include <sstream>
49#include <thread>
50#include <unordered_map>
51#include <vector>
52#include <set>
53#include <limits> // For MaxTreeSizeRAII. Revert when #6640 will be solved.
54
55using namespace ROOT::Detail::RDF;
56using namespace ROOT::Internal::RDF;
57
58namespace {
59/// A helper function that returns all RDF code that is currently scheduled for just-in-time compilation.
60/// This allows different RLoopManager instances to share these data.
61/// We want RLoopManagers to be able to add their code to a global "code to execute via cling",
62/// so that, lazily, we can jit everything that's needed by all RDFs in one go, which is potentially
63/// much faster than jitting each RLoopManager's code separately.
64static std::string &GetCodeToJit()
65{
66 static std::string code;
67 return code;
68}
69
70static bool ContainsLeaf(const std::set<TLeaf *> &leaves, TLeaf *leaf)
71{
72 return (leaves.find(leaf) != leaves.end());
73}
74
75///////////////////////////////////////////////////////////////////////////////
76/// This overload does not check whether the leaf/branch is already in bNamesReg. In case this is a friend leaf/branch,
77/// `allowDuplicates` controls whether we add both `friendname.bname` and `bname` or just the shorter version.
78static void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
79 const std::string &friendName, bool allowDuplicates)
80{
81 if (!friendName.empty()) {
82 // In case of a friend tree, users might prepend its name/alias to the branch names
83 const auto friendBName = friendName + "." + branchName;
84 if (bNamesReg.insert(friendBName).second)
85 bNames.push_back(friendBName);
86 }
87
88 if (allowDuplicates || friendName.empty()) {
89 if (bNamesReg.insert(branchName).second)
90 bNames.push_back(branchName);
91 }
92}
93
94///////////////////////////////////////////////////////////////////////////////
95/// This overload makes sure that the TLeaf has not been already inserted.
96static void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
97 const std::string &friendName, std::set<TLeaf *> &foundLeaves, TLeaf *leaf,
98 bool allowDuplicates)
99{
100 const bool canAdd = allowDuplicates ? true : !ContainsLeaf(foundLeaves, leaf);
101 if (!canAdd) {
102 return;
103 }
104
105 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
106
107 foundLeaves.insert(leaf);
108}
109
110static void ExploreBranch(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames, TBranch *b,
111 std::string prefix, std::string &friendName, bool allowDuplicates)
112{
113 for (auto sb : *b->GetListOfBranches()) {
114 TBranch *subBranch = static_cast<TBranch *>(sb);
115 auto subBranchName = std::string(subBranch->GetName());
116 auto fullName = prefix + subBranchName;
117
118 std::string newPrefix;
119 if (!prefix.empty())
120 newPrefix = fullName + ".";
121
122 ExploreBranch(t, bNamesReg, bNames, subBranch, newPrefix, friendName, allowDuplicates);
123
124 auto branchDirectlyFromTree = t.GetBranch(fullName.c_str());
125 if (!branchDirectlyFromTree)
126 branchDirectlyFromTree = t.FindBranch(fullName.c_str()); // try harder
127 if (branchDirectlyFromTree)
128 InsertBranchName(bNamesReg, bNames, std::string(branchDirectlyFromTree->GetFullName()), friendName,
129 allowDuplicates);
130
131 if (t.GetBranch(subBranchName.c_str()))
132 InsertBranchName(bNamesReg, bNames, subBranchName, friendName, allowDuplicates);
133 }
134}
135
136static void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
137 std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
138{
139 std::set<TLeaf *> foundLeaves;
140 if (!analysedTrees.insert(&t).second) {
141 return;
142 }
143
144 const auto branches = t.GetListOfBranches();
145 // Getting the branches here triggered the read of the first file of the chain if t is a chain.
146 // We check if a tree has been successfully read, otherwise we throw (see ROOT-9984) to avoid further
147 // operations
148 if (!t.GetTree()) {
149 std::string err("GetBranchNames: error in opening the tree ");
150 err += t.GetName();
151 throw std::runtime_error(err);
152 }
153 if (branches) {
154 for (auto b : *branches) {
155 TBranch *branch = static_cast<TBranch *>(b);
156 const auto branchName = std::string(branch->GetName());
157 if (branch->IsA() == TBranch::Class()) {
158 // Leaf list
159 auto listOfLeaves = branch->GetListOfLeaves();
160 if (listOfLeaves->GetEntriesUnsafe() == 1) {
161 auto leaf = static_cast<TLeaf *>(listOfLeaves->UncheckedAt(0));
162 InsertBranchName(bNamesReg, bNames, branchName, friendName, foundLeaves, leaf, allowDuplicates);
163 }
164
165 for (auto leaf : *listOfLeaves) {
166 auto castLeaf = static_cast<TLeaf *>(leaf);
167 const auto leafName = std::string(leaf->GetName());
168 const auto fullName = branchName + "." + leafName;
169 InsertBranchName(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
170 }
171 } else if (branch->IsA() == TBranchObject::Class()) {
172 // TBranchObject
173 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName, allowDuplicates);
174 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
175 } else {
176 // TBranchElement
177 // Check if there is explicit or implicit dot in the name
178
179 bool dotIsImplied = false;
180 auto be = dynamic_cast<TBranchElement *>(b);
181 if (!be)
182 throw std::runtime_error("GetBranchNames: unsupported branch type");
183 // TClonesArray (3) and STL collection (4)
184 if (be->GetType() == 3 || be->GetType() == 4)
185 dotIsImplied = true;
186
187 if (dotIsImplied || branchName.back() == '.')
188 ExploreBranch(t, bNamesReg, bNames, branch, "", friendName, allowDuplicates);
189 else
190 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName, allowDuplicates);
191
192 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
193 }
194 }
195 }
196
197 // The list of friends needs to be accessed via GetTree()->GetListOfFriends()
198 // (and not via GetListOfFriends() directly), otherwise when `t` is a TChain we
199 // might not recover the list correctly (https://github.com/root-project/root/issues/6741).
200 auto friendTrees = t.GetTree()->GetListOfFriends();
201
202 if (!friendTrees)
203 return;
204
205 for (auto friendTreeObj : *friendTrees) {
206 auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
207
208 std::string frName;
209 auto alias = t.GetFriendAlias(friendTree);
210 if (alias != nullptr)
211 frName = std::string(alias);
212 else
213 frName = std::string(friendTree->GetName());
214
215 GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
216 }
217}
218
219static void ThrowIfNSlotsChanged(unsigned int nSlots)
220{
221 const auto currentSlots = RDFInternal::GetNSlots();
222 if (currentSlots != nSlots) {
223 std::string msg = "RLoopManager::Run: when the RDataFrame was constructed the number of slots required was " +
224 std::to_string(nSlots) + ", but when starting the event loop it was " +
225 std::to_string(currentSlots) + ".";
226 if (currentSlots > nSlots)
227 msg += " Maybe EnableImplicitMT() was called after the RDataFrame was constructed?";
228 else
229 msg += " Maybe DisableImplicitMT() was called after the RDataFrame was constructed?";
230 throw std::runtime_error(msg);
231 }
232}
233
234/**
235\struct MaxTreeSizeRAII
236\brief Scope-bound change of `TTree::fgMaxTreeSize`.
237
238This RAII object stores the current value result of `TTree::GetMaxTreeSize`,
239changes it to maximum at construction time and restores it back at destruction
240time. Needed for issue #6523 and should be reverted when #6640 will be solved.
241*/
242struct MaxTreeSizeRAII {
243 Long64_t fOldMaxTreeSize;
244
245 MaxTreeSizeRAII() : fOldMaxTreeSize(TTree::GetMaxTreeSize())
246 {
247 TTree::SetMaxTreeSize(std::numeric_limits<Long64_t>::max());
248 }
249
250 ~MaxTreeSizeRAII() { TTree::SetMaxTreeSize(fOldMaxTreeSize); }
251};
252
253struct DatasetLogInfo {
254 std::string fDataSet;
255 ULong64_t fRangeStart;
256 ULong64_t fRangeEnd;
257 unsigned int fSlot;
258};
259
260std::string LogRangeProcessing(const DatasetLogInfo &info)
261{
262 std::stringstream msg;
263 msg << "Processing " << info.fDataSet << ": entry range [" << info.fRangeStart << "," << info.fRangeEnd - 1
264 << "], using slot " << info.fSlot << " in thread " << std::this_thread::get_id() << '.';
265 return msg.str();
266}
267
268DatasetLogInfo TreeDatasetLogInfo(const TTreeReader &r, unsigned int slot)
269{
270 const auto tree = r.GetTree();
271 const auto chain = dynamic_cast<TChain *>(tree);
272 std::string what;
273 if (chain) {
274 auto files = chain->GetListOfFiles();
275 std::vector<std::string> treeNames;
276 std::vector<std::string> fileNames;
277 for (TObject *f : *files) {
278 treeNames.emplace_back(f->GetName());
279 fileNames.emplace_back(f->GetTitle());
280 }
281 what = "trees {";
282 for (const auto &t : treeNames) {
283 what += t + ",";
284 }
285 what.back() = '}';
286 what += " in files {";
287 for (const auto &f : fileNames) {
288 what += f + ",";
289 }
290 what.back() = '}';
291 } else {
292 assert(tree != nullptr); // to make clang-tidy happy
293 const auto treeName = tree->GetName();
294 what = std::string("tree \"") + treeName + "\"";
295 const auto file = tree->GetCurrentFile();
296 if (file)
297 what += std::string(" in file \"") + file->GetName() + "\"";
298 }
299 const auto entryRange = r.GetEntriesRange();
300 const ULong64_t end = entryRange.second == -1ll ? tree->GetEntries() : entryRange.second;
301 return {std::move(what), static_cast<ULong64_t>(entryRange.first), end, slot};
302}
303
304} // anonymous namespace
305
306namespace ROOT {
307namespace Detail {
308namespace RDF {
309
310/// A RAII object that calls RLoopManager::CleanUpTask at destruction
313 unsigned int fArg;
315
316 RCallCleanUpTask(RLoopManager &lm, unsigned int arg = 0u, TTreeReader *reader = nullptr)
317 : fLoopManager(lm), fArg(arg), fReader(reader)
318 {
319 }
321};
322
323} // namespace RDF
324} // namespace Detail
325} // namespace ROOT
326
327///////////////////////////////////////////////////////////////////////////////
328/// Get all the branches names, including the ones of the friend trees
330{
331 std::set<std::string> bNamesSet;
332 ColumnNames_t bNames;
333 std::set<TTree *> analysedTrees;
334 std::string emptyFrName = "";
335 GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
336 return bNames;
337}
338
339RLoopManager::RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
340 : fTree(std::shared_ptr<TTree>(tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
341 fNSlots(RDFInternal::GetNSlots()),
342 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles),
343 fNewSampleNotifier(fNSlots), fSampleInfos(fNSlots)
344{
345}
346
348 : fNEmptyEntries(nEmptyEntries), fNSlots(RDFInternal::GetNSlots()),
349 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles), fNewSampleNotifier(fNSlots),
350 fSampleInfos(fNSlots)
351{
352}
353
354RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
355 : fDefaultColumns(defaultBranches), fNSlots(RDFInternal::GetNSlots()),
356 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
357 fDataSource(std::move(ds)), fNewSampleNotifier(fNSlots), fSampleInfos(fNSlots)
358{
359 fDataSource->SetNSlots(fNSlots);
360}
361
362struct RSlotRAII {
364 unsigned int fSlot;
365 RSlotRAII(RSlotStack &slotStack) : fSlotStack(slotStack), fSlot(slotStack.GetSlot()) {}
366 ~RSlotRAII() { fSlotStack.ReturnSlot(fSlot); }
367};
368
369/// Run event loop with no source files, in parallel.
371{
372#ifdef R__USE_IMT
373 RSlotStack slotStack(fNSlots);
374 // Working with an empty tree.
375 // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
376 const auto nEntriesPerSlot = fNEmptyEntries / (fNSlots * 2);
377 auto remainder = fNEmptyEntries % (fNSlots * 2);
378 std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
379 ULong64_t start = 0;
380 while (start < fNEmptyEntries) {
381 ULong64_t end = start + nEntriesPerSlot;
382 if (remainder > 0) {
383 ++end;
384 --remainder;
385 }
386 entryRanges.emplace_back(start, end);
387 start = end;
388 }
389
390 // Each task will generate a subrange of entries
391 auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
392 RSlotRAII slotRAII(slotStack);
393 auto slot = slotRAII.fSlot;
394 RCallCleanUpTask cleanup(*this, slot);
395 InitNodeSlots(nullptr, slot);
396 R__LOG_INFO(RDFLogChannel()) << LogRangeProcessing({"an empty source", range.first, range.second, slot});
397 try {
398 UpdateSampleInfo(slot, range);
399 for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
400 RunAndCheckFilters(slot, currEntry);
401 }
402 } catch (...) {
403 // Error might throw in experiment frameworks like CMSSW
404 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
405 throw;
406 }
407 };
408
410 pool.Foreach(genFunction, entryRanges);
411
412#endif // not implemented otherwise
413}
414
415/// Run event loop with no source files, in sequence.
417{
418 InitNodeSlots(nullptr, 0);
419 R__LOG_INFO(RDFLogChannel()) << LogRangeProcessing({"an empty source", 0, fNEmptyEntries, 0u});
420 RCallCleanUpTask cleanup(*this);
421 try {
422 UpdateSampleInfo(/*slot*/0, {0, fNEmptyEntries});
423 for (ULong64_t currEntry = 0; currEntry < fNEmptyEntries && fNStopsReceived < fNChildren; ++currEntry) {
424 RunAndCheckFilters(0, currEntry);
425 }
426 } catch (...) {
427 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
428 throw;
429 }
430}
431
432/// Run event loop over one or multiple ROOT files, in parallel.
434{
435#ifdef R__USE_IMT
436 RSlotStack slotStack(fNSlots);
437 const auto &entryList = fTree->GetEntryList() ? *fTree->GetEntryList() : TEntryList();
438 auto tp = std::make_unique<ROOT::TTreeProcessorMT>(*fTree, entryList, fNSlots);
439
440 std::atomic<ULong64_t> entryCount(0ull);
441
442 tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
443 RSlotRAII slotRAII(slotStack);
444 auto slot = slotRAII.fSlot;
445 RCallCleanUpTask cleanup(*this, slot, &r);
446 InitNodeSlots(&r, slot);
447 R__LOG_INFO(RDFLogChannel()) << LogRangeProcessing(TreeDatasetLogInfo(r, slot));
448 const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
449 const auto nEntries = entryRange.second - entryRange.first;
450 auto count = entryCount.fetch_add(nEntries);
451 try {
452 // recursive call to check filters and conditionally execute actions
453 while (r.Next()) {
454 if (fNewSampleNotifier.CheckFlag(slot)) {
455 UpdateSampleInfo(slot, r);
456 }
457 RunAndCheckFilters(slot, count++);
458 }
459 } catch (...) {
460 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
461 throw;
462 }
463 // fNStopsReceived < fNChildren is always true at the moment as we don't support event loop early quitting in
464 // multi-thread runs, but it costs nothing to be safe and future-proof in case we add support for that later.
465 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
466 // something went wrong in the TTreeReader event loop
467 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
468 std::to_string(r.GetEntryStatus()));
469 }
470 });
471#endif // no-op otherwise (will not be called)
472}
473
474/// Run event loop over one or multiple ROOT files, in sequence.
476{
477 TTreeReader r(fTree.get(), fTree->GetEntryList());
478 if (0 == fTree->GetEntriesFast())
479 return;
480 RCallCleanUpTask cleanup(*this, 0u, &r);
481 InitNodeSlots(&r, 0);
482 R__LOG_INFO(RDFLogChannel()) << LogRangeProcessing(TreeDatasetLogInfo(r, 0u));
483
484 // recursive call to check filters and conditionally execute actions
485 // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
486 try {
487 while (r.Next() && fNStopsReceived < fNChildren) {
489 UpdateSampleInfo(/*slot*/0, r);
490 }
491 RunAndCheckFilters(0, r.GetCurrentEntry());
492 }
493 } catch (...) {
494 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
495 throw;
496 }
497 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
498 // something went wrong in the TTreeReader event loop
499 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
500 std::to_string(r.GetEntryStatus()));
501 }
502}
503
504/// Run event loop over data accessed through a DataSource, in sequence.
506{
507 assert(fDataSource != nullptr);
508 fDataSource->CallInitialize();
509 auto ranges = fDataSource->GetEntryRanges();
510 while (!ranges.empty() && fNStopsReceived < fNChildren) {
511 InitNodeSlots(nullptr, 0u);
512 fDataSource->InitSlot(0u, 0ull);
513 RCallCleanUpTask cleanup(*this);
514 try {
515 for (const auto &range : ranges) {
516 const auto start = range.first;
517 const auto end = range.second;
518 R__LOG_INFO(RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, 0u});
519 for (auto entry = start; entry < end && fNStopsReceived < fNChildren; ++entry) {
520 if (fDataSource->SetEntry(0u, entry)) {
521 RunAndCheckFilters(0u, entry);
522 }
523 }
524 }
525 } catch (...) {
526 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
527 throw;
528 }
529 fDataSource->CallFinalizeSlot(0u);
530 ranges = fDataSource->GetEntryRanges();
531 }
532 fDataSource->CallFinalize();
533}
534
535/// Run event loop over data accessed through a DataSource, in parallel.
537{
538#ifdef R__USE_IMT
539 assert(fDataSource != nullptr);
540 RSlotStack slotStack(fNSlots);
542
543 // Each task works on a subrange of entries
544 auto runOnRange = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
545 RSlotRAII slotRAII(slotStack);
546 const auto slot = slotRAII.fSlot;
547 InitNodeSlots(nullptr, slot);
548 RCallCleanUpTask cleanup(*this, slot);
549 fDataSource->InitSlot(slot, range.first);
550 const auto start = range.first;
551 const auto end = range.second;
552 R__LOG_INFO(RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, slot});
553 try {
554 for (auto entry = start; entry < end; ++entry) {
555 if (fDataSource->SetEntry(slot, entry)) {
556 RunAndCheckFilters(slot, entry);
557 }
558 }
559 } catch (...) {
560 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
561 throw;
562 }
563 fDataSource->CallFinalizeSlot(slot);
564 };
565
566 fDataSource->CallInitialize();
567 auto ranges = fDataSource->GetEntryRanges();
568 while (!ranges.empty()) {
569 pool.Foreach(runOnRange, ranges);
570 ranges = fDataSource->GetEntryRanges();
571 }
572 fDataSource->CallFinalize();
573#endif // not implemented otherwise (never called)
574}
575
576/// Execute actions and make sure named filters are called for each event.
577/// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
578void RLoopManager::RunAndCheckFilters(unsigned int slot, Long64_t entry)
579{
580 // data-block callbacks run before the rest of the graph
581 if (fNewSampleNotifier.CheckFlag(slot)) {
582 for (auto &callback : fSampleCallbacks) {
583 callback(slot, fSampleInfos[slot]);
584 }
586 }
587
588 for (auto &actionPtr : fBookedActions)
589 actionPtr->Run(slot, entry);
590 for (auto &namedFilterPtr : fBookedNamedFilters)
591 namedFilterPtr->CheckFilters(slot, entry);
592 for (auto &callback : fCallbacks)
593 callback(slot);
594}
595
596/// Build TTreeReaderValues for all nodes
597/// This method loops over all filters, actions and other booked objects and
598/// calls their `InitSlot` method, to get them ready for running a task.
600{
601 SetupSampleCallbacks(r, slot);
602 for (auto &ptr : fBookedActions)
603 ptr->InitSlot(r, slot);
604 for (auto &ptr : fBookedFilters)
605 ptr->InitSlot(r, slot);
606 for (auto &ptr : fBookedDefines)
607 ptr->InitSlot(r, slot);
608 for (auto &ptr : fBookedVariations)
609 ptr->InitSlot(r, slot);
610
611 for (auto &callback : fCallbacksOnce)
612 callback(slot);
613}
614
616 if (r != nullptr) {
617 // we need to set a notifier so that we run the callbacks every time we switch to a new TTree
618 // `PrependLink` inserts this notifier into the TTree/TChain's linked list of notifiers
619 fNewSampleNotifier.GetChainNotifyLink(slot).PrependLink(*r->GetTree());
620 }
621 // Whatever the data source, initially set the "new data block" flag:
622 // - for TChains, this ensures that we don't skip the first data block because
623 // the correct tree is already loaded
624 // - for RDataSources and empty sources, which currently don't have data blocks, this
625 // ensures that we run once per task
627}
628
629void RLoopManager::UpdateSampleInfo(unsigned int slot, const std::pair<ULong64_t, ULong64_t> &range) {
631 "Empty source, range: {" + std::to_string(range.first) + ", " + std::to_string(range.second) + "}", range);
632}
633
635 // one GetTree to retrieve the TChain, another to retrieve the underlying TTree
636 auto *tree = r.GetTree()->GetTree();
637 R__ASSERT(tree != nullptr);
638 const std::string treename = ROOT::Internal::TreeUtils::GetTreeFullPaths(*tree)[0];
639 auto *file = tree->GetCurrentFile();
640 const std::string fname = file != nullptr ? file->GetName() : "#inmemorytree#";
641
642
643 std::pair<Long64_t, Long64_t> range = r.GetEntriesRange();
644 R__ASSERT(range.first >= 0);
645 if (range.second == -1) {
646 range.second = tree->GetEntries(); // convert '-1', i.e. 'until the end', to the actual entry number
647 }
648
649 fSampleInfos[slot] = RSampleInfo(fname + "/" + treename, range);
650}
651
652/// Initialize all nodes of the functional graph before running the event loop.
653/// This method is called once per event-loop and performs generic initialization
654/// operations that do not depend on the specific processing slot (i.e. operations
655/// that are common for all threads).
657{
659 for (auto &filter : fBookedFilters)
660 filter->InitNode();
661 for (auto &range : fBookedRanges)
662 range->InitNode();
663 for (auto &ptr : fBookedActions)
664 ptr->Initialize();
665}
666
667/// Perform clean-up operations. To be called at the end of each event loop.
669{
670 fMustRunNamedFilters = false;
671
672 // forget RActions and detach TResultProxies
673 for (auto &ptr : fBookedActions)
674 ptr->Finalize();
675
676 fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
677 fBookedActions.clear();
678
679 // reset children counts
680 fNChildren = 0;
681 fNStopsReceived = 0;
682 for (auto &ptr : fBookedFilters)
683 ptr->ResetChildrenCount();
684 for (auto &ptr : fBookedRanges)
685 ptr->ResetChildrenCount();
686
687 fCallbacks.clear();
688 fCallbacksOnce.clear();
689 fSampleCallbacks.clear();
690}
691
692/// Perform clean-up operations. To be called at the end of each task execution.
693void RLoopManager::CleanUpTask(TTreeReader *r, unsigned int slot)
694{
695 if (r != nullptr)
696 fNewSampleNotifier.GetChainNotifyLink(slot).RemoveLink(*r->GetTree());
697 for (auto &ptr : fBookedActions)
698 ptr->FinalizeSlot(slot);
699 for (auto &ptr : fBookedFilters)
700 ptr->FinalizeSlot(slot);
701 for (auto &ptr : fBookedDefines)
702 ptr->FinalizeSlot(slot);
703}
704
705/// Add RDF nodes that require just-in-time compilation to the computation graph.
706/// This method also clears the contents of GetCodeToJit().
708{
709 // TODO this should be a read lock unless we find GetCodeToJit non-empty
711
712 const std::string code = std::move(GetCodeToJit());
713 if (code.empty()) {
714 R__LOG_INFO(RDFLogChannel()) << "Nothing to jit and execute.";
715 return;
716 }
717
719 s.Start();
720 RDFInternal::InterpreterCalc(code, "RLoopManager::Run");
721 s.Stop();
722 R__LOG_INFO(RDFLogChannel()) << "Just-in-time compilation phase completed"
723 << (s.RealTime() > 1e-3 ? " in " + std::to_string(s.RealTime()) + " seconds." : ".");
724}
725
726/// Trigger counting of number of children nodes for each node of the functional graph.
727/// This is done once before starting the event loop. Each action sends an `increase children count` signal
728/// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
729/// children counter. Each node only propagates the signal once, even if it receives it multiple times.
730/// Named filters also send an `increase children count` signal, just like actions, as they always execute during
731/// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
733{
734 for (auto &actionPtr : fBookedActions)
735 actionPtr->TriggerChildrenCount();
736 for (auto &namedFilterPtr : fBookedNamedFilters)
737 namedFilterPtr->TriggerChildrenCount();
738}
739
740/// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
741/// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
743{
744 // Change value of TTree::GetMaxTreeSize only for this scope. Revert when #6640 will be solved.
745 MaxTreeSizeRAII ctxtmts;
746
747 R__LOG_INFO(RDFLogChannel()) << "Starting event loop number " << fNRuns << '.';
748
749 ThrowIfNSlotsChanged(GetNSlots());
750
751 Jit();
752
753 InitNodes();
754
756 s.Start();
757 switch (fLoopType) {
761 case ELoopType::kNoFiles: RunEmptySource(); break;
764 }
765 s.Stop();
766
767 CleanUpNodes();
768
769 fNRuns++;
770
771 R__LOG_INFO(RDFLogChannel()) << "Finished event loop number " << fNRuns - 1 << " (" << s.CpuTime() << "s CPU, "
772 << s.RealTime() << "s elapsed).";
773}
774
775/// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
777{
778 return fDefaultColumns;
779}
780
782{
783 return fTree.get();
784}
785
787{
788 fBookedActions.emplace_back(actionPtr);
789}
790
792{
795}
796
798{
799 fBookedFilters.emplace_back(filterPtr);
800 if (filterPtr->HasName()) {
801 fBookedNamedFilters.emplace_back(filterPtr);
803 }
804}
805
807{
810}
811
813{
814 fBookedRanges.emplace_back(rangePtr);
815}
816
818{
820}
821
823{
824 fBookedDefines.emplace_back(ptr);
825}
826
828{
830}
831
833{
834 fBookedVariations.emplace_back(v);
835}
836
838{
840}
841
842// dummy call, end of recursive chain of calls
844{
845 return true;
846}
847
848/// Call `FillReport` on all booked filters
850{
851 for (const auto &fPtr : fBookedNamedFilters)
852 fPtr->FillReport(rep);
853}
854
855void RLoopManager::ToJitExec(const std::string &code) const
856{
858 GetCodeToJit().append(code);
859}
860
861void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
862{
863 if (everyNEvents == 0ull)
864 fCallbacksOnce.emplace_back(std::move(f), fNSlots);
865 else
866 fCallbacks.emplace_back(everyNEvents, std::move(f), fNSlots);
867}
868
869std::vector<std::string> RLoopManager::GetFiltersNames()
870{
871 std::vector<std::string> filters;
872 for (auto &filter : fBookedFilters) {
873 auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
874 filters.push_back(name);
875 }
876 return filters;
877}
878
879std::vector<RNodeBase *> RLoopManager::GetGraphEdges() const
880{
881 std::vector<RNodeBase *> nodes(fBookedFilters.size() + fBookedRanges.size());
882 auto it = std::copy(fBookedFilters.begin(), fBookedFilters.end(), nodes.begin());
883 std::copy(fBookedRanges.begin(), fBookedRanges.end(), it);
884 return nodes;
885}
886
887std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions() const
888{
889 std::vector<RDFInternal::RActionBase *> actions(fBookedActions.size() + fRunActions.size());
890 auto it = std::copy(fBookedActions.begin(), fBookedActions.end(), actions.begin());
891 std::copy(fRunActions.begin(), fRunActions.end(), it);
892 return actions;
893}
894
895std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph(
896 std::unordered_map<void *, std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode>> &visitedMap)
897{
898 // If there is already a node for the RLoopManager return it. If there is not, return a new one.
899 auto duplicateRLoopManagerIt = visitedMap.find((void *)this);
900 if (duplicateRLoopManagerIt != visitedMap.end())
901 return duplicateRLoopManagerIt->second;
902
903 std::string name;
904 if (fDataSource) {
905 name = fDataSource->GetLabel();
906 } else if (fTree) {
907 name = fTree->GetName();
908 } else {
909 name = "Empty source<BR/>Entries: " + std::to_string(fNEmptyEntries);
910 }
911 auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(
913 visitedMap[(void *)this] = thisNode;
914 return thisNode;
915}
916
917////////////////////////////////////////////////////////////////////////////
918/// Return all valid TTree::Branch names (caching results for subsequent calls).
919/// Never use fBranchNames directy, always request it through this method.
921{
922 if (fValidBranchNames.empty() && fTree) {
923 fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
924 }
925 return fValidBranchNames;
926}
927
928bool RLoopManager::HasDSValuePtrs(const std::string &col) const
929{
930 return fDSValuePtrMap.find(col) != fDSValuePtrMap.end();
931}
932
933void RLoopManager::AddDSValuePtrs(const std::string &col, const std::vector<void *> ptrs)
934{
935 fDSValuePtrMap[col] = ptrs;
936}
937
939{
940 if (callback)
941 fSampleCallbacks.emplace_back(std::move(callback));
942}
#define R__LOG_INFO(...)
Definition: RLogger.hxx:364
#define f(i)
Definition: RSha256.hxx:104
#define e(i)
Definition: RSha256.hxx:103
long long Long64_t
Definition: RtypesCore.h:80
unsigned long long ULong64_t
Definition: RtypesCore.h:81
#define R__ASSERT(e)
Definition: TError.h:118
const char * filters[]
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 Float_t Float_t b
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
char name[80]
Definition: TGX11.cxx:110
R__EXTERN TVirtualMutex * gROOTMutex
Definition: TROOT.h:63
#define R__LOCKGUARD(mutex)
The head node of a RDF computation graph.
void UpdateSampleInfo(unsigned int slot, const std::pair< ULong64_t, ULong64_t > &range)
RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
unsigned int fNRuns
Number of event loops run.
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.
void AddDSValuePtrs(const std::string &col, const std::vector< void * > ptrs)
const ColumnNames_t & GetBranchNames()
Return all valid TTree::Branch names (caching results for subsequent calls).
void ToJitExec(const std::string &) const
void AddSampleCallback(ROOT::RDF::SampleCallback_t &&callback)
std::vector< RDFInternal::RActionBase * > GetAllActions() const
Return all actions, either booked or already run.
std::vector< ROOT::RDF::RSampleInfo > fSampleInfos
std::shared_ptr< TTree > fTree
Shared pointer to the input TTree.
std::vector< RDefineBase * > fBookedDefines
void RunTreeReader()
Run event loop over one or multiple ROOT files, in sequence.
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< ROOT::RDF::SampleCallback_t > fSampleCallbacks
Registered callbacks to call at the beginning of each "data block".
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.
std::vector< RFilterBase * > fBookedFilters
std::vector< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
std::vector< RDFInternal::RCallback > fCallbacks
Registered callbacks.
const ELoopType fLoopType
The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
void SetupSampleCallbacks(TTreeReader *r, unsigned int slot)
ColumnNames_t fValidBranchNames
Cache of the tree/chain branch names. Never access directy, always use GetBranchNames().
void CleanUpTask(TTreeReader *r, unsigned int slot)
Perform clean-up operations. To be called at the end of each task execution.
std::map< std::string, std::vector< void * > > fDSValuePtrMap
Registry of per-slot value pointers for booked data-source columns.
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.
unsigned int GetNSlots() const
void RunDataSourceMT()
Run event loop over data accessed through a DataSource, in parallel.
bool HasDSValuePtrs(const std::string &col) const
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.
RDFInternal::RNewSampleNotifier fNewSampleNotifier
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...
std::vector< RDFInternal::ROneTimeCallback > fCallbacksOnce
Registered callbacks to invoke just once before running the loop.
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.
std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > GetGraph(std::unordered_map< void *, std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > > &visitedMap)
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:47
unsigned int fNChildren
Number of nodes of the functional graph hanging from this object.
Definition: RNodeBase.hxx:46
bool CheckFlag(unsigned int slot) const
TNotifyLink< RNewSampleFlag > & GetChainNotifyLink(unsigned int slot)
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 type includes all parts of RVariation that do not depend on the callable signature.
This type represents a sample identifier, to be used in conjunction with RDataFrame features such as ...
Definition: RSampleInfo.hxx:32
This class provides a simple interface to execute the same task multiple times in parallel threads,...
void Foreach(F func, unsigned nTimes, unsigned nChunks=0)
Execute a function without arguments several times in parallel, dividing the execution in nChunks.
A Branch for the case of an object.
static TClass * Class()
A TTree is a list of TBranches.
Definition: TBranch.h:89
static TClass * Class()
TClass * IsA() const override
Definition: TBranch.h:291
TObjArray * GetListOfLeaves()
Definition: TBranch.h:243
A chain is a collection of files containing TTree objects.
Definition: TChain.h:33
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.
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:57
const char * GetName() const override
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
Stopwatch class.
Definition: TStopwatch.h:28
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition: TTreeReader.h:44
@ kEntryBeyondEnd
last entry loop has reached its end
Definition: TTreeReader.h:139
A TTree represents a columnar dataset.
Definition: TTree.h:79
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:4809
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:5256
static void SetMaxTreeSize(Long64_t maxsize=100000000000LL)
Set the maximum size in bytes of a Tree file (static function).
Definition: TTree.cxx:9149
virtual TObjArray * GetListOfBranches()
Definition: TTree.h:484
virtual TTree * GetTree() const
Definition: TTree.h:513
virtual TList * GetListOfFriends() const
Definition: TTree.h:486
virtual const char * GetFriendAlias(TTree *) const
If the 'tree' is a friend, this method returns its alias name.
Definition: TTree.cxx:5999
std::vector< std::string > GetTreeFullPaths(const TTree &tree)
Retrieve the full path(s) to a TTree or the trees in a TChain.
RVec< PromoteTypes< T0, T1 > > remainder(const T0 &x, const RVec< T1 > &v)
Definition: RVec.hxx:1742
ROOT::Experimental::RLogChannel & RDFLogChannel()
Definition: RDFUtils.cxx:37
std::vector< std::string > GetBranchNames(TTree &t, bool allowDuplicates=true)
Get all the branches names, including the ones of the friend trees.
unsigned int GetNSlots()
Definition: RDFUtils.cxx:285
void Erase(const T &that, std::vector< T > &v)
Erase that element from vector v
Definition: Utils.hxx:187
Long64_t InterpreterCalc(const std::string &code, const std::string &context="")
Jit code in the interpreter with TInterpreter::Calc, throw in case of errors.
Definition: RDFUtils.cxx:329
std::vector< std::string > ColumnNames_t
Definition: Utils.hxx:35
std::function< void(unsigned int, const ROOT::RDF::RSampleInfo &)> SampleCallback_t
The type of a data-block callback, registered with a RDataFrame computation graph via e....
Definition: RSampleInfo.hxx:84
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:167
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition: TROOT.cxx:558
static constexpr double s
Definition: file.py:1
Definition: tree.py:1
static const char * what
Definition: stlLoader.cc:6
A RAII object that calls RLoopManager::CleanUpTask at destruction.
RCallCleanUpTask(RLoopManager &lm, unsigned int arg=0u, TTreeReader *reader=nullptr)
RSlotRAII(RSlotStack &slotStack)
RSlotStack & fSlotStack
unsigned int fSlot