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