Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
imt001_parBranchProcessing.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_multicore
3 /// \notebook
4 /// Demonstrate how to activate and use the implicit parallelisation of TTree::GetEntry.
5 /// Such parallelisation creates one task per top-level branch of the tree being read.
6 /// In this example, most of the branches are floating point numbers, which are very fast to read.
7 /// This parallelisation can be used, though, on bigger trees with many (complex) branches, which
8 /// are more likely to benefit from speedup gains.
9 ///
10 /// \macro_code
11 ///
12 /// \date 26/09/2016
13 /// \author Enric Tejedor
14 
15 int imt001_parBranchProcessing()
16 {
17  // First enable implicit multi-threading globally, so that the implicit parallelisation is on.
18  // The parameter of the call specifies the number of threads to use.
19  int nthreads = 4;
20  ROOT::EnableImplicitMT(nthreads);
21 
22  // Open the file containing the tree
23  auto file = TFile::Open("http://root.cern.ch/files/h1/dstarmb.root");
24 
25  // Get the tree
26  auto tree = file->Get<TTree>("h42");
27 
28  const auto nEntries = tree->GetEntries();
29 
30  // Read the branches in parallel.
31  // Note that the interface does not change, the parallelisation is internal
32  for (auto i : ROOT::TSeqUL(nEntries)) {
33  tree->GetEntry(i); // parallel read
34  }
35 
36  // IMT parallelisation can be disabled for a specific tree
37  tree->SetImplicitMT(false);
38 
39  // If now GetEntry is invoked on the tree, the reading is sequential
40  for (auto i : ROOT::TSeqUL(nEntries)) {
41  tree->GetEntry(i); // sequential read
42  }
43 
44  // Parallel reading can be re-enabled
45  tree->SetImplicitMT(true);
46 
47  // IMT can be also disabled globally.
48  // As a result, no tree will run GetEntry in parallel
49  ROOT::DisableImplicitMT();
50 
51  // This is still sequential: the global flag is disabled, even if the
52  // flag for this particular tree is enabled
53  for (auto i : ROOT::TSeqUL(nEntries)) {
54  tree->GetEntry(i); // sequential read
55  }
56 
57  return 0;
58 }