Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
mt102_readNtuplesFillHistosAndFit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_multicore
3 /// \notebook
4 /// Read n-tuples in distinct workers, fill histograms, merge them and fit.
5 /// Knowing that other facilities like TProcessExecutor might be more adequate for
6 /// this operation, this tutorial complements mc101, reading and merging.
7 /// We convey another message with this tutorial: the synergy of ROOT and
8 /// STL algorithms is possible.
9 ///
10 /// \macro_output
11 /// \macro_code
12 ///
13 /// \date January 2016
14 /// \author Danilo Piparo
15 
16 Int_t mt102_readNtuplesFillHistosAndFit()
17 {
18 
19  // No nuisance for batch execution
20  gROOT->SetBatch();
21 
22  // Perform the operation sequentially ---------------------------------------
23  TChain inputChain("multiCore");
24  inputChain.Add("mt101_multiCore_*.root");
25  TH1F outHisto("outHisto", "Random Numbers", 128, -4, 4);
26  inputChain.Draw("r >> outHisto");
27  outHisto.Fit("gaus");
28 
29  // We now go MT! ------------------------------------------------------------
30 
31  // The first, fundamental operation to be performed in order to make ROOT
32  // thread-aware.
33  ROOT::EnableThreadSafety();
34 
35  // We adapt our parallelisation to the number of input files
36  const auto nFiles = inputChain.GetListOfFiles()->GetEntries();
37 
38  // We define the histograms we'll fill
39  std::vector<TH1F> histograms;
40  auto workerIDs = ROOT::TSeqI(nFiles);
41  histograms.reserve(nFiles);
42  for (auto workerID : workerIDs) {
43  histograms.emplace_back(TH1F(Form("outHisto_%u", workerID), "Random Numbers", 128, -4, 4));
44  }
45 
46  // We define our work item
47  auto workItem = [&histograms](UInt_t workerID) {
48  TFile f(Form("mt101_multiCore_%u.root", workerID));
49  auto ntuple = f.Get<TNtuple>("multiCore");
50  auto &histo = histograms.at(workerID);
51  for (auto index : ROOT::TSeqL(ntuple->GetEntriesFast())) {
52  ntuple->GetEntry(index);
53  histo.Fill(ntuple->GetArgs()[0]);
54  }
55  };
56 
57  TH1F sumHistogram("SumHisto", "Random Numbers", 128, -4, 4);
58 
59  // Create the collection which will hold the threads, our "pool"
60  std::vector<std::thread> workers;
61 
62  // Spawn workers
63  // Fill the "pool" with workers
64  for (auto workerID : workerIDs) {
65  workers.emplace_back(workItem, workerID);
66  }
67 
68  // Now join them
69  for (auto &&worker : workers)
70  worker.join();
71 
72  // And reduce with a simple lambda
73  std::for_each(std::begin(histograms), std::end(histograms),
74  [&sumHistogram](const TH1F &h) { sumHistogram.Add(&h); });
75 
76  sumHistogram.Fit("gaus", 0);
77 
78  return 0;
79 }