Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RLoopManager.cxx
Go to the documentation of this file.
1/*************************************************************************
2 * Copyright (C) 1995-2021, Rene Brun and Fons Rademakers. *
3 * All rights reserved. *
4 * *
5 * For the licensing terms see $ROOTSYS/LICENSE. *
6 * For the list of contributors see $ROOTSYS/README/CREDITS. *
7 *************************************************************************/
8
9#include "RConfigure.h" // R__USE_IMT
10#include "ROOT/RDataSource.hxx"
12#include "ROOT/InternalTreeUtils.hxx" // GetTreeFullPaths
19#include "ROOT/RLogger.hxx"
20#include "RtypesCore.h" // Long64_t
21#include "TStopwatch.h"
22#include "TBranchElement.h"
23#include "TBranchObject.h"
24#include "TChain.h"
25#include "TEntryList.h"
26#include "TFile.h"
27#include "TFriendElement.h"
28#include "TROOT.h" // IsImplicitMTEnabled
29#include "TTreeReader.h"
30#include "TTree.h" // For MaxTreeSizeRAII. Revert when #6640 will be solved.
31
32#ifdef R__USE_IMT
35#include "ROOT/RSlotStack.hxx"
36#endif
37
38#include <algorithm>
39#include <atomic>
40#include <cassert>
41#include <functional>
42#include <iostream>
43#include <memory>
44#include <stdexcept>
45#include <string>
46#include <sstream>
47#include <thread>
48#include <unordered_map>
49#include <vector>
50#include <set>
51#include <limits> // For MaxTreeSizeRAII. Revert when #6640 will be solved.
52
53using namespace ROOT::Detail::RDF;
54using namespace ROOT::Internal::RDF;
55
56namespace {
57/// A helper function that returns all RDF code that is currently scheduled for just-in-time compilation.
58/// This allows different RLoopManager instances to share these data.
59/// We want RLoopManagers to be able to add their code to a global "code to execute via cling",
60/// so that, lazily, we can jit everything that's needed by all RDFs in one go, which is potentially
61/// much faster than jitting each RLoopManager's code separately.
62static std::string &GetCodeToJit()
63{
64 static std::string code;
65 return code;
66}
67
68static bool ContainsLeaf(const std::set<TLeaf *> &leaves, TLeaf *leaf)
69{
70 return (leaves.find(leaf) != leaves.end());
71}
72
73///////////////////////////////////////////////////////////////////////////////
74/// This overload does not check whether the leaf/branch is already in bNamesReg. In case this is a friend leaf/branch,
75/// `allowDuplicates` controls whether we add both `friendname.bname` and `bname` or just the shorter version.
76static void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
77 const std::string &friendName, bool allowDuplicates)
78{
79 if (!friendName.empty()) {
80 // In case of a friend tree, users might prepend its name/alias to the branch names
81 const auto friendBName = friendName + "." + branchName;
82 if (bNamesReg.insert(friendBName).second)
83 bNames.push_back(friendBName);
84 }
85
86 if (allowDuplicates || friendName.empty()) {
87 if (bNamesReg.insert(branchName).second)
88 bNames.push_back(branchName);
89 }
90}
91
92///////////////////////////////////////////////////////////////////////////////
93/// This overload makes sure that the TLeaf has not been already inserted.
94static void InsertBranchName(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
95 const std::string &friendName, std::set<TLeaf *> &foundLeaves, TLeaf *leaf,
96 bool allowDuplicates)
97{
98 const bool canAdd = allowDuplicates ? true : !ContainsLeaf(foundLeaves, leaf);
99 if (!canAdd) {
100 return;
101 }
102
103 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
104
105 foundLeaves.insert(leaf);
106}
107
108static void ExploreBranch(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames, TBranch *b,
109 std::string prefix, std::string &friendName, bool allowDuplicates)
110{
111 for (auto sb : *b->GetListOfBranches()) {
112 TBranch *subBranch = static_cast<TBranch *>(sb);
113 auto subBranchName = std::string(subBranch->GetName());
114 auto fullName = prefix + subBranchName;
115
116 std::string newPrefix;
117 if (!prefix.empty())
118 newPrefix = fullName + ".";
119
120 ExploreBranch(t, bNamesReg, bNames, subBranch, newPrefix, friendName, allowDuplicates);
121
122 auto branchDirectlyFromTree = t.GetBranch(fullName.c_str());
123 if (!branchDirectlyFromTree)
124 branchDirectlyFromTree = t.FindBranch(fullName.c_str()); // try harder
125 if (branchDirectlyFromTree)
126 InsertBranchName(bNamesReg, bNames, std::string(branchDirectlyFromTree->GetFullName()), friendName,
127 allowDuplicates);
128
129 if (bNamesReg.find(subBranchName) == bNamesReg.end() && t.GetBranch(subBranchName.c_str()))
130 InsertBranchName(bNamesReg, bNames, subBranchName, friendName, allowDuplicates);
131 }
132}
133
134static void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
135 std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
136{
137 std::set<TLeaf *> foundLeaves;
138 if (!analysedTrees.insert(&t).second) {
139 return;
140 }
141
142 const auto branches = t.GetListOfBranches();
143 // Getting the branches here triggered the read of the first file of the chain if t is a chain.
144 // We check if a tree has been successfully read, otherwise we throw (see ROOT-9984) to avoid further
145 // operations
146 if (!t.GetTree()) {
147 std::string err("GetBranchNames: error in opening the tree ");
148 err += t.GetName();
149 throw std::runtime_error(err);
150 }
151 if (branches) {
152 for (auto b : *branches) {
153 TBranch *branch = static_cast<TBranch *>(b);
154 const auto branchName = std::string(branch->GetName());
155 if (branch->IsA() == TBranch::Class()) {
156 // Leaf list
157 auto listOfLeaves = branch->GetListOfLeaves();
158 if (listOfLeaves->GetEntriesUnsafe() == 1) {
159 auto leaf = static_cast<TLeaf *>(listOfLeaves->UncheckedAt(0));
160 InsertBranchName(bNamesReg, bNames, branchName, friendName, foundLeaves, leaf, allowDuplicates);
161 }
162
163 for (auto leaf : *listOfLeaves) {
164 auto castLeaf = static_cast<TLeaf *>(leaf);
165 const auto leafName = std::string(leaf->GetName());
166 const auto fullName = branchName + "." + leafName;
167 InsertBranchName(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
168 }
169 } else if (branch->IsA() == TBranchObject::Class()) {
170 // TBranchObject
171 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName, allowDuplicates);
172 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
173 } else {
174 // TBranchElement
175 // Check if there is explicit or implicit dot in the name
176
177 bool dotIsImplied = false;
178 auto be = dynamic_cast<TBranchElement *>(b);
179 if (!be)
180 throw std::runtime_error("GetBranchNames: unsupported branch type");
181 // TClonesArray (3) and STL collection (4)
182 if (be->GetType() == 3 || be->GetType() == 4)
183 dotIsImplied = true;
184
185 if (dotIsImplied || branchName.back() == '.')
186 ExploreBranch(t, bNamesReg, bNames, branch, "", friendName, allowDuplicates);
187 else
188 ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName, allowDuplicates);
189
190 InsertBranchName(bNamesReg, bNames, branchName, friendName, allowDuplicates);
191 }
192 }
193 }
194
195 // The list of friends needs to be accessed via GetTree()->GetListOfFriends()
196 // (and not via GetListOfFriends() directly), otherwise when `t` is a TChain we
197 // might not recover the list correctly (https://github.com/root-project/root/issues/6741).
198 auto friendTrees = t.GetTree()->GetListOfFriends();
199
200 if (!friendTrees)
201 return;
202
203 for (auto friendTreeObj : *friendTrees) {
204 auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
205
206 std::string frName;
207 auto alias = t.GetFriendAlias(friendTree);
208 if (alias != nullptr)
209 frName = std::string(alias);
210 else
211 frName = std::string(friendTree->GetName());
212
213 GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
214 }
215}
216
217static void ThrowIfNSlotsChanged(unsigned int nSlots)
218{
219 const auto currentSlots = RDFInternal::GetNSlots();
220 if (currentSlots != nSlots) {
221 std::string msg = "RLoopManager::Run: when the RDataFrame was constructed the number of slots required was " +
222 std::to_string(nSlots) + ", but when starting the event loop it was " +
223 std::to_string(currentSlots) + ".";
224 if (currentSlots > nSlots)
225 msg += " Maybe EnableImplicitMT() was called after the RDataFrame was constructed?";
226 else
227 msg += " Maybe DisableImplicitMT() was called after the RDataFrame was constructed?";
228 throw std::runtime_error(msg);
229 }
230}
231
232/**
233\struct MaxTreeSizeRAII
234\brief Scope-bound change of `TTree::fgMaxTreeSize`.
235
236This RAII object stores the current value result of `TTree::GetMaxTreeSize`,
237changes it to maximum at construction time and restores it back at destruction
238time. Needed for issue #6523 and should be reverted when #6640 will be solved.
239*/
240struct MaxTreeSizeRAII {
241 Long64_t fOldMaxTreeSize;
242
243 MaxTreeSizeRAII() : fOldMaxTreeSize(TTree::GetMaxTreeSize())
244 {
245 TTree::SetMaxTreeSize(std::numeric_limits<Long64_t>::max());
246 }
247
248 ~MaxTreeSizeRAII() { TTree::SetMaxTreeSize(fOldMaxTreeSize); }
249};
250
251struct DatasetLogInfo {
252 std::string fDataSet;
253 ULong64_t fRangeStart;
254 ULong64_t fRangeEnd;
255 unsigned int fSlot;
256};
257
258std::string LogRangeProcessing(const DatasetLogInfo &info)
259{
260 std::stringstream msg;
261 msg << "Processing " << info.fDataSet << ": entry range [" << info.fRangeStart << "," << info.fRangeEnd - 1
262 << "], using slot " << info.fSlot << " in thread " << std::this_thread::get_id() << '.';
263 return msg.str();
264}
265
266DatasetLogInfo TreeDatasetLogInfo(const TTreeReader &r, unsigned int slot)
267{
268 const auto tree = r.GetTree();
269 const auto chain = dynamic_cast<TChain *>(tree);
270 std::string what;
271 if (chain) {
272 auto files = chain->GetListOfFiles();
273 std::vector<std::string> treeNames;
274 std::vector<std::string> fileNames;
275 for (TObject *f : *files) {
276 treeNames.emplace_back(f->GetName());
277 fileNames.emplace_back(f->GetTitle());
278 }
279 what = "trees {";
280 for (const auto &t : treeNames) {
281 what += t + ",";
282 }
283 what.back() = '}';
284 what += " in files {";
285 for (const auto &f : fileNames) {
286 what += f + ",";
287 }
288 what.back() = '}';
289 } else {
290 assert(tree != nullptr); // to make clang-tidy happy
291 const auto treeName = tree->GetName();
292 what = std::string("tree \"") + treeName + "\"";
293 const auto file = tree->GetCurrentFile();
294 if (file)
295 what += std::string(" in file \"") + file->GetName() + "\"";
296 }
297 const auto entryRange = r.GetEntriesRange();
298 const ULong64_t end = entryRange.second == -1ll ? tree->GetEntries() : entryRange.second;
299 return {std::move(what), static_cast<ULong64_t>(entryRange.first), end, slot};
300}
301
302static auto MakeDatasetColReadersKey(const std::string &colName, const std::type_info &ti)
303{
304 // We use a combination of column name and column type name as the key because in some cases we might end up
305 // with concrete readers that use different types for the same column, e.g. std::vector and RVec here:
306 // df.Sum<vector<int>>("stdVectorBranch");
307 // df.Sum<RVecI>("stdVectorBranch");
308 return colName + ':' + ti.name();
309}
310} // anonymous namespace
311
312namespace ROOT {
313namespace Detail {
314namespace RDF {
315
316/// A RAII object that calls RLoopManager::CleanUpTask at destruction
319 unsigned int fArg;
321
322 RCallCleanUpTask(RLoopManager &lm, unsigned int arg = 0u, TTreeReader *reader = nullptr)
323 : fLoopManager(lm), fArg(arg), fReader(reader)
324 {
325 }
327};
328
329} // namespace RDF
330} // namespace Detail
331} // namespace ROOT
332
333///////////////////////////////////////////////////////////////////////////////
334/// Get all the branches names, including the ones of the friend trees
336{
337 std::set<std::string> bNamesSet;
338 ColumnNames_t bNames;
339 std::set<TTree *> analysedTrees;
340 std::string emptyFrName = "";
341 GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
342 return bNames;
343}
344
346 : fTree(std::shared_ptr<TTree>(tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
347 fNSlots(RDFInternal::GetNSlots()),
348 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles),
349 fNewSampleNotifier(fNSlots), fSampleInfos(fNSlots), fDatasetColumnReaders(fNSlots)
350{
351}
352
354 : fNEmptyEntries(nEmptyEntries), fNSlots(RDFInternal::GetNSlots()),
355 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles), fNewSampleNotifier(fNSlots),
356 fSampleInfos(fNSlots), fDatasetColumnReaders(fNSlots)
357{
358}
359
360RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
361 : fDefaultColumns(defaultBranches), fNSlots(RDFInternal::GetNSlots()),
362 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
363 fDataSource(std::move(ds)), fNewSampleNotifier(fNSlots), fSampleInfos(fNSlots), fDatasetColumnReaders(fNSlots)
364{
365 fDataSource->SetNSlots(fNSlots);
366}
367
369 : fBeginEntry(spec.GetEntryRangeBegin()),
370 fEndEntry(spec.GetEntryRangeEnd()),
371 fSamples(spec.MoveOutSamples()),
372 fNSlots(RDFInternal::GetNSlots()),
373 fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles),
374 fNewSampleNotifier(fNSlots),
375 fSampleInfos(fNSlots),
376 fDatasetColumnReaders(fNSlots)
377{
378 auto chain = std::make_shared<TChain>("");
379 for (auto &sample : fSamples) {
380 const auto &trees = sample.GetTreeNames();
381 const auto &files = sample.GetFileNameGlobs();
382 for (auto i = 0u; i < files.size(); ++i) {
383 // We need to use `<filename>?#<treename>` as an argument to TChain::Add
384 // (see https://github.com/root-project/root/pull/8820 for why)
385 const auto fullpath = files[i] + "?#" + trees[i];
386 chain->Add(fullpath.c_str());
387 // ...but instead we use `<filename>/<treename>` as a sample ID (cannot
388 // change this easily because of backward compatibility: the sample ID
389 // is exposed to users via RSampleInfo and DefinePerSample).
390 const auto sampleId = files[i] + '/' + trees[i];
391 fSampleMap.insert({sampleId, &sample});
392 }
393 }
394
395 SetTree(std::move(chain));
396
397 // Create the friends from the list of friends
398 const auto &friendInfo = spec.GetFriendInfo();
399 const auto &friendNames = friendInfo.fFriendNames;
400 const auto &friendFileNames = friendInfo.fFriendFileNames;
401 const auto &friendChainSubNames = friendInfo.fFriendChainSubNames;
402 const auto nFriends = friendNames.size();
403
404 for (auto i = 0u; i < nFriends; ++i) {
405 const auto &thisFriendName = friendNames[i].first;
406 const auto &thisFriendAlias = friendNames[i].second;
407 const auto &thisFriendFiles = friendFileNames[i];
408 const auto &thisFriendChainSubNames = friendChainSubNames[i];
409
410 // Build a friend chain
411 auto frChain = std::make_unique<TChain>(thisFriendName.c_str());
412 const auto nFileNames = friendFileNames[i].size();
413 if (thisFriendChainSubNames.empty()) {
414 // If there are no chain subnames, the friend was a TTree. It's safe
415 // to add to the chain the filename directly.
416 for (auto j = 0u; j < nFileNames; ++j) {
417 frChain->Add(thisFriendFiles[j].c_str());
418 }
419 } else {
420 // Otherwise, the new friend chain needs to be built using the nomenclature
421 // "filename#?treename" as argument to `TChain::Add`.
422 // See https://github.com/root-project/root/pull/8820 for why "filename#?treename"
423 // is better than "filename#?treename".
424 for (auto j = 0u; j < nFileNames; ++j)
425 frChain->Add((thisFriendFiles[j] + "?#" + thisFriendChainSubNames[j]).c_str());
426 }
427
428 // Make it friends with the main chain
429 fTree->AddFriend(frChain.get(), thisFriendAlias.c_str());
430 // the friend trees must have same lifetime as fTree
431 fFriends.emplace_back(std::move(frChain));
432 }
433}
434
435/// Run event loop with no source files, in parallel.
437{
438#ifdef R__USE_IMT
440 // Working with an empty tree.
441 // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
442 const auto nEntriesPerSlot = fNEmptyEntries / (fNSlots * 2);
443 auto remainder = fNEmptyEntries % (fNSlots * 2);
444 std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
445 ULong64_t start = 0;
446 while (start < fNEmptyEntries) {
447 ULong64_t end = start + nEntriesPerSlot;
448 if (remainder > 0) {
449 ++end;
450 --remainder;
451 }
452 entryRanges.emplace_back(start, end);
453 start = end;
454 }
455
456 // Each task will generate a subrange of entries
457 auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
458 ROOT::Internal::RSlotStackRAII slotRAII(slotStack);
459 auto slot = slotRAII.fSlot;
460 RCallCleanUpTask cleanup(*this, slot);
461 InitNodeSlots(nullptr, slot);
462 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({"an empty source", range.first, range.second, slot});
463 try {
464 UpdateSampleInfo(slot, range);
465 for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
466 RunAndCheckFilters(slot, currEntry);
467 }
468 } catch (...) {
469 // Error might throw in experiment frameworks like CMSSW
470 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
471 throw;
472 }
473 };
474
476 pool.Foreach(genFunction, entryRanges);
477
478#endif // not implemented otherwise
479}
480
481/// Run event loop with no source files, in sequence.
483{
484 InitNodeSlots(nullptr, 0);
485 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({"an empty source", 0, fNEmptyEntries, 0u});
487 try {
488 UpdateSampleInfo(/*slot*/0, {0, fNEmptyEntries});
489 for (ULong64_t currEntry = 0; currEntry < fNEmptyEntries && fNStopsReceived < fNChildren; ++currEntry) {
490 RunAndCheckFilters(0, currEntry);
491 }
492 } catch (...) {
493 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
494 throw;
495 }
496}
497
498/// Run event loop over one or multiple ROOT files, in parallel.
500{
501#ifdef R__USE_IMT
502 if (fEndEntry == fBeginEntry) // empty range => no work needed
503 return;
505 const auto &entryList = fTree->GetEntryList() ? *fTree->GetEntryList() : TEntryList();
506 auto tp = (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
507 ? std::make_unique<ROOT::TTreeProcessorMT>(*fTree, fNSlots, std::make_pair(fBeginEntry, fEndEntry))
508 : std::make_unique<ROOT::TTreeProcessorMT>(*fTree, entryList, fNSlots);
509
510 std::atomic<ULong64_t> entryCount(0ull);
511
512 tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
513 ROOT::Internal::RSlotStackRAII slotRAII(slotStack);
514 auto slot = slotRAII.fSlot;
515 RCallCleanUpTask cleanup(*this, slot, &r);
516 InitNodeSlots(&r, slot);
517 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing(TreeDatasetLogInfo(r, slot));
518 const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
519 const auto nEntries = entryRange.second - entryRange.first;
520 auto count = entryCount.fetch_add(nEntries);
521 try {
522 // recursive call to check filters and conditionally execute actions
523 while (r.Next()) {
524 if (fNewSampleNotifier.CheckFlag(slot)) {
525 UpdateSampleInfo(slot, r);
526 }
527 RunAndCheckFilters(slot, count++);
528 }
529 } catch (...) {
530 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
531 throw;
532 }
533 // fNStopsReceived < fNChildren is always true at the moment as we don't support event loop early quitting in
534 // multi-thread runs, but it costs nothing to be safe and future-proof in case we add support for that later.
535 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
536 // something went wrong in the TTreeReader event loop
537 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
538 std::to_string(r.GetEntryStatus()));
539 }
540 });
541#endif // no-op otherwise (will not be called)
542}
543
544/// Run event loop over one or multiple ROOT files, in sequence.
546{
547 TTreeReader r(fTree.get(), fTree->GetEntryList());
548 if (0 == fTree->GetEntriesFast() || fBeginEntry == fEndEntry)
549 return;
550 // Apply the range if there is any
551 // In case of a chain with a total of N entries, calling SetEntriesRange(N + 1, ...) does not error out
552 // This is a bug, reported here: https://github.com/root-project/root/issues/10774
553 if (fBeginEntry != 0 || fEndEntry != std::numeric_limits<Long64_t>::max())
554 if (r.SetEntriesRange(fBeginEntry, fEndEntry) != TTreeReader::kEntryValid)
555 throw std::logic_error("Something went wrong in initializing the TTreeReader.");
556
557 RCallCleanUpTask cleanup(*this, 0u, &r);
558 InitNodeSlots(&r, 0);
559 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing(TreeDatasetLogInfo(r, 0u));
560
561 // recursive call to check filters and conditionally execute actions
562 // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
563 try {
564 while (r.Next() && fNStopsReceived < fNChildren) {
566 UpdateSampleInfo(/*slot*/0, r);
567 }
568 RunAndCheckFilters(0, r.GetCurrentEntry());
569 }
570 } catch (...) {
571 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
572 throw;
573 }
574 if (r.GetEntryStatus() != TTreeReader::kEntryBeyondEnd && fNStopsReceived < fNChildren) {
575 // something went wrong in the TTreeReader event loop
576 throw std::runtime_error("An error was encountered while processing the data. TTreeReader status code is: " +
577 std::to_string(r.GetEntryStatus()));
578 }
579}
580
581/// Run event loop over data accessed through a DataSource, in sequence.
583{
584 assert(fDataSource != nullptr);
585 fDataSource->CallInitialize();
586 auto ranges = fDataSource->GetEntryRanges();
587 while (!ranges.empty() && fNStopsReceived < fNChildren) {
588 InitNodeSlots(nullptr, 0u);
589 fDataSource->InitSlot(0u, 0ull);
591 try {
592 for (const auto &range : ranges) {
593 const auto start = range.first;
594 const auto end = range.second;
595 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, 0u});
596 for (auto entry = start; entry < end && fNStopsReceived < fNChildren; ++entry) {
597 if (fDataSource->SetEntry(0u, entry)) {
598 RunAndCheckFilters(0u, entry);
599 }
600 }
601 }
602 } catch (...) {
603 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
604 throw;
605 }
606 fDataSource->CallFinalizeSlot(0u);
607 ranges = fDataSource->GetEntryRanges();
608 }
609 fDataSource->CallFinalize();
610}
611
612/// Run event loop over data accessed through a DataSource, in parallel.
614{
615#ifdef R__USE_IMT
616 assert(fDataSource != nullptr);
619
620 // Each task works on a subrange of entries
621 auto runOnRange = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
622 ROOT::Internal::RSlotStackRAII slotRAII(slotStack);
623 const auto slot = slotRAII.fSlot;
624 InitNodeSlots(nullptr, slot);
625 RCallCleanUpTask cleanup(*this, slot);
626 fDataSource->InitSlot(slot, range.first);
627 const auto start = range.first;
628 const auto end = range.second;
629 R__LOG_DEBUG(0, RDFLogChannel()) << LogRangeProcessing({fDataSource->GetLabel(), start, end, slot});
630 try {
631 for (auto entry = start; entry < end; ++entry) {
632 if (fDataSource->SetEntry(slot, entry)) {
633 RunAndCheckFilters(slot, entry);
634 }
635 }
636 } catch (...) {
637 std::cerr << "RDataFrame::Run: event loop was interrupted\n";
638 throw;
639 }
640 fDataSource->CallFinalizeSlot(slot);
641 };
642
643 fDataSource->CallInitialize();
644 auto ranges = fDataSource->GetEntryRanges();
645 while (!ranges.empty()) {
646 pool.Foreach(runOnRange, ranges);
647 ranges = fDataSource->GetEntryRanges();
648 }
649 fDataSource->CallFinalize();
650#endif // not implemented otherwise (never called)
651}
652
653/// Execute actions and make sure named filters are called for each event.
654/// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
655void RLoopManager::RunAndCheckFilters(unsigned int slot, Long64_t entry)
656{
657 // data-block callbacks run before the rest of the graph
658 if (fNewSampleNotifier.CheckFlag(slot)) {
659 for (auto &callback : fSampleCallbacks)
660 callback.second(slot, fSampleInfos[slot]);
662 }
663
664 for (auto &actionPtr : fBookedActions)
665 actionPtr->Run(slot, entry);
666 for (auto &namedFilterPtr : fBookedNamedFilters)
667 namedFilterPtr->CheckFilters(slot, entry);
668 for (auto &callback : fCallbacks)
669 callback(slot);
670}
671
672/// Build TTreeReaderValues for all nodes
673/// This method loops over all filters, actions and other booked objects and
674/// calls their `InitSlot` method, to get them ready for running a task.
676{
677 SetupSampleCallbacks(r, slot);
678 for (auto &ptr : fBookedActions)
679 ptr->InitSlot(r, slot);
680 for (auto &ptr : fBookedFilters)
681 ptr->InitSlot(r, slot);
682 for (auto &ptr : fBookedDefines)
683 ptr->InitSlot(r, slot);
684 for (auto &ptr : fBookedVariations)
685 ptr->InitSlot(r, slot);
686
687 for (auto &callback : fCallbacksOnce)
688 callback(slot);
689}
690
692 if (r != nullptr) {
693 // we need to set a notifier so that we run the callbacks every time we switch to a new TTree
694 // `PrependLink` inserts this notifier into the TTree/TChain's linked list of notifiers
695 fNewSampleNotifier.GetChainNotifyLink(slot).PrependLink(*r->GetTree());
696 }
697 // Whatever the data source, initially set the "new data block" flag:
698 // - for TChains, this ensures that we don't skip the first data block because
699 // the correct tree is already loaded
700 // - for RDataSources and empty sources, which currently don't have data blocks, this
701 // ensures that we run once per task
703}
704
705void RLoopManager::UpdateSampleInfo(unsigned int slot, const std::pair<ULong64_t, ULong64_t> &range) {
707 "Empty source, range: {" + std::to_string(range.first) + ", " + std::to_string(range.second) + "}", range);
708}
709
711 // one GetTree to retrieve the TChain, another to retrieve the underlying TTree
712 auto *tree = r.GetTree()->GetTree();
713 R__ASSERT(tree != nullptr);
714 const std::string treename = ROOT::Internal::TreeUtils::GetTreeFullPaths(*tree)[0];
715 auto *file = tree->GetCurrentFile();
716 const std::string fname = file != nullptr ? file->GetName() : "#inmemorytree#";
717
718 std::pair<Long64_t, Long64_t> range = r.GetEntriesRange();
719 R__ASSERT(range.first >= 0);
720 if (range.second == -1) {
721 range.second = tree->GetEntries(); // convert '-1', i.e. 'until the end', to the actual entry number
722 }
723 const std::string &id = fname + '/' + treename;
724 fSampleInfos[slot] = fSampleMap.empty() ? RSampleInfo(id, range) : RSampleInfo(id, range, fSampleMap[id]);
725}
726
727/// Initialize all nodes of the functional graph before running the event loop.
728/// This method is called once per event-loop and performs generic initialization
729/// operations that do not depend on the specific processing slot (i.e. operations
730/// that are common for all threads).
732{
734 for (auto &filter : fBookedFilters)
735 filter->InitNode();
736 for (auto &range : fBookedRanges)
737 range->InitNode();
738 for (auto &ptr : fBookedActions)
739 ptr->Initialize();
740}
741
742/// Perform clean-up operations. To be called at the end of each event loop.
744{
745 fMustRunNamedFilters = false;
746
747 // forget RActions and detach TResultProxies
748 for (auto &ptr : fBookedActions)
749 ptr->Finalize();
750
751 fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
752 fBookedActions.clear();
753
754 // reset children counts
755 fNChildren = 0;
756 fNStopsReceived = 0;
757 for (auto &ptr : fBookedFilters)
758 ptr->ResetChildrenCount();
759 for (auto &ptr : fBookedRanges)
760 ptr->ResetChildrenCount();
761
762 fCallbacks.clear();
763 fCallbacksOnce.clear();
764}
765
766/// Perform clean-up operations. To be called at the end of each task execution.
767void RLoopManager::CleanUpTask(TTreeReader *r, unsigned int slot)
768{
769 if (r != nullptr)
770 fNewSampleNotifier.GetChainNotifyLink(slot).RemoveLink(*r->GetTree());
771 for (auto &ptr : fBookedActions)
772 ptr->FinalizeSlot(slot);
773 for (auto &ptr : fBookedFilters)
774 ptr->FinalizeSlot(slot);
775 for (auto &ptr : fBookedDefines)
776 ptr->FinalizeSlot(slot);
777
779 // we are reading from a tree/chain and we need to re-create the RTreeColumnReaders at every task
780 // because the TTreeReader object changes at every task
781 for (auto &v : fDatasetColumnReaders[slot])
782 v.second.reset();
783 }
784}
785
786/// Add RDF nodes that require just-in-time compilation to the computation graph.
787/// This method also clears the contents of GetCodeToJit().
789{
790 // TODO this should be a read lock unless we find GetCodeToJit non-empty
792
793 const std::string code = std::move(GetCodeToJit());
794 if (code.empty()) {
795 R__LOG_INFO(RDFLogChannel()) << "Nothing to jit and execute.";
796 return;
797 }
798
799 TStopwatch s;
800 s.Start();
801 RDFInternal::InterpreterCalc(code, "RLoopManager::Run");
802 s.Stop();
803 R__LOG_INFO(RDFLogChannel()) << "Just-in-time compilation phase completed"
804 << (s.RealTime() > 1e-3 ? " in " + std::to_string(s.RealTime()) + " seconds."
805 : " in less than 1ms.");
806}
807
808/// Trigger counting of number of children nodes for each node of the functional graph.
809/// This is done once before starting the event loop. Each action sends an `increase children count` signal
810/// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
811/// children counter. Each node only propagates the signal once, even if it receives it multiple times.
812/// Named filters also send an `increase children count` signal, just like actions, as they always execute during
813/// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
815{
816 for (auto &actionPtr : fBookedActions)
817 actionPtr->TriggerChildrenCount();
818 for (auto &namedFilterPtr : fBookedNamedFilters)
819 namedFilterPtr->TriggerChildrenCount();
820}
821
822/// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
823/// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
824/// The jitting phase is skipped if the `jit` parameter is `false` (unsafe, use with care).
825void RLoopManager::Run(bool jit)
826{
827 // Change value of TTree::GetMaxTreeSize only for this scope. Revert when #6640 will be solved.
828 MaxTreeSizeRAII ctxtmts;
829
830 R__LOG_INFO(RDFLogChannel()) << "Starting event loop number " << fNRuns << '.';
831
832 ThrowIfNSlotsChanged(GetNSlots());
833
834 if (jit)
835 Jit();
836
837 InitNodes();
838
839 TStopwatch s;
840 s.Start();
841 switch (fLoopType) {
845 case ELoopType::kNoFiles: RunEmptySource(); break;
848 }
849 s.Stop();
850
851 CleanUpNodes();
852
853 fNRuns++;
854
855 R__LOG_INFO(RDFLogChannel()) << "Finished event loop number " << fNRuns - 1 << " (" << s.CpuTime() << "s CPU, "
856 << s.RealTime() << "s elapsed).";
857}
858
859/// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
861{
862 return fDefaultColumns;
863}
864
866{
867 return fTree.get();
868}
869
871{
872 fBookedActions.emplace_back(actionPtr);
873 AddSampleCallback(actionPtr, actionPtr->GetSampleCallback());
874}
875
877{
880 fSampleCallbacks.erase(actionPtr);
881}
882
884{
885 fBookedFilters.emplace_back(filterPtr);
886 if (filterPtr->HasName()) {
887 fBookedNamedFilters.emplace_back(filterPtr);
889 }
890}
891
893{
896}
897
899{
900 fBookedRanges.emplace_back(rangePtr);
901}
902
904{
906}
907
909{
910 fBookedDefines.emplace_back(ptr);
911}
912
914{
916 fSampleCallbacks.erase(ptr);
917}
918
920{
921 fBookedVariations.emplace_back(v);
922}
923
925{
927}
928
929// dummy call, end of recursive chain of calls
931{
932 return true;
933}
934
935/// Call `FillReport` on all booked filters
937{
938 for (const auto &fPtr : fBookedNamedFilters)
939 fPtr->FillReport(rep);
940}
941
942void RLoopManager::SetTree(std::shared_ptr<TTree> tree)
943{
944 fTree = std::move(tree);
945
946 TChain *ch = nullptr;
947 if ((ch = dynamic_cast<TChain *>(fTree.get())))
949}
950
951void RLoopManager::ToJitExec(const std::string &code) const
952{
954 GetCodeToJit().append(code);
955}
956
957void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
958{
959 if (everyNEvents == 0ull)
960 fCallbacksOnce.emplace_back(std::move(f), fNSlots);
961 else
962 fCallbacks.emplace_back(everyNEvents, std::move(f), fNSlots);
963}
964
965std::vector<std::string> RLoopManager::GetFiltersNames()
966{
967 std::vector<std::string> filters;
968 for (auto &filter : fBookedFilters) {
969 auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
970 filters.push_back(name);
971 }
972 return filters;
973}
974
975std::vector<RNodeBase *> RLoopManager::GetGraphEdges() const
976{
977 std::vector<RNodeBase *> nodes(fBookedFilters.size() + fBookedRanges.size());
978 auto it = std::copy(fBookedFilters.begin(), fBookedFilters.end(), nodes.begin());
979 std::copy(fBookedRanges.begin(), fBookedRanges.end(), it);
980 return nodes;
981}
982
983std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions() const
984{
985 std::vector<RDFInternal::RActionBase *> actions(fBookedActions.size() + fRunActions.size());
986 auto it = std::copy(fBookedActions.begin(), fBookedActions.end(), actions.begin());
987 std::copy(fRunActions.begin(), fRunActions.end(), it);
988 return actions;
989}
990
991std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph(
992 std::unordered_map<void *, std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode>> &visitedMap)
993{
994 // If there is already a node for the RLoopManager return it. If there is not, return a new one.
995 auto duplicateRLoopManagerIt = visitedMap.find((void *)this);
996 if (duplicateRLoopManagerIt != visitedMap.end())
997 return duplicateRLoopManagerIt->second;
998
999 std::string name;
1000 if (fDataSource) {
1001 name = fDataSource->GetLabel();
1002 } else if (fTree) {
1003 name = fTree->GetName();
1004 } else {
1005 name = "Empty source\\nEntries: " + std::to_string(fNEmptyEntries);
1006 }
1007 auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(
1009 visitedMap[(void *)this] = thisNode;
1010 return thisNode;
1011}
1012
1013////////////////////////////////////////////////////////////////////////////
1014/// Return all valid TTree::Branch names (caching results for subsequent calls).
1015/// Never use fBranchNames directy, always request it through this method.
1017{
1018 if (fValidBranchNames.empty() && fTree) {
1019 fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
1020 }
1021 return fValidBranchNames;
1022}
1023
1024/// Return true if AddDataSourceColumnReaders was called for column name col.
1025bool RLoopManager::HasDataSourceColumnReaders(const std::string &col, const std::type_info &ti) const
1026{
1027 const auto key = MakeDatasetColReadersKey(col, ti);
1028 assert(fDataSource != nullptr);
1029 // since data source column readers are always added for all slots at the same time,
1030 // if the reader is present for slot 0 we have it for all other slots as well.
1031 return fDatasetColumnReaders[0].find(key) != fDatasetColumnReaders[0].end();
1032}
1033
1035 std::vector<std::unique_ptr<RColumnReaderBase>> &&readers,
1036 const std::type_info &ti)
1037{
1038 const auto key = MakeDatasetColReadersKey(col, ti);
1039 assert(fDataSource != nullptr && !HasDataSourceColumnReaders(col, ti));
1040 assert(readers.size() == fNSlots);
1041
1042 for (auto slot = 0u; slot < fNSlots; ++slot) {
1043 fDatasetColumnReaders[slot][key] = std::move(readers[slot]);
1044 }
1045}
1046
1047// Differently from AddDataSourceColumnReaders, this can be called from multiple threads concurrently
1048/// \brief Register a new RTreeColumnReader with this RLoopManager.
1049/// \return A shared pointer to the inserted column reader.
1050RColumnReaderBase *RLoopManager::AddTreeColumnReader(unsigned int slot, const std::string &col,
1051 std::unique_ptr<RColumnReaderBase> &&reader,
1052 const std::type_info &ti)
1053{
1054 auto &readers = fDatasetColumnReaders[slot];
1055 const auto key = MakeDatasetColReadersKey(col, ti);
1056 // if a reader for this column and this slot was already there, we are doing something wrong
1057 assert(readers.find(key) == readers.end() || readers[key] == nullptr);
1058 auto *rptr = reader.get();
1059 readers[key] = std::move(reader);
1060 return rptr;
1061}
1062
1064RLoopManager::GetDatasetColumnReader(unsigned int slot, const std::string &col, const std::type_info &ti) const
1065{
1066 const auto key = MakeDatasetColReadersKey(col, ti);
1067 auto it = fDatasetColumnReaders[slot].find(key);
1068 if (it != fDatasetColumnReaders[slot].end())
1069 return it->second.get();
1070 else
1071 return nullptr;
1072}
1073
1075{
1076 if (callback)
1077 fSampleCallbacks.insert({nodePtr, std::move(callback)});
1078}
#define R__LOG_DEBUG(DEBUGLEVEL,...)
Definition RLogger.hxx:365
#define R__LOG_INFO(...)
Definition RLogger.hxx:364
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define e(i)
Definition RSha256.hxx:103
static void cleanup()
long long Long64_t
Definition RtypesCore.h:80
unsigned long long ULong64_t
Definition RtypesCore.h:81
#define R__ASSERT(e)
Definition TError.h:117
const char * filters[]
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
char name[80]
Definition TGX11.cxx:110
R__EXTERN TVirtualMutex * gROOTMutex
Definition TROOT.h:63
#define R__LOCKGUARD(mutex)
The head node of a RDF computation graph.
void UpdateSampleInfo(unsigned int slot, const std::pair< ULong64_t, ULong64_t > &range)
RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
unsigned int fNRuns
Number of event loops run.
bool CheckFilters(unsigned int, Long64_t) final
void EvalChildrenCounts()
Trigger counting of number of children nodes for each node of the functional graph.
void CleanUpNodes()
Perform clean-up operations. To be called at the end of each event loop.
void RunEmptySource()
Run event loop with no source files, in sequence.
void 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 SetTree(std::shared_ptr< TTree > tree)
std::shared_ptr< TTree > fTree
Shared pointer to the input TTree.
std::vector< std::unique_ptr< TTree > > fFriends
Friends of the fTree. Only used if we constructed fTree ourselves.
std::vector< RDefineBase * > fBookedDefines
void RunTreeReader()
Run event loop over one or multiple ROOT files, in sequence.
ROOT::Internal::TreeUtils::RNoCleanupNotifier fNoCleanupNotifier
std::vector< RDFInternal::RActionBase * > fRunActions
Non-owning pointers to actions already run.
RColumnReaderBase * GetDatasetColumnReader(unsigned int slot, const std::string &col, const std::type_info &ti) const
std::vector< RRangeBase * > fBookedRanges
std::vector< ROOT::RDF::Experimental::RSample > fSamples
Samples need to survive throughout the whole event loop, hence stored as an attribute.
std::vector< std::string > ColumnNames_t
void RunAndCheckFilters(unsigned int slot, Long64_t entry)
Execute actions and make sure named filters are called for each event.
std::vector< RFilterBase * > fBookedFilters
void Run(bool jit=true)
Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
std::unordered_map< void *, ROOT::RDF::SampleCallback_t > fSampleCallbacks
Registered callbacks to call at the beginning of each "data block".
std::vector< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
std::vector< RDFInternal::RCallback > fCallbacks
Registered callbacks.
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< 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. Null if no data-source.
RDFInternal::RNewSampleNotifier fNewSampleNotifier
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.
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).
A dataset specification for RDataFrame.
This type represents a sample identifier, to be used in conjunction with RDataFrame features such as ...
This class provides a simple interface to execute the same task multiple times in parallel threads,...
void Foreach(F func, unsigned nTimes, unsigned nChunks=0)
Execute a function without arguments several times in parallel, dividing the execution in nChunks.
A Branch for the case of an object.
static TClass * Class()
A TTree is a list of TBranches.
Definition TBranch.h:89
static TClass * Class()
TClass * IsA() const override
Definition TBranch.h:291
TObjArray * GetListOfLeaves()
Definition TBranch.h:243
A chain is a collection of files containing TTree objects.
Definition TChain.h:33
TObjArray * GetListOfFiles() const
Definition TChain.h:111
A List of entry numbers in a TTree or TChain.
Definition TEntryList.h:26
A TFriendElement TF describes a TTree object TF in a file.
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition TLeaf.h:57
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
Stopwatch class.
Definition TStopwatch.h:28
Double_t RealTime()
Stop the stopwatch (if it is running) and return the realtime (in seconds) passed between the start a...
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Double_t CpuTime()
Stop the stopwatch (if it is running) and return the cputime (in seconds) passed between the start an...
void Stop()
Stop the stopwatch.
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Definition TTreeReader.h:44
@ kEntryBeyondEnd
last entry loop has reached its end
@ kEntryValid
data read okay
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual TBranch * FindBranch(const char *name)
Return the branch that correspond to the path 'branchname', which can include the name of the tree or...
Definition TTree.cxx:4832
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition TTree.cxx:5285
static void SetMaxTreeSize(Long64_t maxsize=100000000000LL)
Set the maximum size in bytes of a Tree file (static function).
Definition TTree.cxx:9178
virtual TObjArray * GetListOfBranches()
Definition TTree.h:485
virtual TTree * GetTree() const
Definition TTree.h:514
virtual TList * GetListOfFriends() const
Definition TTree.h:487
virtual const char * GetFriendAlias(TTree *) const
If the 'tree' is a friend, this method returns its alias name.
Definition TTree.cxx:6023
ROOT::Experimental::RLogChannel & RDFLogChannel()
Definition RDFUtils.cxx:37
std::vector< std::string > GetBranchNames(TTree &t, bool allowDuplicates=true)
Get all the branches names, including the ones of the friend trees.
unsigned int GetNSlots()
Definition RDFUtils.cxx:283
void Erase(const T &that, std::vector< T > &v)
Erase that element from vector v
Definition Utils.hxx:187
Long64_t InterpreterCalc(const std::string &code, const std::string &context="")
Jit code in the interpreter with TInterpreter::Calc, throw in case of errors.
Definition RDFUtils.cxx:327
std::vector< std::string > GetTreeFullPaths(const TTree &tree)
std::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 a RDataFrame computation graph via e....
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
Definition TROOT.cxx:558
Definition file.py:1
Definition tree.py:1
static const char * what
Definition stlLoader.cc:6
A RAII object that calls RLoopManager::CleanUpTask at destruction.
RCallCleanUpTask(RLoopManager &lm, unsigned int arg=0u, TTreeReader *reader=nullptr)
A RAII object to pop and push slot numbers from a RSlotStack object.