Logo ROOT   6.07/09
Reference Guide
mt102_readNtuplesFillHistosAndFit.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_multicore
3 /// Read n-tuples in distinct workers, fill histograms, merge them and fit.
4 /// Knowing that other facilities like TProcessExecutor might be more adequate for
5 /// this operation, this tutorial complements mc101, reading and merging.
6 /// We convey another message with this tutorial: the synergy of ROOT and
7 /// STL algorithms is possible.
8 ///
9 /// \macro_output
10 /// \macro_code
11 ///
12 /// \author Danilo Piparo
13 /// \date January 2016
14 
15 Int_t mt102_readNtuplesFillHistosAndFit()
16 {
17 
18  // No nuisance for batch execution
19  gROOT->SetBatch();
20 
21  // Perform the operation sequentially ---------------------------------------
22  TChain inputChain("multiCore");
23  inputChain.Add("mt101_multiCore_*.root");
24  TH1F outHisto("outHisto", "Random Numbers", 128, -4, 4);
25  inputChain.Draw("r >> outHisto");
26  outHisto.Fit("gaus");
27 
28  // We now go MT! ------------------------------------------------------------
29 
30  // The first, fundamental operation to be performed in order to make ROOT
31  // thread-aware.
33 
34  // We adapt our parallelisation to the number of input files
35  const auto nFiles = inputChain.GetListOfFiles()->GetEntries();
36 
37  // We define the histograms we'll fill
38  std::vector<TH1F> histograms;
39  auto workerIDs = ROOT::TSeqI(nFiles);
40  histograms.reserve(nFiles);
41  for (auto workerID : workerIDs){
42  histograms.emplace_back(TH1F(Form("outHisto_%u", workerID), "Random Numbers", 128, -4, 4));
43  }
44 
45  // We define our work item
46  auto workItem = [&histograms](UInt_t workerID) {
47  TFile f(Form("mt101_multiCore_%u.root", workerID));
48  TNtuple *ntuple = nullptr;
49  f.GetObject("multiCore", ntuple);
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) worker.join();
70 
71  // And reduce with a simple lambda
72  std::for_each(std::begin(histograms), std::end(histograms),
73  [&sumHistogram](const TH1F & h) {
74  sumHistogram.Add(&h);
75  });
76 
77  sumHistogram.Fit("gaus",0);
78 
79  return 0;
80 
81 }
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:302
TH1 * h
Definition: legend2.C:5
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:50
#define gROOT
Definition: TROOT.h:364
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:5210
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
int Int_t
Definition: RtypesCore.h:41
Float_t * GetArgs() const
Definition: TNtuple.h:58
A simple TTree restricted to a list of float variables only.
Definition: TNtuple.h:30
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
double f(double x)
void EnableThreadSafety()
Enables the global mutex to make ROOT thread safe/aware.
Definition: TROOT.cxx:498
A pseudo container class which is a generator of indices.
Definition: TSeq.hxx:66
virtual Long64_t GetEntriesFast() const
Definition: TTree.h:394
A chain is a collection of files containg TTree objects.
Definition: TChain.h:35
TSeq< int > TSeqI
Definition: TSeq.hxx:164