Logo ROOT   6.14/05
Reference Guide
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.
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  TNtuple *ntuple = nullptr;
50  f.GetObject("multiCore", ntuple);
51  auto &histo = histograms.at(workerID);
52  for (auto index : ROOT::TSeqL(ntuple->GetEntriesFast())) {
53  ntuple->GetEntry(index);
54  histo.Fill(ntuple->GetArgs()[0]);
55  }
56  };
57 
58  TH1F sumHistogram("SumHisto", "Random Numbers", 128, -4, 4);
59 
60  // Create the collection which will hold the threads, our "pool"
61  std::vector<std::thread> workers;
62 
63  // Spawn workers
64  // Fill the "pool" with workers
65  for (auto workerID : workerIDs) {
66  workers.emplace_back(workItem, workerID);
67  }
68 
69  // Now join them
70  for (auto &&worker : workers)
71  worker.join();
72 
73  // And reduce with a simple lambda
74  std::for_each(std::begin(histograms), std::end(histograms),
75  [&sumHistogram](const TH1F &h) { sumHistogram.Add(&h); });
76 
77  sumHistogram.Fit("gaus", 0);
78 
79  return 0;
80 }
virtual Long64_t GetEntriesFast() const
Definition: TTree.h:386
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:285
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:47
#define gROOT
Definition: TROOT.h:410
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all branches of entry and return total number of bytes read.
Definition: TTree.cxx:5363
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:567
#define f(i)
Definition: RSha256.hxx:104
int Int_t
Definition: RtypesCore.h:41
A simple TTree restricted to a list of float variables only.
Definition: TNtuple.h:28
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
#define h(i)
Definition: RSha256.hxx:106
Float_t * GetArgs() const
Definition: TNtuple.h:56
void EnableThreadSafety()
Enables the global mutex to make ROOT thread safe/aware.
Definition: TROOT.cxx:546
A pseudo container class which is a generator of indices.
Definition: TSeq.hxx:66
A chain is a collection of files containing TTree objects.
Definition: TChain.h:33
TSeq< int > TSeqI
Definition: TSeq.hxx:194