34   return (leaves.find(leaf) != leaves.end());
 
   41                std::string &friendName)
 
   44   if (!friendName.empty()) {
 
   46      auto friendBName = friendName + 
"." + 
branchName;
 
   47      if (bNamesReg.insert(friendBName).second)
 
   48         bNames.push_back(friendBName);
 
   58                std::string &friendName, std::set<TLeaf *> &foundLeaves, 
TLeaf *leaf, 
bool allowDuplicates)
 
   60   const bool canAdd = allowDuplicates ? true : !
ContainsLeaf(foundLeaves, leaf);
 
   67   foundLeaves.insert(leaf);
 
   71                   std::string &friendName)
 
   73   for (
auto sb : *
b->GetListOfBranches()) {
 
   75      auto subBranchName = std::string(subBranch->
GetName());
 
   76      auto fullName = prefix + subBranchName;
 
   78      std::string newPrefix;
 
   80         newPrefix = fullName + 
".";
 
   82      ExploreBranch(t, bNamesReg, bNames, subBranch, newPrefix, friendName);
 
   84      if (t.GetBranch(fullName.c_str())) {
 
   85         UpdateList(bNamesReg, bNames, fullName, friendName);
 
   87      } 
else if (t.GetBranch(subBranchName.c_str())) {
 
   88         UpdateList(bNamesReg, bNames, subBranchName, friendName);
 
   94                        std::set<TTree *> &analysedTrees, std::string &friendName, 
bool allowDuplicates)
 
   96   std::set<TLeaf *> foundLeaves;
 
   97   if (!analysedTrees.insert(&t).second) {
 
  101   const auto branches = t.GetListOfBranches();
 
  103      std::string prefix = 
"";
 
  110            if (listOfLeaves->GetEntries() == 1) {
 
  111               auto leaf = 
static_cast<TLeaf *
>(listOfLeaves->At(0));
 
  112               auto leafName = std::string(leaf->GetName());
 
  118            for (
auto leaf : *listOfLeaves) {
 
  119               auto castLeaf = 
static_cast<TLeaf *
>(leaf);
 
  120               auto leafName = std::string(leaf->GetName());
 
  122               UpdateList(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
 
  132            bool dotIsImplied = 
false;
 
  135               throw std::runtime_error(
"GetBranchNames: unsupported branch type");
 
  137            if (be->GetType() == 3 || be->GetType() == 4)
 
  150   auto friendTrees = t.GetListOfFriends();
 
  155   for (
auto friendTreeObj : *friendTrees) {
 
  159      auto alias = t.GetFriendAlias(friendTree);
 
  160      if (alias != 
nullptr)
 
  161         frName = std::string(alias);
 
  163         frName = std::string(friendTree->GetName());
 
  165      GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
 
  173   std::set<std::string> bNamesSet;
 
  175   std::set<TTree *> analysedTrees;
 
  176   std::string emptyFrName = 
"";
 
  177   GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
 
  183   : fTree(
std::shared_ptr<TTree>(
tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
 
  198     fDataSource(
std::move(ds))
 
  212   std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
 
  220      entryRanges.emplace_back(start, end);
 
  225   auto genFunction = [
this, &slotStack](
const std::pair<ULong64_t, ULong64_t> &range) {
 
  226      auto slot = slotStack.
GetSlot();
 
  228      for (
auto currEntry = range.first; currEntry < range.second; ++currEntry) {
 
  236   pool.
Foreach(genFunction, entryRanges);
 
  255   auto tp = std::make_unique<ROOT::TTreeProcessorMT>(*
fTree);
 
  257   std::atomic<ULong64_t> entryCount(0ull);
 
  259   tp->Process([
this, &slotStack, &entryCount](
TTreeReader &
r) -> 
void {
 
  260      auto slot = slotStack.
GetSlot();
 
  262      const auto entryRange = 
r.GetEntriesRange(); 
 
  263      const auto nEntries = entryRange.second - entryRange.first;
 
  264      auto count = entryCount.fetch_add(nEntries);
 
  279   if (0 == 
fTree->GetEntriesFast())
 
  296   while (!ranges.empty()) {
 
  299      for (
const auto &range : ranges) {
 
  300         auto end = range.second;
 
  301         for (
auto entry = range.first; entry < end; ++entry) {
 
  322   auto runOnRange = [
this, &slotStack](
const std::pair<ULong64_t, ULong64_t> &range) {
 
  323      const auto slot = slotStack.
GetSlot();
 
  326      const auto end = range.second;
 
  327      for (
auto entry = range.first; entry < end; ++entry) {
 
  339   while (!ranges.empty()) {
 
  340      pool.
Foreach(runOnRange, ranges);
 
  352      actionPtr->Run(slot, entry);
 
  354      namedFilterPtr->CheckFilters(slot, entry);
 
  367      ptr->InitSlot(
r, slot);
 
  369      ptr->InitSlot(
r, slot);
 
  407      ptr->ResetChildrenCount();
 
  409      ptr->ResetChildrenCount();
 
  419      ptr->FinalizeSlot(slot);
 
  421      ptr->ClearTask(slot);
 
  427   auto error = TInterpreter::EErrorCode::kNoError;
 
  429   if (TInterpreter::EErrorCode::kNoError != error) {
 
  430      std::string exceptionText =
 
  431         "An error occurred while jitting. The lines above might indicate the cause of the crash\n";
 
  432      throw std::runtime_error(exceptionText.c_str());
 
  446      actionPtr->TriggerChildrenCount();
 
  448      namedFilterPtr->TriggerChildrenCount();
 
  453   static unsigned int id = 0;
 
  536      fPtr->FillReport(rep);
 
  541   if (everyNEvents == 0ull)
 
  549   std::vector<std::string> 
filters;
 
  551      auto name = (filter->HasName() ? filter->GetName() : 
"Unnamed Filter");
 
  559   std::vector<RDFInternal::RActionBase *> actions;
 
  576   auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(
name);
 
  578   thisNode->SetCounter(0);
 
void GetBranchNamesImpl(TTree &t, std::set< std::string > &bNamesReg, ColumnNames_t &bNames, std::set< TTree * > &analysedTrees, std::string &friendName, bool allowDuplicates)
void UpdateList(std::set< std::string > &bNamesReg, ColumnNames_t &bNames, std::string &branchName, std::string &friendName)
This overload does not perform any check on the duplicates.
void ExploreBranch(TTree &t, std::set< std::string > &bNamesReg, ColumnNames_t &bNames, TBranch *b, std::string prefix, std::string &friendName)
bool ContainsLeaf(const std::set< TLeaf * > &leaves, TLeaf *leaf)
unsigned long long ULong64_t
RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
bool CheckFilters(unsigned int, Long64_t) final
void EvalChildrenCounts()
Trigger counting of number of children nodes for each node of the functional graph.
void CleanUpNodes()
Perform clean-up operations. To be called at the end of each event loop.
void RunEmptySource()
Run event loop with no source files, in sequence.
void Report(ROOT::RDF::RCutFlowReport &rep) const final
Call FillReport on all booked filters.
std::vector< RFilterBase * > fBookedNamedFilters
Contains a subset of fBookedFilters, i.e. only the named filters.
void RunEmptySourceMT()
Run event loop with no source files, in parallel.
const ColumnNames_t & GetBranchNames()
Return all valid TTree::Branch names (caching results for subsequent calls).
bool fMustRunNamedFilters
std::shared_ptr< TTree > fTree
Shared pointer to the input TTree.
const ULong64_t fNEmptyEntries
void RunTreeReader()
Run event loop over one or multiple ROOT files, in sequence.
std::vector< RDFInternal::RActionBase * > GetAllActions()
For all the actions, either booked or run.
void CleanUpTask(unsigned int slot)
Perform clean-up operations. To be called at the end of each task execution.
std::vector< RDFInternal::RActionBase * > fRunActions
Non-owning pointers to actions already run.
void Run()
Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
std::vector< RRangeBase * > fBookedRanges
std::vector< TCallback > fCallbacks
Registered callbacks.
void RunAndCheckFilters(unsigned int slot, Long64_t entry)
Execute actions and make sure named filters are called for each event.
std::vector< RFilterBase * > fBookedFilters
std::vector< RDFInternal::RActionBase * > fBookedActions
Non-owning pointers to actions to be run.
std::shared_ptr< ROOT::Internal::RDF::GraphDrawing::GraphNode > GetGraph()
const ELoopType fLoopType
The kind of event loop that is going to be run (e.g. on ROOT files, on no files)
ColumnNames_t fValidBranchNames
Cache of the tree/chain branch names. Never access directy, always use GetBranchNames().
std::vector< RCustomColumnBase * > fCustomColumns
Non-owning container of all custom columns created so far.
const unsigned int fNSlots
const ColumnNames_t & GetDefaultColumnNames() const
Return the list of default columns – empty if none was provided when constructing the RDataFrame.
std::string fToJit
code that should be jitted and executed right before the event loop
void BuildJittedNodes()
Jit all actions that required runtime column type inference, and clean the fToJit member variable.
std::vector< TOneTimeCallback > fCallbacksOnce
Registered callbacks to invoke just once before running the loop.
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.
static unsigned int GetNextID()
const ColumnNames_t fDefaultColumns
void Book(RDFInternal::RActionBase *actionPtr)
void InitNodeSlots(TTreeReader *r, unsigned int slot)
Build TTreeReaderValues for all nodes This method loops over all filters, actions and other booked ob...
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 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.
unsigned int fNStopsReceived
Number of times that a children node signaled to stop processing entries.
unsigned int fNChildren
Number of nodes of the functional graph hanging from this object.
This is an helper class to allow to pick a slot resorting to a map indexed by thread ids.
void ReturnSlot(unsigned int slotNumber)
This class provides a simple interface to execute the same task multiple times in parallel,...
void Foreach(F func, unsigned nTimes, unsigned nChunks=0)
Execute func (with no arguments) nTimes in parallel.
A Branch for the case of an object.
A TTree is a list of TBranches.
TObjArray * GetListOfLeaves()
A TFriendElement TF describes a TTree object TF in a file.
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
virtual const char * GetName() const
Returns name of object.
A simple, robust and fast interface to read values from ROOT colmnar datasets such as TTree,...
ColumnNames_t GetBranchNames(TTree &t, bool allowDuplicates=true)
Get all the branches names, including the ones of the friend trees.
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Namespace for new ROOT classes and functions.
Bool_t IsImplicitMTEnabled()
Returns true if the implicit multi-threading in ROOT is enabled.
ROOT::Detail::RDF::ColumnNames_t ColumnNames_t