Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RLoopManager.cxx
Go to the documentation of this file.
1 #include "RConfigure.h" // R__USE_IMT
9 #include "RtypesCore.h" // Long64_t
10 #include "TBranchElement.h"
11 #include "TBranchObject.h"
12 #include "TEntryList.h"
13 #include "TError.h"
14 #include "TInterpreter.h"
15 #include "TROOT.h" // IsImplicitMTEnabled
16 #include "TTreeReader.h"
17 
18 #ifdef R__USE_IMT
19 #include "ROOT/TThreadExecutor.hxx"
20 #endif
21 
22 #include <atomic>
23 #include <functional>
24 #include <memory>
25 #include <stdexcept>
26 #include <string>
27 #include <vector>
28 
29 using namespace ROOT::Detail::RDF;
30 using namespace ROOT::Internal::RDF;
31 
32 bool ContainsLeaf(const std::set<TLeaf *> &leaves, TLeaf *leaf)
33 {
34  return (leaves.find(leaf) != leaves.end());
35 }
36 
37 ///////////////////////////////////////////////////////////////////////////////
38 /// This overload does not perform any check on the duplicates.
39 /// It is used for TBranch objects.
40 void UpdateList(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
41  const std::string &friendName)
42 {
43 
44  if (!friendName.empty()) {
45  // In case of a friend tree, users might prepend its name/alias to the branch names
46  const auto friendBName = friendName + "." + branchName;
47  if (bNamesReg.insert(friendBName).second)
48  bNames.push_back(friendBName);
49  }
50 
51  if (bNamesReg.insert(branchName).second)
52  bNames.push_back(branchName);
53 }
54 
55 ///////////////////////////////////////////////////////////////////////////////
56 /// This overloads makes sure that the TLeaf has not been already inserted.
57 void UpdateList(std::set<std::string> &bNamesReg, ColumnNames_t &bNames, const std::string &branchName,
58  const std::string &friendName, std::set<TLeaf *> &foundLeaves, TLeaf *leaf, bool allowDuplicates)
59 {
60  const bool canAdd = allowDuplicates ? true : !ContainsLeaf(foundLeaves, leaf);
61  if (!canAdd) {
62  return;
63  }
64 
65  UpdateList(bNamesReg, bNames, branchName, friendName);
66 
67  foundLeaves.insert(leaf);
68 }
69 
70 void ExploreBranch(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames, TBranch *b, std::string prefix,
71  std::string &friendName)
72 {
73  for (auto sb : *b->GetListOfBranches()) {
74  TBranch *subBranch = static_cast<TBranch *>(sb);
75  auto subBranchName = std::string(subBranch->GetName());
76  auto fullName = prefix + subBranchName;
77 
78  std::string newPrefix;
79  if (!prefix.empty())
80  newPrefix = fullName + ".";
81 
82  ExploreBranch(t, bNamesReg, bNames, subBranch, newPrefix, friendName);
83 
84  if (t.GetBranch(fullName.c_str())) {
85  UpdateList(bNamesReg, bNames, fullName, friendName);
86 
87  } else if (t.GetBranch(subBranchName.c_str())) {
88  UpdateList(bNamesReg, bNames, subBranchName, friendName);
89  }
90  }
91 }
92 
93 void GetBranchNamesImpl(TTree &t, std::set<std::string> &bNamesReg, ColumnNames_t &bNames,
94  std::set<TTree *> &analysedTrees, std::string &friendName, bool allowDuplicates)
95 {
96  std::set<TLeaf *> foundLeaves;
97  if (!analysedTrees.insert(&t).second) {
98  return;
99  }
100 
101  const auto branches = t.GetListOfBranches();
102  // Getting the branches here triggered the read of the first file of the chain if t is a chain.
103  // We check if a tree has been successfully read, otherwise we throw (see ROOT-9984) to avoid further
104  // operations
105  if (!t.GetTree()) {
106  std::string err("GetBranchNames: error in opening the tree ");
107  err += t.GetName();
108  throw std::runtime_error(err);
109  }
110  if (branches) {
111  for (auto b : *branches) {
112  TBranch *branch = static_cast<TBranch *>(b);
113  const auto branchName = std::string(branch->GetName());
114  if (branch->IsA() == TBranch::Class()) {
115  // Leaf list
116  auto listOfLeaves = branch->GetListOfLeaves();
117  if (listOfLeaves->GetEntries() == 1) {
118  auto leaf = static_cast<TLeaf *>(listOfLeaves->At(0));
119  const auto leafName = std::string(leaf->GetName());
120  if (leafName == branchName) {
121  UpdateList(bNamesReg, bNames, branchName, friendName, foundLeaves, leaf, allowDuplicates);
122  }
123  }
124 
125  for (auto leaf : *listOfLeaves) {
126  auto castLeaf = static_cast<TLeaf *>(leaf);
127  const auto leafName = std::string(leaf->GetName());
128  const auto fullName = branchName + "." + leafName;
129  UpdateList(bNamesReg, bNames, fullName, friendName, foundLeaves, castLeaf, allowDuplicates);
130  }
131  } else if (branch->IsA() == TBranchObject::Class()) {
132  // TBranchObject
133  ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName);
134  UpdateList(bNamesReg, bNames, branchName, friendName);
135  } else {
136  // TBranchElement
137  // Check if there is explicit or implicit dot in the name
138 
139  bool dotIsImplied = false;
140  auto be = dynamic_cast<TBranchElement *>(b);
141  if (!be)
142  throw std::runtime_error("GetBranchNames: unsupported branch type");
143  // TClonesArray (3) and STL collection (4)
144  if (be->GetType() == 3 || be->GetType() == 4)
145  dotIsImplied = true;
146 
147  if (dotIsImplied || branchName.back() == '.')
148  ExploreBranch(t, bNamesReg, bNames, branch, "", friendName);
149  else
150  ExploreBranch(t, bNamesReg, bNames, branch, branchName + ".", friendName);
151 
152  UpdateList(bNamesReg, bNames, branchName, friendName);
153  }
154  }
155  }
156 
157  auto friendTrees = t.GetListOfFriends();
158 
159  if (!friendTrees)
160  return;
161 
162  for (auto friendTreeObj : *friendTrees) {
163  auto friendTree = ((TFriendElement *)friendTreeObj)->GetTree();
164 
165  std::string frName;
166  auto alias = t.GetFriendAlias(friendTree);
167  if (alias != nullptr)
168  frName = std::string(alias);
169  else
170  frName = std::string(friendTree->GetName());
171 
172  GetBranchNamesImpl(*friendTree, bNamesReg, bNames, analysedTrees, frName, allowDuplicates);
173  }
174 }
175 
176 ///////////////////////////////////////////////////////////////////////////////
177 /// Get all the branches names, including the ones of the friend trees
178 ColumnNames_t ROOT::Internal::RDF::GetBranchNames(TTree &t, bool allowDuplicates)
179 {
180  std::set<std::string> bNamesSet;
181  ColumnNames_t bNames;
182  std::set<TTree *> analysedTrees;
183  std::string emptyFrName = "";
184  GetBranchNamesImpl(t, bNamesSet, bNames, analysedTrees, emptyFrName, allowDuplicates);
185  return bNames;
186 }
187 
188 
189 RLoopManager::RLoopManager(TTree *tree, const ColumnNames_t &defaultBranches)
190  : fTree(std::shared_ptr<TTree>(tree, [](TTree *) {})), fDefaultColumns(defaultBranches),
191  fNSlots(RDFInternal::GetNSlots()),
192  fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kROOTFilesMT : ELoopType::kROOTFiles)
193 {
194 }
195 
196 RLoopManager::RLoopManager(ULong64_t nEmptyEntries)
197  : fNEmptyEntries(nEmptyEntries), fNSlots(RDFInternal::GetNSlots()),
198  fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kNoFilesMT : ELoopType::kNoFiles)
199 {
200 }
201 
202 RLoopManager::RLoopManager(std::unique_ptr<RDataSource> ds, const ColumnNames_t &defaultBranches)
203  : fDefaultColumns(defaultBranches), fNSlots(RDFInternal::GetNSlots()),
204  fLoopType(ROOT::IsImplicitMTEnabled() ? ELoopType::kDataSourceMT : ELoopType::kDataSource),
205  fDataSource(std::move(ds))
206 {
207  fDataSource->SetNSlots(fNSlots);
208 }
209 
210 // ROOT-9559: we cannot handle indexed friends
211 void RLoopManager::CheckIndexedFriends()
212 {
213  auto friends = fTree->GetListOfFriends();
214  if (!friends)
215  return;
216  for (auto friendElObj : *friends) {
217  auto friendEl = static_cast<TFriendElement *>(friendElObj);
218  auto friendTree = friendEl->GetTree();
219  if (friendTree && friendTree->GetTreeIndex()) {
220  std::string err = fTree->GetName();
221  err += " has a friend, \"";
222  err += friendTree->GetName();
223  err += "\", which has an index. This is not supported.";
224  throw std::runtime_error(err);
225  }
226  }
227 }
228 
229 /// Run event loop with no source files, in parallel.
230 void RLoopManager::RunEmptySourceMT()
231 {
232 #ifdef R__USE_IMT
233  RSlotStack slotStack(fNSlots);
234  // Working with an empty tree.
235  // Evenly partition the entries according to fNSlots. Produce around 2 tasks per slot.
236  const auto nEntriesPerSlot = fNEmptyEntries / (fNSlots * 2);
237  auto remainder = fNEmptyEntries % (fNSlots * 2);
238  std::vector<std::pair<ULong64_t, ULong64_t>> entryRanges;
239  ULong64_t start = 0;
240  while (start < fNEmptyEntries) {
241  ULong64_t end = start + nEntriesPerSlot;
242  if (remainder > 0) {
243  ++end;
244  --remainder;
245  }
246  entryRanges.emplace_back(start, end);
247  start = end;
248  }
249 
250  // Each task will generate a subrange of entries
251  auto genFunction = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
252  auto slot = slotStack.GetSlot();
253  InitNodeSlots(nullptr, slot);
254  for (auto currEntry = range.first; currEntry < range.second; ++currEntry) {
255  RunAndCheckFilters(slot, currEntry);
256  }
257  CleanUpTask(slot);
258  slotStack.ReturnSlot(slot);
259  };
260 
261  ROOT::TThreadExecutor pool;
262  pool.Foreach(genFunction, entryRanges);
263 
264 #endif // not implemented otherwise
265 }
266 
267 /// Run event loop with no source files, in sequence.
268 void RLoopManager::RunEmptySource()
269 {
270  InitNodeSlots(nullptr, 0);
271  for (ULong64_t currEntry = 0; currEntry < fNEmptyEntries && fNStopsReceived < fNChildren; ++currEntry) {
272  RunAndCheckFilters(0, currEntry);
273  }
274  CleanUpTask(0u);
275 }
276 
277 /// Run event loop over one or multiple ROOT files, in parallel.
278 void RLoopManager::RunTreeProcessorMT()
279 {
280 #ifdef R__USE_IMT
281  CheckIndexedFriends();
282  RSlotStack slotStack(fNSlots);
283  const auto &entryList = fTree->GetEntryList() ? *fTree->GetEntryList() : TEntryList();
284  auto tp = std::make_unique<ROOT::TTreeProcessorMT>(*fTree, entryList);
285 
286  std::atomic<ULong64_t> entryCount(0ull);
287 
288  tp->Process([this, &slotStack, &entryCount](TTreeReader &r) -> void {
289  auto slot = slotStack.GetSlot();
290  InitNodeSlots(&r, slot);
291  const auto entryRange = r.GetEntriesRange(); // we trust TTreeProcessorMT to call SetEntriesRange
292  const auto nEntries = entryRange.second - entryRange.first;
293  auto count = entryCount.fetch_add(nEntries);
294  // recursive call to check filters and conditionally execute actions
295  while (r.Next()) {
296  RunAndCheckFilters(slot, count++);
297  }
298  CleanUpTask(slot);
299  slotStack.ReturnSlot(slot);
300  });
301 #endif // no-op otherwise (will not be called)
302 }
303 
304 /// Run event loop over one or multiple ROOT files, in sequence.
305 void RLoopManager::RunTreeReader()
306 {
307  CheckIndexedFriends();
308  TTreeReader r(fTree.get(), fTree->GetEntryList());
309  if (0 == fTree->GetEntriesFast())
310  return;
311  InitNodeSlots(&r, 0);
312 
313  // recursive call to check filters and conditionally execute actions
314  // in the non-MT case processing can be stopped early by ranges, hence the check on fNStopsReceived
315  while (r.Next() && fNStopsReceived < fNChildren) {
316  RunAndCheckFilters(0, r.GetCurrentEntry());
317  }
318  CleanUpTask(0u);
319 }
320 
321 /// Run event loop over data accessed through a DataSource, in sequence.
322 void RLoopManager::RunDataSource()
323 {
324  R__ASSERT(fDataSource != nullptr);
325  fDataSource->Initialise();
326  auto ranges = fDataSource->GetEntryRanges();
327  while (!ranges.empty()) {
328  InitNodeSlots(nullptr, 0u);
329  fDataSource->InitSlot(0u, 0ull);
330  for (const auto &range : ranges) {
331  auto end = range.second;
332  for (auto entry = range.first; entry < end; ++entry) {
333  if (fDataSource->SetEntry(0u, entry)) {
334  RunAndCheckFilters(0u, entry);
335  }
336  }
337  }
338  CleanUpTask(0u);
339  fDataSource->FinaliseSlot(0u);
340  ranges = fDataSource->GetEntryRanges();
341  }
342  fDataSource->Finalise();
343 }
344 
345 /// Run event loop over data accessed through a DataSource, in parallel.
346 void RLoopManager::RunDataSourceMT()
347 {
348 #ifdef R__USE_IMT
349  R__ASSERT(fDataSource != nullptr);
350  RSlotStack slotStack(fNSlots);
351  ROOT::TThreadExecutor pool;
352 
353  // Each task works on a subrange of entries
354  auto runOnRange = [this, &slotStack](const std::pair<ULong64_t, ULong64_t> &range) {
355  const auto slot = slotStack.GetSlot();
356  InitNodeSlots(nullptr, slot);
357  fDataSource->InitSlot(slot, range.first);
358  const auto end = range.second;
359  for (auto entry = range.first; entry < end; ++entry) {
360  if (fDataSource->SetEntry(slot, entry)) {
361  RunAndCheckFilters(slot, entry);
362  }
363  }
364  CleanUpTask(slot);
365  fDataSource->FinaliseSlot(slot);
366  slotStack.ReturnSlot(slot);
367  };
368 
369  fDataSource->Initialise();
370  auto ranges = fDataSource->GetEntryRanges();
371  while (!ranges.empty()) {
372  pool.Foreach(runOnRange, ranges);
373  ranges = fDataSource->GetEntryRanges();
374  }
375  fDataSource->Finalise();
376 #endif // not implemented otherwise (never called)
377 }
378 
379 /// Execute actions and make sure named filters are called for each event.
380 /// Named filters must be called even if the analysis logic would not require it, lest they report confusing results.
381 void RLoopManager::RunAndCheckFilters(unsigned int slot, Long64_t entry)
382 {
383  for (auto &actionPtr : fBookedActions)
384  actionPtr->Run(slot, entry);
385  for (auto &namedFilterPtr : fBookedNamedFilters)
386  namedFilterPtr->CheckFilters(slot, entry);
387  for (auto &callback : fCallbacks)
388  callback(slot);
389 }
390 
391 /// Build TTreeReaderValues for all nodes
392 /// This method loops over all filters, actions and other booked objects and
393 /// calls their `InitRDFValues` methods. It is called once per node per slot, before
394 /// running the event loop. It also informs each node of the TTreeReader that
395 /// a particular slot will be using.
396 void RLoopManager::InitNodeSlots(TTreeReader *r, unsigned int slot)
397 {
398  for (auto &ptr : fBookedActions)
399  ptr->InitSlot(r, slot);
400  for (auto &ptr : fBookedFilters)
401  ptr->InitSlot(r, slot);
402  for (auto &callback : fCallbacksOnce)
403  callback(slot);
404 }
405 
406 /// Initialize all nodes of the functional graph before running the event loop.
407 /// This method is called once per event-loop and performs generic initialization
408 /// operations that do not depend on the specific processing slot (i.e. operations
409 /// that are common for all threads).
410 void RLoopManager::InitNodes()
411 {
412  EvalChildrenCounts();
413  for (auto column : fCustomColumns)
414  column->InitNode();
415  for (auto &filter : fBookedFilters)
416  filter->InitNode();
417  for (auto &range : fBookedRanges)
418  range->InitNode();
419  for (auto &ptr : fBookedActions)
420  ptr->Initialize();
421 }
422 
423 /// Perform clean-up operations. To be called at the end of each event loop.
424 void RLoopManager::CleanUpNodes()
425 {
426  fMustRunNamedFilters = false;
427 
428  // forget RActions and detach TResultProxies
429  for (auto &ptr : fBookedActions)
430  ptr->Finalize();
431 
432  fRunActions.insert(fRunActions.begin(), fBookedActions.begin(), fBookedActions.end());
433  fBookedActions.clear();
434 
435  // reset children counts
436  fNChildren = 0;
437  fNStopsReceived = 0;
438  for (auto &ptr : fBookedFilters)
439  ptr->ResetChildrenCount();
440  for (auto &ptr : fBookedRanges)
441  ptr->ResetChildrenCount();
442 
443  fCallbacks.clear();
444  fCallbacksOnce.clear();
445 }
446 
447 /// Perform clean-up operations. To be called at the end of each task execution.
448 void RLoopManager::CleanUpTask(unsigned int slot)
449 {
450  for (auto &ptr : fBookedActions)
451  ptr->FinalizeSlot(slot);
452  for (auto &ptr : fBookedFilters)
453  ptr->ClearTask(slot);
454 }
455 
456 /// Declare to the interpreter type aliases and other entities required by RDF jitted nodes.
457 /// This method clears the `fToJitDeclare` member variable.
458 void RLoopManager::JitDeclarations()
459 {
460  if (fToJitDeclare.empty())
461  return;
462 
463  RDFInternal::InterpreterDeclare(fToJitDeclare);
464  fToJitDeclare.clear();
465 }
466 
467 /// Add RDF nodes that require just-in-time compilation to the computation graph.
468 /// This method also invokes JitDeclarations() if needed, and clears the `fToJitExec` member variable.
469 void RLoopManager::Jit()
470 {
471  if (fToJitExec.empty())
472  return;
473 
474  JitDeclarations();
475  RDFInternal::InterpreterCalc(fToJitExec, "RLoopManager::Run");
476  fToJitExec.clear();
477 }
478 
479 /// Trigger counting of number of children nodes for each node of the functional graph.
480 /// This is done once before starting the event loop. Each action sends an `increase children count` signal
481 /// upstream, which is propagated until RLoopManager. Each time a node receives the signal, in increments its
482 /// children counter. Each node only propagates the signal once, even if it receives it multiple times.
483 /// Named filters also send an `increase children count` signal, just like actions, as they always execute during
484 /// the event loop so the graph branch they belong to must count as active even if it does not end in an action.
485 void RLoopManager::EvalChildrenCounts()
486 {
487  for (auto &actionPtr : fBookedActions)
488  actionPtr->TriggerChildrenCount();
489  for (auto &namedFilterPtr : fBookedNamedFilters)
490  namedFilterPtr->TriggerChildrenCount();
491 }
492 
493 unsigned int RLoopManager::GetNextID()
494 {
495  static unsigned int id = 0;
496  ++id;
497  return id;
498 }
499 
500 /// Start the event loop with a different mechanism depending on IMT/no IMT, data source/no data source.
501 /// Also perform a few setup and clean-up operations (jit actions if necessary, clear booked actions after the loop...).
502 void RLoopManager::Run()
503 {
504  Jit();
505 
506  InitNodes();
507 
508  switch (fLoopType) {
509  case ELoopType::kNoFilesMT: RunEmptySourceMT(); break;
510  case ELoopType::kROOTFilesMT: RunTreeProcessorMT(); break;
511  case ELoopType::kDataSourceMT: RunDataSourceMT(); break;
512  case ELoopType::kNoFiles: RunEmptySource(); break;
513  case ELoopType::kROOTFiles: RunTreeReader(); break;
514  case ELoopType::kDataSource: RunDataSource(); break;
515  }
516 
517  CleanUpNodes();
518 }
519 
520 /// Return the list of default columns -- empty if none was provided when constructing the RDataFrame
521 const ColumnNames_t &RLoopManager::GetDefaultColumnNames() const
522 {
523  return fDefaultColumns;
524 }
525 
526 TTree *RLoopManager::GetTree() const
527 {
528  return fTree.get();
529 }
530 
531 void RLoopManager::Book(RDFInternal::RActionBase *actionPtr)
532 {
533  fBookedActions.emplace_back(actionPtr);
534 }
535 
536 void RLoopManager::Deregister(RDFInternal::RActionBase *actionPtr)
537 {
538  RDFInternal::Erase(actionPtr, fRunActions);
539  RDFInternal::Erase(actionPtr, fBookedActions);
540 }
541 
542 void RLoopManager::Book(RFilterBase *filterPtr)
543 {
544  fBookedFilters.emplace_back(filterPtr);
545  if (filterPtr->HasName()) {
546  fBookedNamedFilters.emplace_back(filterPtr);
547  fMustRunNamedFilters = true;
548  }
549 }
550 
551 void RLoopManager::Deregister(RFilterBase *filterPtr)
552 {
553  RDFInternal::Erase(filterPtr, fBookedFilters);
554  RDFInternal::Erase(filterPtr, fBookedNamedFilters);
555 }
556 
557 void RLoopManager::Book(RRangeBase *rangePtr)
558 {
559  fBookedRanges.emplace_back(rangePtr);
560 }
561 
562 void RLoopManager::Deregister(RRangeBase *rangePtr)
563 {
564  RDFInternal::Erase(rangePtr, fBookedRanges);
565 }
566 
567 // dummy call, end of recursive chain of calls
568 bool RLoopManager::CheckFilters(unsigned int, Long64_t)
569 {
570  return true;
571 }
572 
573 /// Call `FillReport` on all booked filters
574 void RLoopManager::Report(ROOT::RDF::RCutFlowReport &rep) const
575 {
576  for (const auto &fPtr : fBookedNamedFilters)
577  fPtr->FillReport(rep);
578 }
579 
580 void RLoopManager::RegisterCallback(ULong64_t everyNEvents, std::function<void(unsigned int)> &&f)
581 {
582  if (everyNEvents == 0ull)
583  fCallbacksOnce.emplace_back(std::move(f), fNSlots);
584  else
585  fCallbacks.emplace_back(everyNEvents, std::move(f), fNSlots);
586 }
587 
588 std::vector<std::string> RLoopManager::GetFiltersNames()
589 {
590  std::vector<std::string> filters;
591  for (auto &filter : fBookedFilters) {
592  auto name = (filter->HasName() ? filter->GetName() : "Unnamed Filter");
593  filters.push_back(name);
594  }
595  return filters;
596 }
597 
598 std::vector<RDFInternal::RActionBase *> RLoopManager::GetAllActions()
599 {
600  std::vector<RDFInternal::RActionBase *> actions;
601  actions.insert(actions.begin(), fBookedActions.begin(), fBookedActions.end());
602  actions.insert(actions.begin(), fRunActions.begin(), fRunActions.end());
603  return actions;
604 }
605 
606 std::shared_ptr<ROOT::Internal::RDF::GraphDrawing::GraphNode> RLoopManager::GetGraph()
607 {
608  std::string name;
609  if (fDataSource) {
610  name = fDataSource->GetLabel();
611  } else if (fTree) {
612  name = fTree->GetName();
613  } else {
614  name = std::to_string(fNEmptyEntries);
615  }
616 
617  auto thisNode = std::make_shared<ROOT::Internal::RDF::GraphDrawing::GraphNode>(name);
618  thisNode->SetRoot();
619  thisNode->SetCounter(0);
620  return thisNode;
621 }
622 
623 ////////////////////////////////////////////////////////////////////////////
624 /// Return all valid TTree::Branch names (caching results for subsequent calls).
625 /// Never use fBranchNames directy, always request it through this method.
626 const ColumnNames_t &RLoopManager::GetBranchNames()
627 {
628  if (fValidBranchNames.empty() && fTree) {
629  fValidBranchNames = RDFInternal::GetBranchNames(*fTree, /*allowRepetitions=*/true);
630  }
631  return fValidBranchNames;
632 }