ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
mp102_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 /// We express parallelism with multiprocessing as it is done with multithreading
5 /// in mt102_readNtuplesFillHistosAndFit.
6 ///
7 /// \macro_code
8 ///
9 /// \author Danilo Piparo
10 
11 // Measure time in a scope
12 class TimerRAII {
13  TStopwatch fTimer;
14  std::string fMeta;
15 public:
16  TimerRAII(const char *meta): fMeta(meta) {
17  fTimer.Start();
18  }
19  ~TimerRAII() {
20  fTimer.Stop();
21  std::cout << fMeta << " - real time elapsed " << fTimer.RealTime() << "s" << std::endl;
22  }
23 };
24 
25 Int_t mp102_readNtuplesFillHistosAndFit()
26 {
27 
28  // No nuisance for batch execution
29  gROOT->SetBatch();
30 
31  // Perform the operation sequentially ---------------------------------------
32  TChain inputChain("multiCore");
33  inputChain.Add("mp101_multiCore_*.root");
34  TH1F outHisto("outHisto", "Random Numbers", 128, -4, 4);
35  {
36  TimerRAII t("Sequential read and fit");
37  inputChain.Draw("r >> outHisto");
38  outHisto.Fit("gaus");
39  }
40 
41  // We now go MP! ------------------------------------------------------------
42  // TProcPool offers an interface to directly process trees and chains without
43  // the need for the user to go through the low level implementation of a
44  // map-reduce.
45 
46  // We adapt our parallelisation to the number of input files
47  const auto nFiles = inputChain.GetListOfFiles()->GetEntries();
48 
49 
50  // This is the function invoked during the processing of the trees.
51  auto workItem = [](TTreeReader & reader) {
52  TTreeReaderValue<Float_t> randomRV(reader, "r");
53  auto partialHisto = new TH1F("outHistoMP", "Random Numbers", 128, -4, 4);
54  while (reader.Next()) {
55  partialHisto->Fill(*randomRV);
56  }
57  return partialHisto;
58  };
59 
60  // Create the pool of processes
61  TProcPool workers(nFiles);
62 
63  // Process the TChain
64  {
65  TimerRAII t("Parallel execution");
66  TH1F *sumHistogram = workers.ProcTree(inputChain, workItem, "multiCore");
67  sumHistogram->Fit("gaus", 0);
68  }
69 
70  return 0;
71 
72 }
TTreeReader is a simple, robust and fast interface to read values from a TTree, TChain or TNtuple...
Definition: TTreeReader.h:48
#define gROOT
Definition: TROOT.h:344
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
int Int_t
Definition: RtypesCore.h:41
This class provides a simple interface to execute the same task multiple times in parallel...
Definition: TProcPool.h:38
TThread * t[5]
Definition: threadsh1.C:13
A chain is a collection of files containg TTree objects.
Definition: TChain.h:35
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3607
Stopwatch class.
Definition: TStopwatch.h:30