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 "ROOT/RNTuple.hxx"
23#include "ROOT/RNTupleDS.hxx"
24#include "RtypesCore.h" // Long64_t
25#include "TStopwatch.h"
26#include "TBranchElement.h"
27#include "TBranchObject.h"
28#include "TChain.h"
29#include "TEntryList.h"
30#include "TFile.h"
31#include "TFriendElement.h"
32#include "TROOT.h" // IsImplicitMTEnabled, gCoreMutex, R__*_LOCKGUARD
33#include "TTreeReader.h"
34#include "TTree.h" // For MaxTreeSizeRAII. Revert when #6640 will be solved.
35
36#include "ROOT/RTTreeDS.hxx"
37
38#ifdef R__USE_IMT
41#include "ROOT/RSlotStack.hxx"
42#endif
43
44#ifdef R__UNIX
45// Functions needed to perform EOS XRootD redirection in ChangeSpec
46#include "TEnv.h"
47#include "TSystem.h"
48#ifndef R__FBSD
49#include <sys/xattr.h>
50#else
51#include <sys/extattr.h>
52#endif
53#ifdef R__MACOSX
54/* On macOS getxattr takes two extra arguments that should be set to 0 */
55#define getxattr(path, name, value, size) getxattr(path, name, value, size, 0u, 0)
56#endif
57#ifdef R__FBSD
58#define getxattr(path, name, value, size) extattr_get_file(path, EXTATTR_NAMESPACE_USER, name, value, size)
59#endif
60#endif
61
62#include <algorithm>
63#include <atomic>
64#include <cassert>
65#include <functional>
66#include <iostream>
67#include <memory>
68#include <stdexcept>
69#include <string>
70#include <sstream>
71#include <thread>
72#include <unordered_map>
73#include <vector>
74#include <set>
75#include <limits> // For MaxTreeSizeRAII. Revert when #6640 will be solved.
76
77using namespace ROOT::Detail::RDF;
78using namespace ROOT::Internal::RDF;
79
80namespace {
81/// A helper function that returns all RDF code that is currently scheduled for just-in-time compilation.
82/// This allows different RLoopManager instances to share these data.
83/// We want RLoopManagers to be able to add their code to a global "code to execute via cling",
84/// so that, lazily, we can jit everything that's needed by all RDFs in one go, which is potentially
85/// much faster than jitting each RLoopManager's code separately.
86std::string &GetCodeToJit()
87{
88 static std::string code;
89 return code;
90}
91
92bool ContainsLeaf(const std::set<TLeaf *> &leaves, TLeaf *leaf)
93{
94 return (leaves.find(leaf) != leaves.end());
95}
96
97///////////////////////////////////////////////////////////////////////////////
98/// This overload does not check whether the leaf/branch is already in bNamesReg. In case this is a friend leaf/branch,
99/// `allowDuplicates` controls whether we add both `friendname.bname` and `bname` or just the shorter version.
100void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
101 const std::string &friendName, bool allowDuplicates)
102{
103 if (!friendName.empty()) {
104 // In case of a friend tree, users might prepend its name/alias to the branch names
105 const auto friendBName = friendName + "." + branchName;
106 if (bNamesReg.insert(friendBName).second)
107 bNames.push_back(friendBName);
108 }
109
110 if (allowDuplicates || friendName.empty()) {
111 if (bNamesReg.insert(branchName).second)
112 bNames.push_back(branchName);
113 }
114}
115
116///////////////////////////////////////////////////////////////////////////////
117/// This overload makes sure that the TLeaf has not been already inserted.
118void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
119 const std::string &friendName, std::set<TLeaf *> &foundLeaves, TLeaf *leaf,
120 bool allowDuplicates)
121{
123 if (!canAdd) {
124 return;
125 }
126
128
129 foundLeaves.insert(leaf);
130}
131
132void ExploreBranch(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames, TBranch *b,
133 std::string prefix, std::string &friendName, bool allowDuplicates)
134{
135 // We want to avoid situations of overlap between the prefix and the
136 // sub-branch name that might happen when the branch is composite, e.g.
137 // prefix=reco_Wdecay2_from_tbar_4vect_NOSYS.fCoordinates.
138 // subBranchName=fCoordinates.fPt
139 // which would lead to a repetition of fCoordinates in the output branch name
140 // Boundary to search for the token before the last dot
141 auto prefixEndingDot = std::string::npos;
142 if (!prefix.empty() && prefix.back() == '.')
143 prefixEndingDot = prefix.size() - 2;
144 std::string lastPrefixToken{};
145 if (auto prefixLastRealDot = prefix.find_last_of('.', prefixEndingDot); prefixLastRealDot != std::string::npos)
147
148 for (auto sb : *b->GetListOfBranches()) {
149 TBranch *subBranch = static_cast<TBranch *>(sb);
150 auto subBranchName = std::string(subBranch->GetName());
151 auto fullName = prefix + subBranchName;
152
153 if (auto subNameFirstDot = subBranchName.find_first_of('.'); subNameFirstDot != std::string::npos) {
154 // Concatenate the prefix to the sub-branch name without overlaps
155 if (!lastPrefixToken.empty() && lastPrefixToken == subBranchName.substr(0, subNameFirstDot))
156 fullName = prefix + subBranchName.substr(subNameFirstDot + 1);
157 }
158
159 std::string newPrefix;
160 if (!prefix.empty())
161 newPrefix = fullName + ".";
162
164
165 auto branchDirectlyFromTree = t.GetBranch(fullName.c_str());
167 branchDirectlyFromTree = t.FindBranch(fullName.c_str()); // try harder
171
172 if (bNamesReg.find(subBranchName) == bNamesReg.end() && t.GetBranch(subBranchName.c_str()))
174 }
175}
176
177void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
178 std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
179{
180 std::set<TLeaf *> foundLeaves;
181 if (!analysedTrees.insert(&t).second) {
182 return;
183 }
184
185 const auto branches = t.GetListOfBranches();
186 // Getting the branches here triggered the read of the first file of the chain if t is a chain.
187 // We check if a tree has been successfully read, otherwise we throw (see ROOT-9984) to avoid further
188 // operations
189 if (!t.GetTree()) {
190 std::string err("GetBranchNames: error in opening the tree ");
191 err += t.GetName();
192 throw std::runtime_error(err);
193 }
194 if (branches) {
195 for (auto b : *branches) {
196 TBranch *branch = static_cast<TBranch *>(b);
197 const auto branchName = std::string(branch->GetName());
198 if (branch->IsA() == TBranch::Class()) {
199 // Leaf list
200 auto listOfLeaves = branch->GetListOfLeaves();
201 if (listOfLeaves->GetEntriesUnsafe() == 1) {
202 auto leaf = static_cast<TLeaf *>(listOfLeaves->UncheckedAt(0));
204 }
205
206 for (auto leaf : *listOfLeaves) {
207 auto castLeaf = static_cast<TLeaf *>(leaf);
208 const auto leafName = std::string(leaf->GetName());
209 const auto fullName = branchName + "." + leafName;
211 }
212 } else if (branch->IsA() == TBranchObject::Class()) {
213 // TBranchObject
216 } else {
217 // TBranchElement
218 // Check if there is explicit or implicit dot in the name
219
220 bool dotIsImplied = false;
221 auto be = dynamic_cast<TBranchElement *>(b);
222 if (!be)
223 throw std::runtime_error("GetBranchNames: unsupported branch type");
224 // TClonesArray (3) and STL collection (4)
225 if (be->GetType() == 3 || be->GetType() == 4)
226 dotIsImplied = true;
227
228 if (dotIsImplied || branchName.back() == '.')
230 else
232
234 }
235 }
236 }
237
238 // The list of friends needs to be accessed via GetTree()->GetListOfFriends()
239 // (and not via GetListOfFriends() directly), otherwise when `t` is a TChain we
240 // might not recover the list correctly (https://github.com/root-project/root/issues/6741).
241 auto friendTrees = t.GetTree()->GetListOfFriends();
242
243 if (!friendTrees)
244 return;
245
246 for (auto friendTreeObj : *friendTrees) {
247 auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
248
249 std::string frName;
251 if (alias != nullptr)
252 frName = std::string(alias);
253 else
254 frName = std::string(friendTree->GetName());
255
257 }
258}
259
260void ThrowIfNSlotsChanged(unsigned int nSlots)
261{
263 if (currentSlots != nSlots) {
264 std::string msg = "RLoopManager::Run: when the RDataFrame was constructed the number of slots required was " +
265 std::to_string(nSlots) + ", but when starting the event loop it was " +
266 std::to_string(currentSlots) + ".";
267 if (currentSlots > nSlots)
268 msg += " Maybe EnableImplicitMT() was called after the RDataFrame was constructed?";
269 else
270 msg += " Maybe DisableImplicitMT() was called after the RDataFrame was constructed?";
271 throw std::runtime_error(msg);
272 }
273}
274
275/**
276\struct MaxTreeSizeRAII
277\brief Scope-bound change of `TTree::fgMaxTreeSize`.
278
279This RAII object stores the current value result of `TTree::GetMaxTreeSize`,
280changes it to maximum at construction time and restores it back at destruction
281time. Needed for issue #6523 and should be reverted when #6640 will be solved.
282*/
283struct MaxTreeSizeRAII {
284 Long64_t fOldMaxTreeSize;
285
286 MaxTreeSizeRAII() : fOldMaxTreeSize(TTree::GetMaxTreeSize())
287 {
288 TTree::SetMaxTreeSize(std::numeric_limits<Long64_t>::max());
289 }
290
291 ~MaxTreeSizeRAII() { TTree::SetMaxTreeSize(fOldMaxTreeSize); }
292};
293
294struct DatasetLogInfo {
295 std::string fDataSet;
296 ULong64_t fRangeStart;
297 ULong64_t fRangeEnd;
298 unsigned int fSlot;
299};
300
301std::string LogRangeProcessing(const DatasetLogInfo &info)
302{
303 std::stringstream msg;
304 msg << "Processing " << info.fDataSet << ": entry range [" << info.fRangeStart << "," << info.fRangeEnd - 1
305 << "], using slot " << info.fSlot << " in thread " << std::this_thread::get_id() << '.';
306 return msg.str();
307}
308
309DatasetLogInfo TreeDatasetLogInfo(const TTreeReader &r, unsigned int slot)
310{
311 const auto r.GetTree();
312 const auto tree);
313 std::string what;
314 if (chain) {
315 auto files = chain->GetListOfFiles();
316 std::vector<std::string> treeNames;
317 std::vector<std::string> fileNames;
318 for (TObject *f : *files) {
319 treeNames.emplace_back(f->GetName());
320 fileNames.emplace_back(f->GetTitle());
321 }
322 what = "trees {";
323 for (const auto &t : treeNames) {
324 what += t + ",";
325 }
326 what.back() = '}';
327 what += " in files {";
328 for (const auto &f : fileNames) {
329 what += f + ",";
330 }
331 what.back() = '}';
332 } else {
333 assert(tree != nullptr); // to make clang-tidy happy
334 const auto tree->GetName();
335 what = std::string("tree \"") + treeName + "\"";
336 const auto file = tree->GetCurrentFile();
337 if (file)
338 what += std::string(" in file \"") + file->GetName() + "\"";
339 }
340 const auto entryRange = r.GetEntriesRange();
341 const entryRange.second;
342 return {std::move(what), static_cast<ULong64_t>(entryRange.first), end, slot};
343}
344
345auto MakeDatasetColReadersKey(std::string_view colName, const std::type_info &ti)
346{
347 // We use a combination of column name and column type name as the key because in some cases we might end up
348 // with concrete readers that use different types for the same column, e.g. std::vector and RVec here:
349 // df.Sum<vector<int>>("stdVectorBranch");
350 // df.Sum<RVecI>("stdVectorBranch");
351 return std::string(colName) + ':' + ti.name();
352}
353} // anonymous namespace
354
355namespace ROOT {
356namespace Detail {
357namespace RDF {
358
359/// A RAII object that calls RLoopManager::CleanUpTask at destruction
371
372} // namespace RDF
373} // namespace Detail
374} // namespace ROOT
375
376///////////////////////////////////////////////////////////////////////////////
377/// Get all the branches names, including the ones of the friend trees
379{
380 std::set<std::string> bNamesSet;
382 std::set<TTree *> analysedTrees;
383 std::string emptyFrName = "";
385 return bNames;
386}
387
388ROOT::Detail::RDF::RLoopManager::RLoopManager(const ROOT::Detail::RDF::ColumnNames_t &defaultColumns)
389 : fDefaultColumns(defaultColumns),
390 fNSlots(RDFInternal::GetNSlots()),
391 fNewSampleNotifier(fNSlots),
392 fSampleInfos(fNSlots),
393 fDatasetColumnReaders(fNSlots)
394{
395}
396
398 : fDefaultColumns(defaultBranches),
399 fNSlots(RDFInternal::GetNSlots()),
400 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
401 fDataSource(std::make_unique<ROOT::Internal::RDF::RTTreeDS>(ROOT::Internal::RDF::MakeAliasedSharedPtr(tree))),
402 fNewSampleNotifier(fNSlots),
403 fSampleInfos(fNSlots),
404 fDatasetColumnReaders(fNSlots)
405{
406 fDataSource->SetNSlots(fNSlots);
407}
408
410 : fEmptyEntryRange(0, nEmptyEntries),
411 fNSlots(RDFInternal::GetNSlots()),
412 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles),
413 fNewSampleNotifier(fNSlots),
414 fSampleInfos(fNSlots),
415 fDatasetColumnReaders(fNSlots)
416{
417}
418
419RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
420 : fDefaultColumns(defaultBranches),
421 fNSlots(RDFInternal::GetNSlots()),
422 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
423 fDataSource(std::move(ds)),
424 fNewSampleNotifier(fNSlots),
425 fSampleInfos(fNSlots),
426 fDatasetColumnReaders(fNSlots)
427{
428 fDataSource->SetNSlots(fNSlots);
429}
430
432 : fNSlots(RDFInternal::GetNSlots()),
433 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
434 fNewSampleNotifier(fNSlots),
435 fSampleInfos(fNSlots),
436 fDatasetColumnReaders(fNSlots)
437{
438 ChangeSpec(std::move(spec));
439}
440
441#ifdef R__UNIX
442namespace {
443std::optional<std::string> GetRedirectedSampleId(std::string_view path, std::string_view datasetName)
444{
445 // Mimick the redirection done in TFile::Open to see if the path points to a FUSE-mounted EOS path.
446 // If so, we create a redirected sample ID with the full xroot URL.
447 TString expandedUrl(path.data());
449 if (gEnv->GetValue("TFile.CrossProtocolRedirects", 1) == 1) {
450 TUrl fileurl(expandedUrl, /* default is file */ kTRUE);
451 if (strcmp(fileurl.GetProtocol(), "file") == 0) {
452 ssize_t len = getxattr(fileurl.GetFile(), "eos.url.xroot", nullptr, 0);
453 if (len > 0) {
454 std::string xurl(len, 0);
455 std::string fileNameFromUrl{fileurl.GetFile()};
456 if (getxattr(fileNameFromUrl.c_str(), "eos.url.xroot", &xurl[0], len) == len) {
457 // Sometimes the `getxattr` call may return an invalid URL due
458 // to the POSIX attribute not being yet completely filled by EOS.
459 if (auto baseName = fileNameFromUrl.substr(fileNameFromUrl.find_last_of("/") + 1);
460 std::equal(baseName.crbegin(), baseName.crend(), xurl.crbegin())) {
461 return xurl + '/' + datasetName.data();
462 }
463 }
464 }
465 }
466 }
467
468 return std::nullopt;
469}
470} // namespace
471#endif
472
473/**
474 * @brief Changes the internal TTree held by the RLoopManager.
475 *
476 * @warning This method may lead to potentially dangerous interactions if used
477 * after the construction of the RDataFrame. Changing the specification makes
478 * sense *if and only if* the schema of the dataset is *unchanged*, i.e. the
479 * new specification refers to exactly the same number of columns, with the
480 * same names and types. The actual use case of this method is moving the
481 * processing of the same RDataFrame to a different range of entries of the
482 * same dataset (which may be stored in a different set of files).
483 *
484 * @param spec The specification of the dataset to be adopted.
485 */
487{
488 // Change the range of entries to be processed
489 fBeginEntry = spec.GetEntryRangeBegin();
490 fEndEntry = spec.GetEntryRangeEnd();
491
492 // Store the samples
493 fSamples = spec.MoveOutSamples();
494 fSampleMap.clear();
495
496 // Create the internal main chain
498 for (auto &sample : fSamples) {
499 const auto &trees = sample.GetTreeNames();
500 const auto &files = sample.GetFileNameGlobs();
501 for (std::size_t i = 0ul; i < files.size(); ++i) {
502 // We need to use `<filename>?#<treename>` as an argument to TChain::Add
503 // (see https://github.com/root-project/root/pull/8820 for why)
504 const auto fullpath = files[i] + "?#" + trees[i];
505 chain->Add(fullpath.c_str());
506 // ...but instead we use `<filename>/<treename>` as a sample ID (cannot
507 // change this easily because of backward compatibility: the sample ID
508 // is exposed to users via RSampleInfo and DefinePerSample).
509 const auto sampleId = files[i] + '/' + trees[i];
510 fSampleMap.insert({sampleId, &sample});
511#ifdef R__UNIX
512 // Also add redirected EOS xroot URL when available
514 fSampleMap.insert({redirectedSampleId.value(), &sample});
515#endif
516 }
517 }
518 fDataSource = std::make_unique<ROOT::Internal::RDF::RTTreeDS>(std::move(chain), spec.GetFriendInfo());
519 fDataSource->SetNSlots(fNSlots);
520 for (unsigned int slot{}; slot < fNSlots; slot++) {
521 for (auto &v : fDatasetColumnReaders[slot])
522 v.second.reset();
523 }
524}
525
526/// Run event loop with no source files, in parallel.
528{
529#ifdef R__USE_IMT
530 std::shared_ptr<ROOT::Internal::RSlotStack> slotStack = SlotStack();
531 // Working with an empty tree.
532 // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
533 const auto nEmptyEntries = GetNEmptyEntries();
534 const auto nEntriesPerSlot = nEmptyEntries / (fNSlots * 2);
535 auto remainder = nEmptyEntries % (fNSlots * 2);
536 std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
537 ULong64_t begin = fEmptyEntryRange.first;
538 while (begin < fEmptyEntryRange.second) {
539 ULong64_t end = begin + nEntriesPerSlot;
540 if (remainder > 0) {
541 ++end;
542 --remainder;
543 }
544 entryRanges.emplace_back(begin, end);
545 begin = end;
546 }
547
548 // Each task will generate a subrange of entries
549 auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
551 auto slot = slotRAII.fSlot;
552 RCallCleanUpTask cleanup(*this, slot);
553 InitNodeSlots(nullptr, slot);
554 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({"an empty source", range.first, range.second, slot});
555 try {
557 for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
559 }
560 } catch (...) {
561 // Error might throw in experiment frameworks like CMSSW
562 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
563 throw;
564 }
565 };
566
568 pool.Foreach(genFunction, entryRanges);
569
570#endif // not implemented otherwise
571}
572
573/// Run event loop with no source files, in sequence.
575{
576 InitNodeSlots(nullptr, 0);
578 {"an empty source", fEmptyEntryRange.first, fEmptyEntryRange.second, 0u});
579 RCallCleanUpTask cleanup(*this);
580 try {
585 }
586 } catch (...) {
587 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
588 throw;
589 }
590}
591
592namespace {
593/// Return true on succesful entry read.
594///
595/// TTreeReader encodes successful reads in the `kEntryValid` enum value, but
596/// there can be other situations where the read is still valid. For now, these
597/// are:
598/// - If there was no match of the current entry in one or more friend trees
599/// according to their respective indexes.
600/// - If there was a missing branch at the start of a new tree in the dataset.
601///
602/// In such situations, although the entry is not complete, the processing
603/// should not be aborted and nodes of the computation graph will take action
604/// accordingly.
606{
607 treeReader.Next();
608 switch (treeReader.GetEntryStatus()) {
609 case TTreeReader::kEntryValid: return true;
610 case TTreeReader::kIndexedFriendNoMatch: return true;
612 default: return false;
613 }
614}
615} // namespace
616
617/// Run event loop over one or multiple ROOT files, in parallel.
619{
620#ifdef R__USE_IMT
621 if (fEndEntry == fBeginEntry) // empty range => no work needed
622 return;
623 std::shared_ptr<ROOT::Internal::RSlotStack> slotStack = SlotStack();
624 const auto &entryList = fTree->GetEntryList() ? *fTree->GetEntryList() : TEntryList();
625 auto tp =
626 (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
627 ? std::make_unique<ROOT::TTreeProcessorMT>(*fTree, fNSlots, std::make_pair(fBeginEntry, fEndEntry),
629 : std::make_unique<ROOT::TTreeProcessorMT>(*fTree, entryList, fNSlots, fSuppressErrorsForMissingBranches);
630
631 std::atomic<ULong64_t> entryCount(0ull);
632
633 tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
635 auto slot = slotRAII.fSlot;
636 RCallCleanUpTask cleanup(*this, slot, &r);
639 const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
640 const auto nEntries = entryRange.second - entryRange.first;
641 auto count = entryCount.fetch_add(nEntries);
642 try {
643 // recursive call to check filters and conditionally execute actions
644 while (validTTreeReaderRead(r)) {
647 }
648 RunAndCheckFilters(slot, count++);
649 }
650 } catch (...) {
651 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
652 throw;
653 }
654 // fNStopsReceived < fNChildren is always true at the moment as we don't support event loop early quitting in
655 // multi-thread runs, but it costs nothing to be safe and future-proof in case we add support for that later.
656 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
657 // something went wrong in the TTreeReader event loop
658 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
659 std::to_string(r.GetEntryStatus()));
660 }
661 });
662
663 auto &&processedEntries = entryCount.load();
664 if (fEndEntry != std::numeric_limits<Long64_t>::max() &&
666 Warning("RDataFrame::Run",
667 "RDataFrame stopped processing after %lld entries, whereas an entry range (begin=%lld,end=%lld) was "
668 "requested. Consider adjusting the end value of the entry range to a maximum of %lld.",
670 }
671#endif // no-op otherwise (will not be called)
672}
673
674/// Run event loop over one or multiple ROOT files, in sequence.
676{
677 TTreeReader r(fTree.get(), fTree->GetEntryList(), /*warnAboutLongerFriends*/ true,
679 if (0 == fTree->GetEntriesFast() || fBeginEntry == fEndEntry)
680 return;
681 // Apply the range if there is any
682 // In case of a chain with a total of N entries, calling SetEntriesRange(N + 1, ...) does not error out
683 // This is a bug, reported here: https://github.com/root-project/root/issues/10774
684 if (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
685 if (r.SetEntriesRange(fBeginEntry, fEndEntry) != TTreeReader::kEntryValid)
686 throw std::logic_error("Something went wrong in initializing the TTreeReader.");
687
688 RCallCleanUpTask cleanup(*this, 0u, &r);
689 InitNodeSlots(&r, 0);
691
692 // recursive call to check filters and conditionally execute actions
693 // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
695 try {
698 UpdateSampleInfo(/*slot*/0, r);
699 }
700 RunAndCheckFilters(0, r.GetCurrentEntry());
702 }
703 } catch (...) {
704 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
705 throw;
706 }
707 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
708 // something went wrong in the TTreeReader event loop
709 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
710 std::to_string(r.GetEntryStatus()));
711 }
712
713 if (fEndEntry != std::numeric_limits<Long64_t>::max() && (fEndEntry - fBeginEntry) > processedEntries) {
714 Warning("RDataFrame::Run",
715 "RDataFrame stopped processing after %lld entries, whereas an entry range (begin=%lld,end=%lld) was "
716 "requested. Consider adjusting the end value of the entry range to a maximum of %lld.",
718 }
719}
720
721namespace {
722struct DSRunRAII {
724 DSRunRAII(ROOT::RDF::RDataSource &ds, const std::set<std::string> &suppressErrorsForMissingColumns) : fDS(ds)
725 {
727 }
728 ~DSRunRAII() { fDS.Finalize(); }
729};
730} // namespace
731
734 unsigned int fSlot;
737 TTreeReader *treeReader = nullptr)
738 : fLM(lm), fSlot(slot), fTreeReader(treeReader)
739 {
740 fLM.InitNodeSlots(fTreeReader, fSlot);
741 fLM.GetDataSource()->InitSlot(fSlot, firstEntry);
742 }
744};
745
746/// Run event loop over data accessed through a DataSource, in sequence.
748{
749 assert(fDataSource != nullptr);
750 // Shortcut if the entry range would result in not reading anything
751 if (fBeginEntry == fEndEntry)
752 return;
753 // Apply global entry range if necessary
754 if (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
755 fDataSource->SetGlobalEntryRange(std::make_pair<std::uint64_t, std::uint64_t>(fBeginEntry, fEndEntry));
756 // Initialize data source and book finalization
758 // Ensure cleanup task is always called at the end. Notably, this also resets the column readers for those data
759 // sources that need it (currently only TTree).
760 RCallCleanUpTask cleanup(*this);
761
762 // Main event loop. We start with an empty vector of ranges because we need to initialize the nodes and the data
763 // source before the first call to GetEntryRanges, since it could trigger reading (currently only happens with
764 // TTree).
765 std::uint64_t processedEntries{};
766 std::vector<std::pair<ULong64_t, ULong64_t>> ranges{};
767 do {
768
770
771 ranges = fDataSource->GetEntryRanges();
772
774
775 try {
776 for (const auto &range : ranges) {
777 const auto start = range.first;
778 const auto end = range.second;
779 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, 0u});
780 for (auto entry = start; entry < end && fNStopsReceived < fNChildren; ++entry) {
781 if (fDataSource->SetEntry(0u, entry)) {
783 }
785 }
786 }
787 } catch (...) {
788 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
789 throw;
790 }
791
792 } while (!ranges.empty() && fNStopsReceived < fNChildren);
793
795
796 if (fEndEntry != std::numeric_limits<Long64_t>::max() &&
797 static_cast<std::uint64_t>(fEndEntry - fBeginEntry) > processedEntries) {
798 std::ostringstream buf{};
799 buf << "RDataFrame stopped processing after ";
800 buf << processedEntries;
801 buf << " entries, whereas an entry range (begin=";
802 buf << fBeginEntry;
803 buf << ",end=";
804 buf << fEndEntry;
805 buf << ") was requested. Consider adjusting the end value of the entry range to a maximum of ";
806 buf << (fBeginEntry + processedEntries);
807 buf << ".";
808 Warning("RDataFrame::Run", "%s", buf.str().c_str());
809 }
810}
811
812/// Run event loop over data accessed through a DataSource, in parallel.
814{
815#ifdef R__USE_IMT
816 assert(fDataSource != nullptr);
817 // Shortcut if the entry range would result in not reading anything
818 if (fBeginEntry == fEndEntry)
819 return;
820 // Apply global entry range if necessary
821 if (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
822 fDataSource->SetGlobalEntryRange(std::make_pair<std::uint64_t, std::uint64_t>(fBeginEntry, fEndEntry));
823
825
827
828#endif // not implemented otherwise (never called)
829}
830
831/// Execute actions and make sure named filters are called for each event.
832/// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
834{
835 // data-block callbacks run before the rest of the graph
837 for (auto &callback : fSampleCallbacks)
838 callback.second(slot, fSampleInfos[slot]);
840 }
841
842 for (auto *actionPtr : fBookedActions)
843 actionPtr->Run(slot, entry);
845 namedFilterPtr->CheckFilters(slot, entry);
846 for (auto &callback : fCallbacksEveryNEvents)
847 callback(slot);
848}
849
850/// Build TTreeReaderValues for all nodes
851/// This method loops over all filters, actions and other booked objects and
852/// calls their `InitSlot` method, to get them ready for running a task.
854{
856 for (auto *ptr : fBookedActions)
857 ptr->InitSlot(r, slot);
858 for (auto *ptr : fBookedFilters)
859 ptr->InitSlot(r, slot);
860 for (auto *ptr : fBookedDefines)
861 ptr->InitSlot(r, slot);
862 for (auto *ptr : fBookedVariations)
863 ptr->InitSlot(r, slot);
864
865 for (auto &callback : fCallbacksOnce)
866 callback(slot);
867}
868
870 if (r != nullptr) {
871 // we need to set a notifier so that we run the callbacks every time we switch to a new TTree
872 // `PrependLink` inserts this notifier into the TTree/TChain's linked list of notifiers
873 fNewSampleNotifier.GetChainNotifyLink(slot).PrependLink(*r->GetTree());
874 }
875 // Whatever the data source, initially set the "new data block" flag:
876 // - for TChains, this ensures that we don't skip the first data block because
877 // the correct tree is already loaded
878 // - for RDataSources and empty sources, which currently don't have data blocks, this
879 // ensures that we run once per task
881}
882
883void RLoopManager::UpdateSampleInfo(unsigned int slot, const std::pair<ULong64_t, ULong64_t> &range) {
885 "Empty source, range: {" + std::to_string(range.first) + ", " + std::to_string(range.second) + "}", range);
886}
887
889 // one GetTree to retrieve the TChain, another to retrieve the underlying TTree
890 auto *tree = r.GetTree()->GetTree();
891 R__ASSERT(tree != nullptr);
892 const std::string treename = ROOT::Internal::TreeUtils::GetTreeFullPaths(*tree)[0];
893 auto *file = tree->GetCurrentFile();
894 const std::string fname = file != nullptr ? file->GetName() : "#inmemorytree#";
895
896 std::pair<Long64_t, Long64_t> range = r.GetEntriesRange();
897 R__ASSERT(range.first >= 0);
898 if (range.second == -1) {
899 range.second = tree->GetEntries(); // convert '-1', i.e. 'until the end', to the actual entry number
900 }
901 // If the tree is stored in a subdirectory, treename will be the full path to it starting with the root directory '/'
902 const std::string &id = fname + (treename.rfind('/', 0) == 0 ? "" : "/") + treename;
903 if (fSampleMap.empty()) {
905 } else {
906 if (fSampleMap.find(id) == fSampleMap.end())
907 throw std::runtime_error("Full sample identifier '" + id + "' cannot be found in the available samples.");
909 }
910}
911
912/// Create a slot stack with the desired number of slots or reuse a shared instance.
913/// When a LoopManager runs in isolation, it will create its own slot stack from the
914/// number of slots. When it runs as part of RunGraphs(), each loop manager will be
915/// assigned a shared slot stack, so dataframe helpers can be shared in a thread-safe
916/// manner.
917std::shared_ptr<ROOT::Internal::RSlotStack> RLoopManager::SlotStack() const
918{
919#ifdef R__USE_IMT
920 if (auto shared = fSlotStack.lock(); shared) {
921 return shared;
922 } else {
923 return std::make_shared<ROOT::Internal::RSlotStack>(fNSlots);
924 }
925#else
926 return nullptr;
927#endif
928}
929
930/// Initialize all nodes of the functional graph before running the event loop.
931/// This method is called once per event-loop and performs generic initialization
932/// operations that do not depend on the specific processing slot (i.e. operations
933/// that are common for all threads).
935{
937 for (auto *filter : fBookedFilters)
938 filter->InitNode();
939 for (auto *range : fBookedRanges)
940 range->InitNode();
941 for (auto *ptr : fBookedActions)
942 ptr->Initialize();
943}
944
945/// Perform clean-up operations. To be called at the end of each event loop.
947{
948 fMustRunNamedFilters = false;
949
950 // forget RActions and detach TResultProxies
951 for (auto *ptr : fBookedActions)
952 ptr->Finalize();
953
954 fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
955 fBookedActions.clear();
956
957 // reset children counts
958 fNChildren = 0;
959 fNStopsReceived = 0;
960 for (auto *ptr : fBookedFilters)
961 ptr->ResetChildrenCount();
962 for (auto *ptr : fBookedRanges)
963 ptr->ResetChildrenCount();
964
966 fCallbacksOnce.clear();
967}
968
969/// Perform clean-up operations. To be called at the end of each task execution.
971{
972 if (r != nullptr)
973 fNewSampleNotifier.GetChainNotifyLink(slot).RemoveLink(*r->GetTree());
974 for (auto *ptr : fBookedActions)
975 ptr->FinalizeSlot(slot);
976 for (auto *ptr : fBookedFilters)
977 ptr->FinalizeSlot(slot);
978 for (auto *ptr : fBookedDefines)
979 ptr->FinalizeSlot(slot);
980
981 if (auto ds = GetDataSource(); ds && ds->GetLabel() == "TTreeDS") {
982 // we are reading from a tree/chain and we need to re-create the RTreeColumnReaders at every task
983 // because the TTreeReader object changes at every task
984 for (auto &v : fDatasetColumnReaders[slot])
985 v.second.reset();
986 }
987}
988
989/// Add RDF nodes that require just-in-time compilation to the computation graph.
990/// This method also clears the contents of GetCodeToJit().
992{
993 {
995 if (GetCodeToJit().empty()) {
996 R__LOG_INFO(RDFLogChannel()) << "Nothing to jit and execute.";
997 return;
998 }
999 }
1000
1001 const std::string code = []() {
1003 return std::move(GetCodeToJit());
1004 }();
1005
1006 TStopwatch s;
1007 s.Start();
1008 RDFInternal::InterpreterCalc(code, "RLoopManager::Run");
1009 s.Stop();
1010 R__LOG_INFO(RDFLogChannel()) << "Just-in-time compilation phase completed"
1011 << (s.RealTime() > 1e-3 ? " in " + std::to_string(s.RealTime()) + " seconds."
1012 : " in less than 1ms.");
1013}
1014
1015/// Trigger counting of number of children nodes for each node of the functional graph.
1016/// This is done once before starting the event loop. Each action sends an `increase children count` signal
1017/// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
1018/// children counter. Each node only propagates the signal once, even if it receives it multiple times.
1019/// Named filters also send an `increase children count` signal, just like actions, as they always execute during
1020/// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
1022{
1023 for (auto *actionPtr : fBookedActions)
1024 actionPtr->TriggerChildrenCount();
1026 namedFilterPtr->TriggerChildrenCount();
1027}
1028
1029/// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
1030/// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
1031/// The jitting phase is skipped if the `jit` parameter is `false` (unsafe, use with care).
1033{
1034 // Change value of TTree::GetMaxTreeSize only for this scope. Revert when #6640 will be solved.
1035 MaxTreeSizeRAII ctxtmts;
1036
1037 R__LOG_INFO(RDFLogChannel()) << "Starting event loop number " << fNRuns << '.';
1038
1040
1041 if (jit)
1042 Jit();
1043
1044 InitNodes();
1045
1046 // Exceptions can occur during the event loop. In order to ensure proper cleanup of nodes
1047 // we use RAII: even in case of an exception, the destructor of the object is invoked and
1048 // all the cleanup takes place.
1049 class NodesCleanerRAII {
1051
1052 public:
1054 ~NodesCleanerRAII() { fRLM.CleanUpNodes(); }
1055 };
1056
1058
1059 TStopwatch s;
1060 s.Start();
1061
1062 switch (fLoopType) {
1064 throw std::runtime_error("RDataFrame: executing the computation graph without a data source, aborting.");
1065 break;
1069 case ELoopType::kNoFiles: RunEmptySource(); break;
1070 case ELoopType::kROOTFiles: RunTreeReader(); break;
1072 }
1073 s.Stop();
1074
1075 fNRuns++;
1076
1077 R__LOG_INFO(RDFLogChannel()) << "Finished event loop number " << fNRuns - 1 << " (" << s.CpuTime() << "s CPU, "
1078 << s.RealTime() << "s elapsed).";
1079}
1080
1081/// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
1086
1088{
1089 // This is currently called in SnapshotTTreeHelper[MT] to retrieve the task-local input TTree. It is not guaranteed
1090 // that the same RLoopManager will always have the same input TTree for its entire lifetime, notably it could be
1091 // changed by ChangeSpec when moving to a different entry range.
1092 if (auto *treeDS = dynamic_cast<ROOT::Internal::RDF::RTTreeDS *>(fDataSource.get())) {
1093 return treeDS->GetTree();
1094 }
1095 return nullptr;
1096}
1097
1103
1110
1112{
1113 fBookedFilters.emplace_back(filterPtr);
1114 if (filterPtr->HasName()) {
1115 fBookedNamedFilters.emplace_back(filterPtr);
1116 fMustRunNamedFilters = true;
1117 }
1118}
1119
1125
1127{
1128 fBookedRanges.emplace_back(rangePtr);
1129}
1130
1135
1137{
1138 fBookedDefines.emplace_back(ptr);
1139}
1140
1146
1151
1156
1157// dummy call, end of recursive chain of calls
1159{
1160 return true;
1161}
1162
1163/// Call `FillReport` on all booked filters
1165{
1166 for (const auto *fPtr : fBookedNamedFilters)
1167 fPtr->FillReport(rep);
1168}
1169
1170void RLoopManager::SetTree(std::shared_ptr<TTree> tree)
1171{
1172 fTree = std::move(tree);
1173
1174 TChain *ch = nullptr;
1175 if ((ch = dynamic_cast<TChain *>(fTree.get())))
1177}
1178
1179void RLoopManager::ToJitExec(const std::string &code) const
1180{
1182 GetCodeToJit().append(code);
1183}
1184
1185void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
1186{
1187 if (everyNEvents == 0ull)
1188 fCallbacksOnce.emplace_back(std::move(f), fNSlots);
1189 else
1190 fCallbacksEveryNEvents.emplace_back(everyNEvents, std::move(f), fNSlots);
1191}
1192
1193std::vector<std::string> RLoopManager::GetFiltersNames()
1194{
1195 std::vector<std::string> filters;
1196 for (auto *filter : fBookedFilters) {
1197 auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
1198 filters.push_back(name);
1199 }
1200 return filters;
1201}
1202
1203std::vector<RNodeBase *> RLoopManager::GetGraphEdges() const
1204{
1205 std::vector<RNodeBase *> nodes(fBookedFilters.size() + fBookedRanges.size());
1206 auto it = std::copy(fBookedFilters.begin(), fBookedFilters.end(), nodes.begin());
1207 std::copy(fBookedRanges.begin(), fBookedRanges.end(), it);
1208 return nodes;
1209}
1210
1211std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions() const
1212{
1213 std::vector<RDFInternal::RActionBase *> actions(fBookedActions.size() + fRunActions.size());
1214 auto it = std::copy(fBookedActions.begin(), fBookedActions.end(), actions.begin());
1215 std::copy(fRunActions.begin(), fRunActions.end(), it);
1216 return actions;
1217}
1218
1219std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph(
1220 std::unordered_map<void *, std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode>> &visitedMap)
1221{
1222 // If there is already a node for the RLoopManager return it. If there is not, return a new one.
1223 auto duplicateRLoopManagerIt = visitedMap.find((void *)this);
1225 return duplicateRLoopManagerIt->second;
1226
1227 std::string name;
1228 if (fDataSource) {
1229 name = fDataSource->GetLabel();
1230 } else if (fTree) {
1231 name = fTree->GetName();
1232 if (name.empty())
1233 name = fTree->ClassName();
1234 } else {
1235 name = "Empty source\\nEntries: " + std::to_string(GetNEmptyEntries());
1236 }
1237 auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(
1239 visitedMap[(void *)this] = thisNode;
1240 return thisNode;
1241}
1242
1243////////////////////////////////////////////////////////////////////////////
1244/// Return all valid TTree::Branch names (caching results for subsequent calls).
1245/// Never use fBranchNames directy, always request it through this method.
1247{
1248 if (fValidBranchNames.empty() && fTree) {
1249 fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
1250 }
1251 if (fValidBranchNames.empty() && fDataSource) {
1252 fValidBranchNames = fDataSource->GetColumnNames();
1253 }
1254 return fValidBranchNames;
1255}
1256
1257/// Return true if AddDataSourceColumnReaders was called for column name col.
1258bool RLoopManager::HasDataSourceColumnReaders(std::string_view col, const std::type_info &ti) const
1259{
1260 const auto key = MakeDatasetColReadersKey(col, ti);
1261 assert(fDataSource != nullptr);
1262 // since data source column readers are always added for all slots at the same time,
1263 // if the reader is present for slot 0 we have it for all other slots as well.
1264 auto it = fDatasetColumnReaders[0].find(key);
1265 return (it != fDatasetColumnReaders[0].end() && it->second);
1266}
1267
1269 std::vector<std::unique_ptr<RColumnReaderBase>> &&readers,
1270 const std::type_info &ti)
1271{
1272 const auto key = MakeDatasetColReadersKey(col, ti);
1273 assert(fDataSource != nullptr && !HasDataSourceColumnReaders(col, ti));
1274 assert(readers.size() == fNSlots);
1275
1276 for (auto slot = 0u; slot < fNSlots; ++slot) {
1277 fDatasetColumnReaders[slot][key] = std::move(readers[slot]);
1278 }
1279}
1280
1281// Differently from AddDataSourceColumnReaders, this can be called from multiple threads concurrently
1282/// \brief Register a new RTreeColumnReader with this RLoopManager.
1283/// \return A shared pointer to the inserted column reader.
1285 std::unique_ptr<RColumnReaderBase> &&reader,
1286 const std::type_info &ti)
1287{
1289 const auto key = MakeDatasetColReadersKey(col, ti);
1290 // if a reader for this column and this slot was already there, we are doing something wrong
1291 assert(readers.find(key) == readers.end() || readers[key] == nullptr);
1292 auto *rptr = reader.get();
1293 readers[key] = std::move(reader);
1294 return rptr;
1295}
1296
1298 const std::type_info &ti, TTreeReader *treeReader)
1299{
1301 const auto key = MakeDatasetColReadersKey(col, ti);
1302 // if a reader for this column and this slot was already there, we are doing something wrong
1303 assert(readers.find(key) == readers.end() || readers[key] == nullptr);
1304 assert(fDataSource && "Missing RDataSource to add column reader.");
1305
1307
1308 return readers[key].get();
1309}
1310
1312RLoopManager::GetDatasetColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti) const
1313{
1314 const auto key = MakeDatasetColReadersKey(col, ti);
1315 if (auto it = fDatasetColumnReaders[slot].find(key); it != fDatasetColumnReaders[slot].end() && it->second)
1316 return it->second.get();
1317 else
1318 return nullptr;
1319}
1320
1322{
1323 if (callback)
1324 fSampleCallbacks.insert({nodePtr, std::move(callback)});
1325}
1326
1327void RLoopManager::SetEmptyEntryRange(std::pair<ULong64_t, ULong64_t> &&newRange)
1328{
1329 fEmptyEntryRange = std::move(newRange);
1330}
1331
1333{
1334 fBeginEntry = begin;
1335 fEndEntry = end;
1336}
1337
1339{
1340 fTTreeLifeline = std::move(lifeline);
1341}
1342
1343/**
1344 * \brief Helper function to open a file (or the first file from a glob).
1345 * This function is used at construction time of an RDataFrame, to check the
1346 * concrete type of the dataset stored inside the file.
1347 */
1348std::unique_ptr<TFile> OpenFileWithSanityChecks(std::string_view fileNameGlob)
1349{
1350 // Follow same logic in TChain::Add to find the correct string to look for globbing:
1351 // - If the extension ".root" is present in the file name, pass along the basename.
1352 // - If not, use the "?" token to delimit the part of the string which represents the basename.
1353 // - Otherwise, pass the full filename.
1354 auto &&baseNameAndQuery = [&fileNameGlob]() {
1355 constexpr std::string_view delim{".root"};
1356 if (auto &&it = std::find_end(fileNameGlob.begin(), fileNameGlob.end(), delim.begin(), delim.end());
1357 it != fileNameGlob.end()) {
1358 auto &&distanceToEndOfDelim = std::distance(fileNameGlob.begin(), it + delim.length());
1359 return std::make_pair(fileNameGlob.substr(0, distanceToEndOfDelim), fileNameGlob.substr(distanceToEndOfDelim));
1360 } else if (auto &&lastQuestionMark = fileNameGlob.find_last_of('?'); lastQuestionMark != std::string_view::npos)
1361 return std::make_pair(fileNameGlob.substr(0, lastQuestionMark), fileNameGlob.substr(lastQuestionMark));
1362 else
1363 return std::make_pair(fileNameGlob, std::string_view{});
1364 }();
1365 // Captured structured bindings variable are only valid since C++20
1366 auto &&baseName = baseNameAndQuery.first;
1367 auto &&query = baseNameAndQuery.second;
1368
1369 const auto nameHasWildcard = [&baseName]() {
1370 constexpr std::array<char, 4> wildCards{'[', ']', '*', '?'}; // Wildcards accepted by TChain::Add
1371 return std::any_of(wildCards.begin(), wildCards.end(),
1372 [&baseName](auto &&wc) { return baseName.find(wc) != std::string_view::npos; });
1373 }();
1374
1375 // Open first file in case of glob, suppose all files in the glob use the same data format
1376 std::string fileToOpen{nameHasWildcard
1377 ? ROOT::Internal::TreeUtils::ExpandGlob(std::string{baseName})[0] + std::string{query}
1378 : fileNameGlob};
1379
1380 ::TDirectory::TContext ctxt; // Avoid changing gDirectory;
1381 std::unique_ptr<TFile> inFile{TFile::Open(fileToOpen.c_str(), "READ_WITHOUT_GLOBALREGISTRATION")};
1382 if (!inFile || inFile->IsZombie())
1383 throw std::invalid_argument("RDataFrame: could not open file \"" + fileToOpen + "\".");
1384
1385 return inFile;
1386}
1387
1388std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1391{
1392 // Introduce the same behaviour as in CreateLMFromFile for consistency.
1393 // Creating an RDataFrame with a non-existing file will throw early rather
1394 // than wait for the start of the graph execution.
1395 if (checkFile) {
1397 }
1398
1399 auto dataSource = std::make_unique<ROOT::Internal::RDF::RTTreeDS>(datasetName, fileNameGlob);
1400 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(dataSource), defaultColumns);
1401 return lm;
1402}
1403
1404std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1405ROOT::Detail::RDF::CreateLMFromTTree(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1406 const std::vector<std::string> &defaultColumns, bool checkFile)
1407{
1408 if (fileNameGlobs.size() == 0)
1409 throw std::invalid_argument("RDataFrame: empty list of input files.");
1410 // Introduce the same behaviour as in CreateLMFromFile for consistency.
1411 // Creating an RDataFrame with a non-existing file will throw early rather
1412 // than wait for the start of the graph execution.
1413 if (checkFile) {
1415 }
1416 auto dataSource = std::make_unique<ROOT::Internal::RDF::RTTreeDS>(datasetName, fileNameGlobs);
1417 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(dataSource), defaultColumns);
1418 return lm;
1419}
1420
1421std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1424{
1425 auto dataSource = std::make_unique<ROOT::RDF::RNTupleDS>(datasetName, fileNameGlob);
1426 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(dataSource), defaultColumns);
1427 return lm;
1428}
1429
1430std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1431ROOT::Detail::RDF::CreateLMFromRNTuple(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1433{
1434 auto dataSource = std::make_unique<ROOT::RDF::RNTupleDS>(datasetName, fileNameGlobs);
1435 auto lm = std::make_shared<ROOT::Detail::RDF::RLoopManager>(std::move(dataSource), defaultColumns);
1436 return lm;
1437}
1438
1439std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1442{
1443
1445
1446 if (inFile->Get<TTree>(datasetName.data())) {
1447 return CreateLMFromTTree(datasetName, fileNameGlob, defaultColumns, /*checkFile=*/false);
1448 } else if (inFile->Get<ROOT::RNTuple>(datasetName.data())) {
1450 }
1451
1452 throw std::invalid_argument("RDataFrame: unsupported data format for dataset \"" + std::string(datasetName) +
1453 "\" in file \"" + inFile->GetName() + "\".");
1454}
1455
1456std::shared_ptr<ROOT::Detail::RDF::RLoopManager>
1457ROOT::Detail::RDF::CreateLMFromFile(std::string_view datasetName, const std::vector<std::string> &fileNameGlobs,
1459{
1460
1461 if (fileNameGlobs.size() == 0)
1462 throw std::invalid_argument("RDataFrame: empty list of input files.");
1463
1465
1466 if (inFile->Get<TTree>(datasetName.data())) {
1467 return CreateLMFromTTree(datasetName, fileNameGlobs, defaultColumns, /*checkFile=*/false);
1468 } else if (inFile->Get<ROOT::RNTuple>(datasetName.data())) {
1470 }
1471
1472 throw std::invalid_argument("RDataFrame: unsupported data format for dataset \"" + std::string(datasetName) +
1473 "\" in file \"" + inFile->GetName() + "\".");
1474}
1475
1476// outlined to pin virtual table
1478
1479void ROOT::Detail::RDF::RLoopManager::SetDataSource(std::unique_ptr<ROOT::RDF::RDataSource> dataSource)
1480{
1481 if (dataSource) {
1482 fDataSource = std::move(dataSource);
1483 fDataSource->SetNSlots(fNSlots);
1484 fLoopType = ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource;
1485 }
1486}
1487
1488void ROOT::Detail::RDF::RLoopManager::DataSourceThreadTask(const std::pair<ULong64_t, ULong64_t> &entryRange,
1490 std::atomic<ULong64_t> &entryCount)
1491{
1492#ifdef R__USE_IMT
1494 const auto &slot = slotRAII.fSlot;
1495
1496 const auto &[start, end] = entryRange;
1497 const auto nEntries = end - start;
1498 entryCount.fetch_add(nEntries);
1499
1500 RCallCleanUpTask cleanup(*this, slot);
1501 RDSRangeRAII _{*this, slot, start};
1502
1503 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, slot});
1504
1505 try {
1506 for (auto entry = start; entry < end; ++entry) {
1507 if (fDataSource->SetEntry(slot, entry)) {
1508 RunAndCheckFilters(slot, entry);
1509 }
1510 }
1511 } catch (...) {
1512 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
1513 throw;
1514 }
1515 fDataSource->FinalizeSlot(slot);
1516#else
1517 (void)entryRange;
1518 (void)slotStack;
1519 (void)entryCount;
1520#endif
1521}
1522
1524 std::atomic<ULong64_t> &entryCount)
1525{
1526#ifdef R__USE_IMT
1528 const auto &slot = slotRAII.fSlot;
1529
1530 const auto entryRange = treeReader.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
1531 const auto &[start, end] = entryRange;
1532 const auto nEntries = end - start;
1533 auto count = entryCount.fetch_add(nEntries);
1534
1535 RDSRangeRAII _{*this, slot, static_cast<ULong64_t>(start), &treeReader};
1536 RCallCleanUpTask cleanup(*this, slot, &treeReader);
1537
1539 {fDataSource->GetLabel(), static_cast<ULong64_t>(start), static_cast<ULong64_t>(end), slot});
1540 try {
1541 // recursive call to check filters and conditionally execute actions
1543 if (fNewSampleNotifier.CheckFlag(slot)) {
1544 UpdateSampleInfo(slot, treeReader);
1545 }
1546 RunAndCheckFilters(slot, count++);
1547 }
1548 } catch (...) {
1549 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
1550 throw;
1551 }
1552 // fNStopsReceived < fNChildren is always true at the moment as we don't support event loop early quitting in
1553 // multi-thread runs, but it costs nothing to be safe and future-proof in case we add support for that later.
1554 if (treeReader.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
1555 // something went wrong in the TTreeReader event loop
1556 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
1557 std::to_string(treeReader.GetEntryStatus()));
1558 }
1559#else
1560 (void)treeReader;
1561 (void)slotStack;
1562 (void)entryCount;
1563#endif
1564}
#define R__LOG_DEBUG(DEBUGLEVEL,...)
Definition RLogger.hxx:360
#define R__LOG_INFO(...)
Definition RLogger.hxx:359
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
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
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
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
char name[80]
Definition TGX11.cxx:110
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
#define R__WRITE_LOCKGUARD(mutex)
#define R__READ_LOCKGUARD(mutex)
#define _(A, B)
Definition cfortran.h:108
The head node of a RDF computation graph.
RColumnReaderBase * AddDataSourceColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti, TTreeReader *treeReader)
void UpdateSampleInfo(unsigned int slot, const std::pair< ULong64_t, ULong64_t > &range)
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.
void AddDataSourceColumnReaders(std::string_view col, std::vector< std::unique_ptr< RColumnReaderBase > > &&readers, const std::type_info &ti)
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
std::set< std::string > fSuppressErrorsForMissingBranches
void ChangeSpec(ROOT::RDF::Experimental::RDatasetSpec &&spec)
Changes the internal TTree held by the RLoopManager.
void SetTree(std::shared_ptr< TTree > tree)
std::weak_ptr< ROOT::Internal::RSlotStack > fSlotStack
Pointer to a shared slot stack in case this instance runs concurrently with others:
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.
void TTreeThreadTask(TTreeReader &treeReader, ROOT::Internal::RSlotStack &slotStack, std::atomic< ULong64_t > &entryCount)
The task run by every thread on an entry range (known by the input TTreeReader), for the TTree data s...
ROOT::Internal::TreeUtils::RNoCleanupNotifier fNoCleanupNotifier
RColumnReaderBase * AddTreeColumnReader(unsigned int slot, std::string_view col, std::unique_ptr< RColumnReaderBase > &&reader, const std::type_info &ti)
Register a new RTreeColumnReader with this RLoopManager.
std::vector< RDFInternal::RActionBase * > fRunActions
Non-owning pointers to actions already run.
RLoopManager(const ColumnNames_t &defaultColumns={})
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< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
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.
RDataSource * GetDataSource() const
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".
RDFInternal::RNewSampleNotifier fNewSampleNotifier
std::pair< ULong64_t, ULong64_t > fEmptyEntryRange
Range of entries created when no data source is specified.
std::unique_ptr< RDataSource > fDataSource
Owning pointer to a data-source object.
void DataSourceThreadTask(const std::pair< ULong64_t, ULong64_t > &entryRange, ROOT::Internal::RSlotStack &slotStack, std::atomic< ULong64_t > &entryCount)
The task run by every thread on the input entry range, for the generic RDataSource.
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 SetDataSource(std::unique_ptr< ROOT::RDF::RDataSource > dataSource)
void RegisterCallback(ULong64_t everyNEvents, std::function< void(unsigned int)> &&f)
void SetTTreeLifeline(std::any lifeline)
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.
RColumnReaderBase * GetDatasetColumnReader(unsigned int slot, std::string_view col, const std::type_info &ti) const
void RunTreeProcessorMT()
Run event loop over one or multiple ROOT files, in parallel.
std::shared_ptr< ROOT::Internal::RSlotStack > SlotStack() const
Create a slot stack with the desired number of slots or reuse a shared instance.
void Deregister(RDFInternal::RActionBase *actionPtr)
ELoopType fLoopType
The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
void InitNodes()
Initialize all nodes of the functional graph before running the event loop.
bool HasDataSourceColumnReaders(std::string_view 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
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 list of N indexes (0 to size - 1).
The dataset specification for RDataFrame.
RDataSource defines an API that RDataFrame can use to read arbitrary data formats.
virtual void Finalize()
Convenience method called after concluding an event-loop.
virtual void InitSlot(unsigned int, ULong64_t)
Convenience method called at the start of the data processing associated to a slot.
virtual void FinalizeSlot(unsigned int)
Convenience method called at the end of the data processing associated to a slot.
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:65
const_iterator begin() const
const_iterator end() const
This class provides a simple interface to execute the same task multiple times in parallel threads,...
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()
A chain is a collection of files containing TTree objects.
Definition TChain.h:33
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
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
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:4131
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:49
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.
Basic string class.
Definition TString.h:139
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition TSystem.cxx:1286
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
A TTree represents a columnar dataset.
Definition TTree.h:84
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:4846
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition TTree.cxx:5299
static void SetMaxTreeSize(Long64_t maxsize=100000000000LL)
Set the maximum size in bytes of a Tree file (static function).
Definition TTree.cxx:9309
virtual TObjArray * GetListOfBranches()
Definition TTree.h:543
virtual TTree * GetTree() const
Definition TTree.h:572
virtual const char * GetFriendAlias(TTree *) const
If the 'tree' is a friend, this method returns its alias name.
Definition TTree.cxx:6042
This class represents a WWW compatible URL.
Definition TUrl.h:33
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::RLogChannel & RDFLogChannel()
Definition RDFUtils.cxx:41
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > CreateLMFromFile(std::string_view datasetName, std::string_view fileNameGlob, const std::vector< std::string > &defaultColumns)
Create an RLoopManager opening a file and checking the data format of the dataset.
std::shared_ptr< ROOT::Detail::RDF::RLoopManager > CreateLMFromRNTuple(std::string_view datasetName, std::string_view fileNameGlob, const std::vector< std::string > &defaultColumns)
Create an RLoopManager that reads an RNTuple.
void RunFinalChecks(const ROOT::RDF::RDataSource &ds, bool nodesLeftNotRun)
Definition RDFUtils.cxx:581
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:305
void CallInitializeWithOpts(ROOT::RDF::RDataSource &ds, const std::set< std::string > &suppressErrorsForMissingColumns)
Definition RDFUtils.cxx:563
void Erase(const T &that, std::vector< T > &v)
Erase that element from vector v
Definition Utils.hxx:192
ROOT::RDF::RSampleInfo CreateSampleInfo(const ROOT::RDF::RDataSource &ds, const std::unordered_map< std::string, ROOT::RDF::Experimental::RSample * > &sampleMap)
Definition RDFUtils.cxx:574
std::unique_ptr< ROOT::Detail::RDF::RColumnReaderBase > CreateColumnReader(ROOT::RDF::RDataSource &ds, unsigned int slot, std::string_view col, const std::type_info &tid, TTreeReader *treeReader)
Definition RDFUtils.cxx:592
void 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:349
void ProcessMT(ROOT::RDF::RDataSource &ds, ROOT::Detail::RDF::RLoopManager &lm)
Definition RDFUtils.cxx:586
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::string > ExpandGlob(const std::string &glob)
Expands input glob into a collection of full paths to files.
auto MakeAliasedSharedPtr(T *rawPtr)
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....
std::vector< std::string > ColumnNames_t
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:595
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)
ROOT::Detail::RDF::RLoopManager & fLM
RDSRangeRAII(ROOT::Detail::RDF::RLoopManager &lm, unsigned int slot, ULong64_t firstEntry, TTreeReader *treeReader=nullptr)
A RAII object to pop and push slot numbers from a RSlotStack object.