Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
19#include "ROOT/RLogger.hxx"
20#include "RtypesCore.h" // Long64_t
21#include "TStopwatch.h"
22#include "TBranchElement.h"
23#include "TBranchObject.h"
24#include "TChain.h"
25#include "TEntryList.h"
26#include "TFile.h"
27#include "TFriendElement.h"
28#include "TROOT.h" // IsImplicitMTEnabled
29#include "TTreeReader.h"
30#include "TTree.h" // For MaxTreeSizeRAII. Revert when #6640 will be solved.
31
32#ifdef R__USE_IMT
35#include "ROOT/RSlotStack.hxx"
36#endif
37
38#ifdef R__HAS_ROOT7
39#include "ROOT/RNTuple.hxx"
40#include "ROOT/RNTupleDS.hxx"
41#endif
42
43#include <algorithm>
44#include <atomic>
45#include <cassert>
46#include <functional>
47#include <iostream>
48#include <memory>
49#include <stdexcept>
50#include <string>
51#include <sstream>
52#include <thread>
53#include <unordered_map>
54#include <vector>
55#include <set>
56#include <limits> // For MaxTreeSizeRAII. Revert when #6640 will be solved.
57
58using namespace ROOT::Detail::RDF;
59using namespace ROOT::Internal::RDF;
60
61namespace {
62/// A helper function that returns all RDF code that is currently scheduled for just-in-time compilation.
63/// This allows different RLoopManager instances to share these data.
64/// We want RLoopManagers to be able to add their code to a global "code to execute via cling",
65/// so that, lazily, we can jit everything that's needed by all RDFs in one go, which is potentially
66/// much faster than jitting each RLoopManager's code separately.
67std::string &GetCodeToJit()
68{
69 static std::string code;
70 return code;
71}
72
73bool ContainsLeaf(const std::set<TLeaf *> &leaves, TLeaf *leaf)
74{
75 return (leaves.find(leaf) != leaves.end());
76}
77
78///////////////////////////////////////////////////////////////////////////////
79/// This overload does not check whether the leaf/branch is already in bNamesReg. In case this is a friend leaf/branch,
80/// `allowDuplicates` controls whether we add both `friendname.bname` and `bname` or just the shorter version.
81void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
82 const std::string &friendName, bool allowDuplicates)
83{
84 if (!friendName.empty()) {
85 // In case of a friend tree, users might prepend its name/alias to the branch names
86 const auto friendBName = friendName + "." + branchName;
87 if (bNamesReg.insert(friendBName).second)
88 bNames.push_back(friendBName);
89 }
90
91 if (allowDuplicates || friendName.empty()) {
92 if (bNamesReg.insert(branchName).second)
93 bNames.push_back(branchName);
94 }
95}
96
97///////////////////////////////////////////////////////////////////////////////
98/// This overload makes sure that the TLeaf has not been already inserted.
99void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
100 const std::string &friendName, std::set<TLeaf *> &foundLeaves, TLeaf *leaf,
101 bool allowDuplicates)
102{
103 const bool canAdd = allowDuplicates ? true : !ContainsLeaf(foundLeaves, leaf);
104 if (!canAdd) {
105 return;
106 }
107
108 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
109
110 foundLeaves.insert(leaf);
111}
112
113void ExploreBranch(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames, TBranch *b,
114 std::string prefix, std::string &friendName, bool allowDuplicates)
115{
116 for (auto sb : *b->GetListOfBranches()) {
117 TBranch *subBranch = static_cast<TBranch *>(sb);
118 auto subBranchName = std::string(subBranch->GetName());
119 auto fullName = prefix + subBranchName;
120
121 std::string newPrefix;
122 if (!prefix.empty())
123 newPrefix = fullName + ".";
124
125 ExploreBranch(t, bNamesReg, bNames, subBranch, newPrefix, friendName, allowDuplicates);
126
127 auto branchDirectlyFromTree = t.GetBranch(fullName.c_str());
128 if (!branchDirectlyFromTree)
129 branchDirectlyFromTree = t.FindBranch(fullName.c_str()); // try harder
130 if (branchDirectlyFromTree)
131 InsertBranchName(bNamesReg, bNames, std::string(branchDirectlyFromTree->GetFullName()), friendName,
132 allowDuplicates);
133
134 if (bNamesReg.find(subBranchName) == bNamesReg.end() && t.GetBranch(subBranchName.c_str()))
135 InsertBranchName(bNamesReg, bNames, subBranchName, friendName, allowDuplicates);
136 }
137}
138
139void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
140 std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
141{
142 std::set<TLeaf *> foundLeaves;
143 if (!analysedTrees.insert(&t).second) {
144 return;
145 }
146
147 const auto branches = t.GetListOfBranches();
148 // Getting the branches here triggered the read of the first file of the chain if t is a chain.
149 // We check if a tree has been successfully read, otherwise we throw (see ROOT-9984) to avoid further
150 // operations
151 if (!t.GetTree()) {
152 std::string err("GetBranchNames: error in opening the tree ");
153 err += t.GetName();
154 throw std::runtime_error(err);
155 }
156 if (branches) {
157 for (auto b : *branches) {
158 TBranch *branch = static_cast<TBranch *>(b);
159 const auto branchName = std::string(branch->GetName());
160 if (branch->IsA() == TBranch::Class()) {
161 // Leaf list
162 auto listOfLeaves = branch->GetListOfLeaves();
163 if (listOfLeaves->GetEntriesUnsafe() == 1) {
164 auto leaf = static_cast<TLeaf *>(listOfLeaves->UncheckedAt(0));
165 InsertBranchName(bNamesReg, bNames, branchName, friendName, foundLeaves, leaf, allowDuplicates);
166 }
167
168 for (auto leaf : *listOfLeaves) {
169 auto castLeaf = static_cast<TLeaf *>(leaf);
170 const auto leafName = std::string(leaf->GetName());
171 const auto fullName = branchName + "." + leafName;
172 InsertBranchName(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
173 }
174 } else if (branch->IsA() == TBranchObject::Class()) {
175 // TBranchObject
176 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName, allowDuplicates);
177 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
178 } else {
179 // TBranchElement
180 // Check if there is explicit or implicit dot in the name
181
182 bool dotIsImplied = false;
183 auto be = dynamic_cast<TBranchElement *>(b);
184 if (!be)
185 throw std::runtime_error("GetBranchNames: unsupported branch type");
186 // TClonesArray (3) and STL collection (4)
187 if (be->GetType() == 3 || be->GetType() == 4)
188 dotIsImplied = true;
189
190 if (dotIsImplied || branchName.back() == '.')
191 ExploreBranch(t, bNamesReg, bNames, branch, "", friendName, allowDuplicates);
192 else
193 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName, allowDuplicates);
194
195 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
196 }
197 }
198 }
199
200 // The list of friends needs to be accessed via GetTree()->GetListOfFriends()
201 // (and not via GetListOfFriends() directly), otherwise when `t` is a TChain we
202 // might not recover the list correctly (https://github.com/root-project/root/issues/6741).
203 auto friendTrees = t.GetTree()->GetListOfFriends();
204
205 if (!friendTrees)
206 return;
207
208 for (auto friendTreeObj : *friendTrees) {
209 auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
210
211 std::string frName;
212 auto alias = t.GetFriendAlias(friendTree);
213 if (alias != nullptr)
214 frName = std::string(alias);
215 else
216 frName = std::string(friendTree->GetName());
217
218 GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
219 }
220}
221
222void ThrowIfNSlotsChanged(unsigned int nSlots)
223{
224 const auto currentSlots = RDFInternal::GetNSlots();
225 if (currentSlots != nSlots) {
226 std::string msg = "RLoopManager::Run: when the RDataFrame was constructed the number of slots required was " +
227 std::to_string(nSlots) + ", but when starting the event loop it was " +
228 std::to_string(currentSlots) + ".";
229 if (currentSlots > nSlots)
230 msg += " Maybe EnableImplicitMT() was called after the RDataFrame was constructed?";
231 else
232 msg += " Maybe DisableImplicitMT() was called after the RDataFrame was constructed?";
233 throw std::runtime_error(msg);
234 }
235}
236
237/**
238\struct MaxTreeSizeRAII
239\brief Scope-bound change of `TTree::fgMaxTreeSize`.
240
241This RAII object stores the current value result of `TTree::GetMaxTreeSize`,
242changes it to maximum at construction time and restores it back at destruction
243time. Needed for issue #6523 and should be reverted when #6640 will be solved.
244*/
245struct MaxTreeSizeRAII {
246 Long64_t fOldMaxTreeSize;
247
248 MaxTreeSizeRAII() : fOldMaxTreeSize(TTree::GetMaxTreeSize())
249 {
250 TTree::SetMaxTreeSize(std::numeric_limits<Long64_t>::max());
251 }
252
253 ~MaxTreeSizeRAII() { TTree::SetMaxTreeSize(fOldMaxTreeSize); }
254};
255
256struct DatasetLogInfo {
257 std::string fDataSet;
258 ULong64_t fRangeStart;
259 ULong64_t fRangeEnd;
260 unsigned int fSlot;
261};
262
263std::string LogRangeProcessing(const DatasetLogInfo &info)
264{
265 std::stringstream msg;
266 msg << "Processing " << info.fDataSet << ": entry range [" << info.fRangeStart << "," << info.fRangeEnd - 1
267 << "], using slot " << info.fSlot << " in thread " << std::this_thread::get_id() << '.';
268 return msg.str();
269}
270
271DatasetLogInfo TreeDatasetLogInfo(const TTreeReader &r, unsigned int slot)
272{
273 const auto tree = r.GetTree();
274 const auto chain = dynamic_cast<TChain *>(tree);
275 std::string what;
276 if (chain) {
277 auto files = chain->GetListOfFiles();
278 std::vector<std::string> treeNames;
279 std::vector<std::string> fileNames;
280 for (TObject *f : *files) {
281 treeNames.emplace_back(f->GetName());
282 fileNames.emplace_back(f->GetTitle());
283 }
284 what = "trees {";
285 for (const auto &t : treeNames) {
286 what += t + ",";
287 }
288 what.back() = '}';
289 what += " in files {";
290 for (const auto &f : fileNames) {
291 what += f + ",";
292 }
293 what.back() = '}';
294 } else {
295 assert(tree != nullptr); // to make clang-tidy happy
296 const auto treeName = tree->GetName();
297 what = std::string("tree \"") + treeName + "\"";
298 const auto file = tree->GetCurrentFile();
299 if (file)
300 what += std::string(" in file \"") + file->GetName() + "\"";
301 }
302 const auto entryRange = r.GetEntriesRange();
303 const ULong64_t end = entryRange.second == -1ll ? tree->GetEntries() : entryRange.second;
304 return {std::move(what), static_cast<ULong64_t>(entryRange.first), end, slot};
305}
306
307auto MakeDatasetColReadersKey(const std::string &colName, const std::type_info &ti)
308{
309 // We use a combination of column name and column type name as the key because in some cases we might end up
310 // with concrete readers that use different types for the same column, e.g. std::vector and RVec here:
311 // df.Sum<vector<int>>("stdVectorBranch");
312 // df.Sum<RVecI>("stdVectorBranch");
313 return colName + ':' + ti.name();
314}
315} // anonymous namespace
316
317namespace ROOT {
318namespace Detail {
319namespace RDF {
320
321/// A RAII object that calls RLoopManager::CleanUpTask at destruction
324 unsigned int fArg;
326
327 RCallCleanUpTask(RLoopManager &lm, unsigned int arg = 0u, TTreeReader *reader = nullptr)
328 : fLoopManager(lm), fArg(arg), fReader(reader)
329 {
330 }
332};
333
334} // namespace RDF
335} // namespace Detail
336} // namespace ROOT
337
338///////////////////////////////////////////////////////////////////////////////
339/// Get all the branches names, including the ones of the friend trees
341{
342 std::set<std::string> bNamesSet;
343 ColumnNames_t bNames;
344 std::set<TTree *> analysedTrees;
345 std::string emptyFrName = "";
346 GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
347 return bNames;
348}
349
350RLoopManager::RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
351 : fTree(std::shared_ptr<TTree>(tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
352 fNSlots(RDFInternal::GetNSlots()),
353 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles),
354 fNewSampleNotifier(fNSlots), fSampleInfos(fNSlots), fDatasetColumnReaders(fNSlots)
355{
356}
357
358RLoopManager::RLoopManager(std::unique_ptr<TTree> tree, const ColumnNames_t &defaultBranches)
359 : fTree(std::move(tree)),
360 fDefaultColumns(defaultBranches),
361 fNSlots(RDFInternal::GetNSlots()),
362 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles),
363 fNewSampleNotifier(fNSlots),
364 fSampleInfos(fNSlots),
365 fDatasetColumnReaders(fNSlots)
366{
367}
368
370 : fEmptyEntryRange(0, nEmptyEntries),
371 fNSlots(RDFInternal::GetNSlots()),
372 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles),
373 fNewSampleNotifier(fNSlots),
374 fSampleInfos(fNSlots),
375 fDatasetColumnReaders(fNSlots)
376{
377}
378
379RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
380 : fDefaultColumns(defaultBranches), fNSlots(RDFInternal::GetNSlots()),
381 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
382 fDataSource(std::move(ds)), fNewSampleNotifier(fNSlots), fSampleInfos(fNSlots), fDatasetColumnReaders(fNSlots)
383{
384 fDataSource->SetNSlots(fNSlots);
385}
386
388 : fNSlots(RDFInternal::GetNSlots()),
389 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles),
390 fNewSampleNotifier(fNSlots),
391 fSampleInfos(fNSlots),
392 fDatasetColumnReaders(fNSlots)
393{
394 ChangeSpec(std::move(spec));
395}
396
397/**
398 * @brief Changes the internal TTree held by the RLoopManager.
399 *
400 * @warning This method may lead to potentially dangerous interactions if used
401 * after the construction of the RDataFrame. Changing the specification makes
402 * sense *if and only if* the schema of the dataset is *unchanged*, i.e. the
403 * new specification refers to exactly the same number of columns, with the
404 * same names and types. The actual use case of this method is moving the
405 * processing of the same RDataFrame to a different range of entries of the
406 * same dataset (which may be stored in a different set of files).
407 *
408 * @param spec The specification of the dataset to be adopted.
409 */
411{
412 // Change the range of entries to be processed
413 fBeginEntry = spec.GetEntryRangeBegin();
414 fEndEntry = spec.GetEntryRangeEnd();
415
416 // Store the samples
417 fSamples = spec.MoveOutSamples();
418 fSampleMap.clear();
419
420 // Create the internal main chain
422 for (auto &sample : fSamples) {
423 const auto &trees = sample.GetTreeNames();
424 const auto &files = sample.GetFileNameGlobs();
425 for (std::size_t i = 0ul; i < files.size(); ++i) {
426 // We need to use `<filename>?#<treename>` as an argument to TChain::Add
427 // (see https://github.com/root-project/root/pull/8820 for why)
428 const auto fullpath = files[i] + "?#" + trees[i];
429 chain->Add(fullpath.c_str());
430 // ...but instead we use `<filename>/<treename>` as a sample ID (cannot
431 // change this easily because of backward compatibility: the sample ID
432 // is exposed to users via RSampleInfo and DefinePerSample).
433 const auto sampleId = files[i] + '/' + trees[i];
434 fSampleMap.insert({sampleId, &sample});
435 }
436 }
437 SetTree(std::move(chain));
438
439 // Create friends from the specification and connect them to the main chain
440 const auto &friendInfo = spec.GetFriendInfo();
442 for (std::size_t i = 0ul; i < fFriends.size(); i++) {
443 const auto &thisFriendAlias = friendInfo.fFriendNames[i].second;
444 fTree->AddFriend(fFriends[i].get(), thisFriendAlias.c_str());
445 }
446}
447
448/// Run event loop with no source files, in parallel.
450{
451#ifdef R__USE_IMT
453 // Working with an empty tree.
454 // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
455 const auto nEmptyEntries = GetNEmptyEntries();
456 const auto nEntriesPerSlot = nEmptyEntries / (fNSlots * 2);
457 auto remainder = nEmptyEntries % (fNSlots * 2);
458 std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
459 ULong64_t begin = fEmptyEntryRange.first;
460 while (begin < fEmptyEntryRange.second) {
461 ULong64_t end = begin + nEntriesPerSlot;
462 if (remainder > 0) {
463 ++end;
464 --remainder;
465 }
466 entryRanges.emplace_back(begin, end);
467 begin = end;
468 }
469
470 // Each task will generate a subrange of entries
471 auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
472 ROOT::Internal::RSlotStackRAII slotRAII(slotStack);
473 auto slot = slotRAII.fSlot;
474 RCallCleanUpTask cleanup(*this, slot);
475 InitNodeSlots(nullptr, slot);
476 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({"an empty source", range.first, range.second, slot});
477 try {
478 UpdateSampleInfo(slot, range);
479 for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
480 RunAndCheckFilters(slot, currEntry);
481 }
482 } catch (...) {
483 // Error might throw in experiment frameworks like CMSSW
484 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
485 throw;
486 }
487 };
488
490 pool.Foreach(genFunction, entryRanges);
491
492#endif // not implemented otherwise
493}
494
495/// Run event loop with no source files, in sequence.
497{
498 InitNodeSlots(nullptr, 0);
499 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing(
500 {"an empty source", fEmptyEntryRange.first, fEmptyEntryRange.second, 0u});
501 RCallCleanUpTask cleanup(*this);
502 try {
504 for (ULong64_t currEntry = fEmptyEntryRange.first;
505 currEntry < fEmptyEntryRange.second && fNStopsReceived < fNChildren; ++currEntry) {
506 RunAndCheckFilters(0, currEntry);
507 }
508 } catch (...) {
509 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
510 throw;
511 }
512}
513
514/// Run event loop over one or multiple ROOT files, in parallel.
516{
517#ifdef R__USE_IMT
518 if (fEndEntry == fBeginEntry) // empty range => no work needed
519 return;
521 const auto &entryList = fTree->GetEntryList() ? *fTree->GetEntryList() : TEntryList();
522 auto tp = (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
523 ? std::make_unique<ROOT::TTreeProcessorMT>(*fTree, fNSlots, std::make_pair(fBeginEntry, fEndEntry))
524 : std::make_unique<ROOT::TTreeProcessorMT>(*fTree, entryList, fNSlots);
525
526 std::atomic<ULong64_t> entryCount(0ull);
527
528 tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
529 ROOT::Internal::RSlotStackRAII slotRAII(slotStack);
530 auto slot = slotRAII.fSlot;
531 RCallCleanUpTask cleanup(*this, slot, &r);
532 InitNodeSlots(&r, slot);
533 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing(TreeDatasetLogInfo(r, slot));
534 const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
535 const auto nEntries = entryRange.second - entryRange.first;
536 auto count = entryCount.fetch_add(nEntries);
537 try {
538 // recursive call to check filters and conditionally execute actions
539 while (r.Next()) {
540 if (fNewSampleNotifier.CheckFlag(slot)) {
541 UpdateSampleInfo(slot, r);
542 }
543 RunAndCheckFilters(slot, count++);
544 }
545 } catch (...) {
546 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
547 throw;
548 }
549 // fNStopsReceived < fNChildren is always true at the moment as we don't support event loop early quitting in
550 // multi-thread runs, but it costs nothing to be safe and future-proof in case we add support for that later.
551 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
552 // something went wrong in the TTreeReader event loop
553 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
554 std::to_string(r.GetEntryStatus()));
555 }
556 });
557#endif // no-op otherwise (will not be called)
558}
559
560/// Run event loop over one or multiple ROOT files, in sequence.
562{
563 TTreeReader r(fTree.get(), fTree->GetEntryList());
564 if (0 == fTree->GetEntriesFast() || fBeginEntry == fEndEntry)
565 return;
566 // Apply the range if there is any
567 // In case of a chain with a total of N entries, calling SetEntriesRange(N + 1, ...) does not error out
568 // This is a bug, reported here: https://github.com/root-project/root/issues/10774
569 if (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
570 if (r.SetEntriesRange(fBeginEntry, fEndEntry) != TTreeReader::kEntryValid)
571 throw std::logic_error("Something went wrong in initializing the TTreeReader.");
572
573 RCallCleanUpTask cleanup(*this, 0u, &r);
574 InitNodeSlots(&r, 0);
575 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing(TreeDatasetLogInfo(r, 0u));
576
577 // recursive call to check filters and conditionally execute actions
578 // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
579 try {
580 while (r.Next() && fNStopsReceived < fNChildren) {
582 UpdateSampleInfo(/*slot*/0, r);
583 }
584 RunAndCheckFilters(0, r.GetCurrentEntry());
585 }
586 } catch (...) {
587 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
588 throw;
589 }
590 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
591 // something went wrong in the TTreeReader event loop
592 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
593 std::to_string(r.GetEntryStatus()));
594 }
595}
596
597/// Run event loop over data accessed through a DataSource, in sequence.
599{
600 assert(fDataSource != nullptr);
601 fDataSource->Initialize();
602 auto ranges = fDataSource->GetEntryRanges();
603 while (!ranges.empty() && fNStopsReceived < fNChildren) {
604 InitNodeSlots(nullptr, 0u);
605 fDataSource->InitSlot(0u, 0ull);
606 RCallCleanUpTask cleanup(*this);
607 try {
608 for (const auto &range : ranges) {
609 const auto start = range.first;
610 const auto end = range.second;
611 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, 0u});
612 for (auto entry = start; entry < end && fNStopsReceived < fNChildren; ++entry) {
613 if (fDataSource->SetEntry(0u, entry)) {
614 RunAndCheckFilters(0u, entry);
615 }
616 }
617 }
618 } catch (...) {
619 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
620 throw;
621 }
622 fDataSource->FinalizeSlot(0u);
623 ranges = fDataSource->GetEntryRanges();
624 }
625 fDataSource->Finalize();
626}
627
628/// Run event loop over data accessed through a DataSource, in parallel.
630{
631#ifdef R__USE_IMT
632 assert(fDataSource != nullptr);
635
636 // Each task works on a subrange of entries
637 auto runOnRange = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
638 ROOT::Internal::RSlotStackRAII slotRAII(slotStack);
639 const auto slot = slotRAII.fSlot;
640 InitNodeSlots(nullptr, slot);
641 RCallCleanUpTask cleanup(*this, slot);
642 fDataSource->InitSlot(slot, range.first);
643 const auto start = range.first;
644 const auto end = range.second;
645 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, slot});
646 try {
647 for (auto entry = start; entry < end; ++entry) {
648 if (fDataSource->SetEntry(slot, entry)) {
649 RunAndCheckFilters(slot, entry);
650 }
651 }
652 } catch (...) {
653 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
654 throw;
655 }
656 fDataSource->FinalizeSlot(slot);
657 };
658
659 fDataSource->Initialize();
660 auto ranges = fDataSource->GetEntryRanges();
661 while (!ranges.empty()) {
662 pool.Foreach(runOnRange, ranges);
663 ranges = fDataSource->GetEntryRanges();
664 }
665 fDataSource->Finalize();
666#endif // not implemented otherwise (never called)
667}
668
669/// Execute actions and make sure named filters are called for each event.
670/// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
671void RLoopManager::RunAndCheckFilters(unsigned int slot, Long64_t entry)
672{
673 // data-block callbacks run before the rest of the graph
674 if (fNewSampleNotifier.CheckFlag(slot)) {
675 for (auto &callback : fSampleCallbacks)
676 callback.second(slot, fSampleInfos[slot]);
678 }
679
680 for (auto *actionPtr : fBookedActions)
681 actionPtr->Run(slot, entry);
682 for (auto *namedFilterPtr : fBookedNamedFilters)
683 namedFilterPtr->CheckFilters(slot, entry);
684 for (auto &callback : fCallbacksEveryNEvents)
685 callback(slot);
686}
687
688/// Build TTreeReaderValues for all nodes
689/// This method loops over all filters, actions and other booked objects and
690/// calls their `InitSlot` method, to get them ready for running a task.
692{
693 SetupSampleCallbacks(r, slot);
694 for (auto *ptr : fBookedActions)
695 ptr->InitSlot(r, slot);
696 for (auto *ptr : fBookedFilters)
697 ptr->InitSlot(r, slot);
698 for (auto *ptr : fBookedDefines)
699 ptr->InitSlot(r, slot);
700 for (auto *ptr : fBookedVariations)
701 ptr->InitSlot(r, slot);
702
703 for (auto &callback : fCallbacksOnce)
704 callback(slot);
705}
706
708 if (r != nullptr) {
709 // we need to set a notifier so that we run the callbacks every time we switch to a new TTree
710 // `PrependLink` inserts this notifier into the TTree/TChain's linked list of notifiers
711 fNewSampleNotifier.GetChainNotifyLink(slot).PrependLink(*r->GetTree());
712 }
713 // Whatever the data source, initially set the "new data block" flag:
714 // - for TChains, this ensures that we don't skip the first data block because
715 // the correct tree is already loaded
716 // - for RDataSources and empty sources, which currently don't have data blocks, this
717 // ensures that we run once per task
719}
720
721void RLoopManager::UpdateSampleInfo(unsigned int slot, const std::pair<ULong64_t, ULong64_t> &range) {
723 "Empty source, range: {" + std::to_string(range.first) + ", " + std::to_string(range.second) + "}", range);
724}
725
727 // one GetTree to retrieve the TChain, another to retrieve the underlying TTree
728 auto *tree = r.GetTree()->GetTree();
729 R__ASSERT(tree != nullptr);
730 const std::string treename = ROOT::Internal::TreeUtils::GetTreeFullPaths(*tree)[0];
731 auto *file = tree->GetCurrentFile();
732 const std::string fname = file != nullptr ? file->GetName() : "#inmemorytree#";
733
734 std::pair<Long64_t, Long64_t> range = r.GetEntriesRange();
735 R__ASSERT(range.first >= 0);
736 if (range.second == -1) {
737 range.second = tree->GetEntries(); // convert '-1', i.e. 'until the end', to the actual entry number
738 }
739 const std::string &id = fname + '/' + treename;
740 fSampleInfos[slot] = fSampleMap.empty() ? RSampleInfo(id, range) : RSampleInfo(id, range, fSampleMap[id]);
741}
742
743/// Initialize all nodes of the functional graph before running the event loop.
744/// This method is called once per event-loop and performs generic initialization
745/// operations that do not depend on the specific processing slot (i.e. operations
746/// that are common for all threads).
748{
750 for (auto *filter : fBookedFilters)
751 filter->InitNode();
752 for (auto *range : fBookedRanges)
753 range->InitNode();
754 for (auto *ptr : fBookedActions)
755 ptr->Initialize();
756}
757
758/// Perform clean-up operations. To be called at the end of each event loop.
760{
761 fMustRunNamedFilters = false;
762
763 // forget RActions and detach TResultProxies
764 for (auto *ptr : fBookedActions)
765 ptr->Finalize();
766
767 fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
768 fBookedActions.clear();
769
770 // reset children counts
771 fNChildren = 0;
772 fNStopsReceived = 0;
773 for (auto *ptr : fBookedFilters)
774 ptr->ResetChildrenCount();
775 for (auto *ptr : fBookedRanges)
776 ptr->ResetChildrenCount();
777
779 fCallbacksOnce.clear();
780}
781
782/// Perform clean-up operations. To be called at the end of each task execution.
783void RLoopManager::CleanUpTask(TTreeReader *r, unsigned int slot)
784{
785 if (r != nullptr)
786 fNewSampleNotifier.GetChainNotifyLink(slot).RemoveLink(*r->GetTree());
787 for (auto *ptr : fBookedActions)
788 ptr->FinalizeSlot(slot);
789 for (auto *ptr : fBookedFilters)
790 ptr->FinalizeSlot(slot);
791 for (auto *ptr : fBookedDefines)
792 ptr->FinalizeSlot(slot);
793
795 // we are reading from a tree/chain and we need to re-create the RTreeColumnReaders at every task
796 // because the TTreeReader object changes at every task
797 for (auto &v : fDatasetColumnReaders[slot])
798 v.second.reset();
799 }
800}
801
802/// Add RDF nodes that require just-in-time compilation to the computation graph.
803/// This method also clears the contents of GetCodeToJit().
805{
806 // TODO this should be a read lock unless we find GetCodeToJit non-empty
808
809 const std::string code = std::move(GetCodeToJit());
810 if (code.empty()) {
811 R__LOG_INFO(RDFLogChannel()) << "Nothing to jit and execute.";
812 return;
813 }
814
815 TStopwatch s;
816 s.Start();
817 RDFInternal::InterpreterCalc(code, "RLoopManager::Run");
818 s.Stop();
819 R__LOG_INFO(RDFLogChannel()) << "Just-in-time compilation phase completed"
820 << (s.RealTime() > 1e-3 ? " in " + std::to_string(s.RealTime()) + " seconds."
821 : " in less than 1ms.");
822}
823
824/// Trigger counting of number of children nodes for each node of the functional graph.
825/// This is done once before starting the event loop. Each action sends an `increase children count` signal
826/// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
827/// children counter. Each node only propagates the signal once, even if it receives it multiple times.
828/// Named filters also send an `increase children count` signal, just like actions, as they always execute during
829/// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
831{
832 for (auto *actionPtr : fBookedActions)
833 actionPtr->TriggerChildrenCount();
834 for (auto *namedFilterPtr : fBookedNamedFilters)
835 namedFilterPtr->TriggerChildrenCount();
836}
837
838/// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
839/// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
840/// The jitting phase is skipped if the `jit` parameter is `false` (unsafe, use with care).
841void RLoopManager::Run(bool jit)
842{
843 // Change value of TTree::GetMaxTreeSize only for this scope. Revert when #6640 will be solved.
844 MaxTreeSizeRAII ctxtmts;
845
846 R__LOG_INFO(RDFLogChannel()) << "Starting event loop number " << fNRuns << '.';
847
848 ThrowIfNSlotsChanged(GetNSlots());
849
850 if (jit)
851 Jit();
852
853 InitNodes();
854
855 // Exceptions can occur during the event loop. In order to ensure proper cleanup of nodes
856 // we use RAII: even in case of an exception, the destructor of the object is invoked and
857 // all the cleanup takes place.
858 class NodesCleanerRAII {
859 RLoopManager &fRLM;
860
861 public:
862 NodesCleanerRAII(RLoopManager &thisRLM) : fRLM(thisRLM) {}
863 ~NodesCleanerRAII() { fRLM.CleanUpNodes(); }
864 };
865
866 NodesCleanerRAII runKeeper(*this);
867
868 TStopwatch s;
869 s.Start();
870
871 switch (fLoopType) {
875 case ELoopType::kNoFiles: RunEmptySource(); break;
878 }
879 s.Stop();
880
881 fNRuns++;
882
883 R__LOG_INFO(RDFLogChannel()) << "Finished event loop number " << fNRuns - 1 << " (" << s.CpuTime() << "s CPU, "
884 << s.RealTime() << "s elapsed).";
885}
886
887/// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
889{
890 return fDefaultColumns;
891}
892
894{
895 return fTree.get();
896}
897
899{
900 fBookedActions.emplace_back(actionPtr);
901 AddSampleCallback(actionPtr, actionPtr->GetSampleCallback());
902}
903
905{
908 fSampleCallbacks.erase(actionPtr);
909}
910
912{
913 fBookedFilters.emplace_back(filterPtr);
914 if (filterPtr->HasName()) {
915 fBookedNamedFilters.emplace_back(filterPtr);
917 }
918}
919
921{
924}
925
927{
928 fBookedRanges.emplace_back(rangePtr);
929}
930
932{
934}
935
937{
938 fBookedDefines.emplace_back(ptr);
939}
940
942{
944 fSampleCallbacks.erase(ptr);
945}
946
948{
949 fBookedVariations.emplace_back(v);
950}
951
953{
955}
956
957// dummy call, end of recursive chain of calls
959{
960 return true;
961}
962
963/// Call `FillReport` on all booked filters
965{
966 for (const auto *fPtr : fBookedNamedFilters)
967 fPtr->FillReport(rep);
968}
969
970void RLoopManager::SetTree(std::shared_ptr<TTree> tree)
971{
972 fTree = std::move(tree);
973
974 TChain *ch = nullptr;
975 if ((ch = dynamic_cast<TChain *>(fTree.get())))
977}
978
979void RLoopManager::ToJitExec(const std::string &code) const
980{
982 GetCodeToJit().append(code);
983}
984
985void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
986{
987 if (everyNEvents == 0ull)
988 fCallbacksOnce.emplace_back(std::move(f), fNSlots);
989 else
990 fCallbacksEveryNEvents.emplace_back(everyNEvents, std::move(f), fNSlots);
991}
992
993std::vector<std::string> RLoopManager::GetFiltersNames()
994{
995 std::vector<std::string> filters;
996 for (auto *filter : fBookedFilters) {
997 auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
998 filters.push_back(name);
999 }
1000 return filters;
1001}
1002
1003std::vector<RNodeBase *> RLoopManager::GetGraphEdges() const
1004{
1005 std::vector<RNodeBase *> nodes(fBookedFilters.size() + fBookedRanges.size());
1006 auto it = std::copy(fBookedFilters.begin(), fBookedFilters.end(), nodes.begin());
1007 std::copy(fBookedRanges.begin(), fBookedRanges.end(), it);
1008 return nodes;
1009}
1010
1011std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions() const
1012{
1013 std::vector<RDFInternal::RActionBase *> actions(fBookedActions.size() + fRunActions.size());
1014 auto it = std::copy(fBookedActions.begin(), fBookedActions.end(), actions.begin());
1015 std::copy(fRunActions.begin(), fRunActions.end(), it);
1016 return actions;
1017}
1018
1019std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph(
1020 std::unordered_map<void *, std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode>> &visitedMap)
1021{
1022 // If there is already a node for the RLoopManager return it. If there is not, return a new one.
1023 auto duplicateRLoopManagerIt = visitedMap.find((void *)this);
1024 if (duplicateRLoopManagerIt != visitedMap.end())
1025 return duplicateRLoopManagerIt->second;
1026
1027 std::string name;
1028 if (fDataSource) {
1029 name = fDataSource->GetLabel();
1030 } else if (fTree) {
1031 name = fTree->GetName();
1032 if (name.empty())
1033 name = fTree->ClassName();
1034 } else {
1035 name = "Empty source\\nEntries: " + std::to_string(GetNEmptyEntries());
1036 }
1037 auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(
1039 visitedMap[(void *)this] = thisNode;
1040 return thisNode;
1041}
1042
1043////////////////////////////////////////////////////////////////////////////
1044/// Return all valid TTree::Branch names (caching results for subsequent calls).
1045/// Never use fBranchNames directy, always request it through this method.
1047{
1048 if (fValidBranchNames.empty() && fTree) {
1049 fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
1050 }
1051 return fValidBranchNames;
1052}
1053
1054/// Return true if AddDataSourceColumnReaders was called for column name col.
1055bool RLoopManager::HasDataSourceColumnReaders(const std::string &col, const std::type_info &ti) const
1056{
1057 const auto key = MakeDatasetColReadersKey(col, ti);
1058 assert(fDataSource != nullptr);
1059 // since data source column readers are always added for all slots at the same time,
1060 // if the reader is present for slot 0 we have it for all other slots as well.
1061 return fDatasetColumnReaders[0].find(key) != fDatasetColumnReaders[0].end();
1062}
1063
1065 std::vector<std::unique_ptr<RColumnReaderBase>> &&readers,
1066 const std::type_info &ti)
1067{
1068 const auto key = MakeDatasetColReadersKey(col, ti);
1069 assert(fDataSource != nullptr && !HasDataSourceColumnReaders(col, ti));
1070 assert(readers.size() == fNSlots);
1071
1072 for (auto slot = 0u; slot < fNSlots; ++slot) {
1073 fDatasetColumnReaders[slot][key] = std::move(readers[slot]);
1074 }
1075}
1076
1077// Differently from AddDataSourceColumnReaders, this can be called from multiple threads concurrently
1078/// \brief Register a new RTreeColumnReader with this RLoopManager.
1079/// \return A shared pointer to the inserted column reader.
1080RColumnReaderBase *RLoopManager::AddTreeColumnReader(unsigned int slot, const std::string &col,
1081 std::unique_ptr<RColumnReaderBase> &&reader,
1082 const std::type_info &ti)
1083{
1084 auto &readers = fDatasetColumnReaders[slot];
1085 const auto key = MakeDatasetColReadersKey(col, ti);
1086 // if a reader for this column and this slot was already there, we are doing something wrong
1087 assert(readers.find(key) == readers.end() || readers[key] == nullptr);
1088 auto *rptr = reader.get();
1089 readers[key] = std::move(reader);
1090 return rptr;
1091}
1092
1094RLoopManager::GetDatasetColumnReader(unsigned int slot, const std::string &col, const std::type_info &ti) const
1095{
1096 const auto key = MakeDatasetColReadersKey(col, ti);
1097 auto it = fDatasetColumnReaders[slot].find(key);
1098 if (it != fDatasetColumnReaders[slot].end())
1099 return it->second.get();
1100 else
1101 return nullptr;
1102}
1103
1105{
1106 if (callback)
1107 fSampleCallbacks.insert({nodePtr, std::move(callback)});
1108}
1109
1110void RLoopManager::SetEmptyEntryRange(std::pair<ULong64_t, ULong64_t> &&newRange)
1111{
1112 fEmptyEntryRange = std::move(newRange);
1113}
1114
1115/**
1116 * \brief Helper function to open a file (or the first file from a glob).
1117 * This function is used at construction time of an RDataFrame, to check the
1118 * concrete type of the dataset stored inside the file.
1119 */
1120std::unique_ptr<TFile> OpenFileWithSanityChecks(std::string_view fileNameGlob)
1121{
1122 bool fileIsGlob = [&fileNameGlob]() {
1123 const std::vector<std::string_view> wildcards = {"[", "]", "*", "?"}; // Wildcards accepted by TChain::Add
1124 return std::any_of(wildcards.begin(), wildcards.end(),
1125 [&fileNameGlob](const auto &wc) { return fileNameGlob.find(wc) != std::string_view::npos; });
1126 }();
1127
1128 // Open first file in case of glob, suppose all files in the glob use the same data format
1129 std::string fileToOpen{fileIsGlob ? ROOT::Internal::TreeUtils::ExpandGlob(std::string{fileNameGlob})[0]
1130 : fileNameGlob};
1131
1132 ::TDirectory::TContext ctxt; // Avoid changing gDirectory;
1133 std::unique_ptr<TFile> inFile{TFile::Open(fileToOpen.c_str(), "READ_WITHOUT_GLOBALREGISTRATION")};
1134 if (!inFile || inFile->IsZombie())
1135 throw std::invalid_argument("RDataFrame: could not open file \"" + fileToOpen + "\".");
1136
1137 return inFile;
1138}
1139
1140std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1141ROOT::Detail::RDF::CreateLMFromTTree(std::string_view datasetName, std::string_view fileNameGlob,
1142 const ROOT::RDF::ColumnNames_t &defaultColumns, bool checkFile)
1143{
1144 // Introduce the same behaviour as in CreateLMFromFile for consistency.
1145 // Creating an RDataFrame with a non-existing file will throw early rather
1146 // than wait for the start of the graph execution.
1147 if (checkFile) {
1148 OpenFileWithSanityChecks(fileNameGlob);
1149 }
1150 std::string datasetNameInt{datasetName};
1151 std::string fileNameGlobInt{fileNameGlob};
1152 auto chain = ROOT::Internal::TreeUtils::MakeChainForMT(datasetNameInt.c_str());
1153 chain->Add(fileNameGlobInt.c_str());
1154 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(chain), defaultColumns);
1155 return lm;
1156}
1157
1158std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1159ROOT::Detail::RDF::CreateLMFromTTree(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1160 const std::vector<std::string> &defaultColumns, bool checkFile)
1161{
1162 // Introduce the same behaviour as in CreateLMFromFile for consistency.
1163 // Creating an RDataFrame with a non-existing file will throw early rather
1164 // than wait for the start of the graph execution.
1165 if (checkFile) {
1166 OpenFileWithSanityChecks(fileNameGlobs[0]);
1167 }
1168 std::string treeNameInt(datasetName);
1169 auto chain = ROOT::Internal::TreeUtils::MakeChainForMT(treeNameInt);
1170 for (auto &f : fileNameGlobs)
1171 chain->Add(f.c_str());
1172 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(chain), defaultColumns);
1173 return lm;
1174}
1175
1176#ifdef R__HAS_ROOT7
1177std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1178ROOT::Detail::RDF::CreateLMFromRNTuple(std::string_view datasetName, std::string_view fileNameGlob,
1179 const ROOT::RDF::ColumnNames_t &defaultColumns)
1180{
1181 auto dataSource = std::make_unique<ROOT::Experimental::RNTupleDS>(datasetName, fileNameGlob);
1182 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(dataSource), defaultColumns);
1183 return lm;
1184}
1185
1186std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1187ROOT::Detail::RDF::CreateLMFromRNTuple(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1188 const ROOT::RDF::ColumnNames_t &defaultColumns)
1189{
1190 auto dataSource = std::make_unique<ROOT::Experimental::RNTupleDS>(datasetName, fileNameGlobs);
1191 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(dataSource), defaultColumns);
1192 return lm;
1193}
1194
1195std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1196ROOT::Detail::RDF::CreateLMFromFile(std::string_view datasetName, std::string_view fileNameGlob,
1197 const ROOT::RDF::ColumnNames_t &defaultColumns)
1198{
1199
1200 auto inFile = OpenFileWithSanityChecks(fileNameGlob);
1201
1202 if (inFile->Get<TTree>(datasetName.data())) {
1203 return CreateLMFromTTree(datasetName, fileNameGlob, defaultColumns, /*checkFile=*/false);
1204 } else if (inFile->Get<ROOT::Experimental::RNTuple>(datasetName.data())) {
1205 return CreateLMFromRNTuple(datasetName, fileNameGlob, defaultColumns);
1206 }
1207
1208 throw std::invalid_argument("RDataFrame: unsupported data format for dataset \"" + std::string(datasetName) +
1209 "\" in file \"" + inFile->GetName() + "\".");
1210}
1211
1212std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1213ROOT::Detail::RDF::CreateLMFromFile(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1214 const ROOT::RDF::ColumnNames_t &defaultColumns)
1215{
1216
1217 auto inFile = OpenFileWithSanityChecks(fileNameGlobs[0]);
1218
1219 if (inFile->Get<TTree>(datasetName.data())) {
1220 return CreateLMFromTTree(datasetName, fileNameGlobs, defaultColumns, /*checkFile=*/false);
1221 } else if (inFile->Get<ROOT::Experimental::RNTuple>(datasetName.data())) {
1222 return CreateLMFromRNTuple(datasetName, fileNameGlobs, defaultColumns);
1223 }
1224
1225 throw std::invalid_argument("RDataFrame: unsupported data format for dataset \"" + std::string(datasetName) +
1226 "\" in file \"" + inFile->GetName() + "\".");
1227}
1228#endif
#define R__LOG_DEBUG(DEBUGLEVEL,...)
Definition RLogger.hxx:365
#define R__LOG_INFO(...)
Definition RLogger.hxx:364
std::unique_ptr< TFile > OpenFileWithSanityChecks(std::string_view fileNameGlob)
Helper function to open a file (or the first file from a glob).
#define b(i)
Definition RSha256.hxx:100
#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 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 SetEmptyEntryRange(std::pair< ULong64_t, ULong64_t > &&newRange)
void Report(ROOT::RDF::RCutFlowReport &rep) const final
Call FillReport on all booked filters.
void AddSampleCallback(void *nodePtr, ROOT::RDF::SampleCallback_t &&callback)
std::vector< RFilterBase * > fBookedNamedFilters
Contains a subset of fBookedFilters, i.e. only the named filters.
void RunEmptySourceMT()
Run event loop with no source files, in parallel.
std::unordered_map< std::string, ROOT::RDF::Experimental::RSample * > fSampleMap
Keys are fname + "/" + treename as RSampleInfo::fID; Values are pointers to the corresponding sample.
std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > GetGraph(std::unordered_map< void *, std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > > &visitedMap) final
const ColumnNames_t & GetBranchNames()
Return all valid TTree::Branch names (caching results for subsequent calls).
void ToJitExec(const std::string &) const
std::vector< RDFInternal::RActionBase * > GetAllActions() const
Return all actions, either booked or already run.
std::vector< ROOT::RDF::RSampleInfo > fSampleInfos
void ChangeSpec(ROOT::RDF::Experimental::RDatasetSpec &&spec)
Changes the internal TTree held by the RLoopManager.
void SetTree(std::shared_ptr< TTree > tree)
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.
ROOT::Internal::TreeUtils::RNoCleanupNotifier fNoCleanupNotifier
std::vector< RDFInternal::RActionBase * > fRunActions
Non-owning pointers to actions already run.
RColumnReaderBase * GetDatasetColumnReader(unsigned int slot, const std::string &col, const std::type_info &ti) const
std::vector< RRangeBase * > fBookedRanges
std::vector< ROOT::RDF::Experimental::RSample > fSamples
Samples need to survive throughout the whole event loop, hence stored as an attribute.
std::vector< std::string > ColumnNames_t
void RunAndCheckFilters(unsigned int slot, Long64_t entry)
Execute actions and make sure named filters are called for each event.
std::vector< RFilterBase * > fBookedFilters
void Run(bool jit=true)
Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
std::unordered_map< void *, ROOT::RDF::SampleCallback_t > fSampleCallbacks
Registered callbacks to call at the beginning of each "data block".
std::vector< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
RColumnReaderBase * AddTreeColumnReader(unsigned int slot, const std::string &col, std::unique_ptr< RColumnReaderBase > &&reader, const std::type_info &ti)
Register a new RTreeColumnReader with this RLoopManager.
const ELoopType fLoopType
The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
void AddDataSourceColumnReaders(const std::string &col, std::vector< std::unique_ptr< RColumnReaderBase > > &&readers, const std::type_info &ti)
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::vector< RDFInternal::RCallback > fCallbacksEveryNEvents
Registered callbacks to be executed every N events.
std::vector< std::unordered_map< std::string, std::unique_ptr< RColumnReaderBase > > > fDatasetColumnReaders
Readers for TTree/RDataSource columns (one per slot), shared by all nodes in the computation graph.
void Register(RDFInternal::RActionBase *actionPtr)
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.
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.
RDFInternal::RNewSampleNotifier fNewSampleNotifier
std::pair< ULong64_t, ULong64_t > fEmptyEntryRange
Range of entries created when no data source is specified.
const ColumnNames_t fDefaultColumns
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.
void Deregister(RDFInternal::RActionBase *actionPtr)
void InitNodes()
Initialize all nodes of the functional graph before running the event loop.
std::vector< std::unique_ptr< TChain > > fFriends
Friends of the fTree. Only used if we constructed fTree ourselves.
bool HasDataSourceColumnReaders(const std::string &col, const std::type_info &ti) const
Return true if AddDataSourceColumnReaders was called for column name col.
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
Representation of an RNTuple data set in a ROOT file.
Definition RNTuple.hxx:61
virtual ROOT::RDF::SampleCallback_t GetSampleCallback()=0
bool CheckFlag(unsigned int slot) const
TNotifyLink< RNewSampleFlag > & GetChainNotifyLink(unsigned int slot)
This type includes all parts of RVariation that do not depend on the callable signature.
A thread-safe stack of N indexes (0 to size - 1).
The dataset specification for RDataFrame.
This type represents a sample identifier, to be used in conjunction with RDataFrame features such as ...
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:93
static TClass * Class()
TClass * IsA() const override
Definition TBranch.h:295
TObjArray * GetListOfLeaves()
Definition TBranch.h:247
A chain is a collection of files containing TTree objects.
Definition TChain.h:33
TObjArray * GetListOfFiles() const
Definition TChain.h:111
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
A List of entry numbers in a TTree or TChain.
Definition TEntryList.h:26
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4082
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:41
Stopwatch class.
Definition TStopwatch.h:28
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
void Stop()
Stop the stopwatch.
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
@ kEntryValid
data read okay
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:4841
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition TTree.cxx:5294
static void SetMaxTreeSize(Long64_t maxsize=100000000000LL)
Set the maximum size in bytes of a Tree file (static function).
Definition TTree.cxx:9197
virtual TObjArray * GetListOfBranches()
Definition TTree.h:488
virtual TTree * GetTree() const
Definition TTree.h:517
virtual TList * GetListOfFriends() const
Definition TTree.h:490
virtual const char * GetFriendAlias(TTree *) const
If the 'tree' is a friend, this method returns its alias name.
Definition TTree.cxx:6032
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > CreateLMFromTTree(std::string_view datasetName, std::string_view fileNameGlob, const std::vector< std::string > &defaultColumns, bool checkFile=true)
Create an RLoopManager that reads a TChain.
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:283
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:327
std::vector< std::string > GetTreeFullPaths(const TTree &tree)
std::unique_ptr< TChain > MakeChainForMT(const std::string &name="", const std::string &title="")
Create a TChain object with options that avoid common causes of thread contention.
std::vector< std::unique_ptr< TChain > > MakeFriends(const ROOT::TreeUtils::RFriendInfo &finfo)
Create friends from the main TTree.
std::vector< std::string > ExpandGlob(const std::string &glob)
Expands input glob into a collection of full paths to files.
std::vector< std::string > ColumnNames_t
std::function< void(unsigned int, const ROOT::RDF::RSampleInfo &)> SampleCallback_t
The type of a data-block callback, registered with an RDataFrame computation graph via e....
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:568
static const char * what
Definition stlLoader.cc:5
A RAII object that calls RLoopManager::CleanUpTask at destruction.
RCallCleanUpTask(RLoopManager &lm, unsigned int arg=0u, TTreeReader *reader=nullptr)
A RAII object to pop and push slot numbers from a RSlotStack object.