Logo ROOT  
Reference Guide
mp105_processEntryList.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_multicore
3/// \notebook -nodraw
4/// Illustrate the usage of the multiproc to process TEntryList with the H1 analysis
5/// example.
6///
7/// \macro_code
8///
9/// \author Gerardo Ganis
10
11#include "TString.h"
12#include "TROOT.h"
13#include "TTree.h"
14#include "TH1F.h"
15#include "TH2F.h"
16#include "TEntryList.h"
17#include "TTreeReader.h"
18#include "TTreeReaderArray.h"
19#include "TTreeReaderValue.h"
20#include "TSystem.h"
21#include "TMath.h"
22#include "TCanvas.h"
23#include "TStyle.h"
24#include "TF1.h"
25#include "TLine.h"
26#include "TPaveStats.h"
27#include "TStopwatch.h"
29
30static std::string tutname = "mp105_processEntryList: ";
31static std::string logfile = "mp105_processEntryList.log";
32static RedirectHandle_t gRH;
33
34std::vector<std::string> files = {"http://root.cern.ch/files/h1/dstarmb.root",
35 "http://root.cern.ch/files/h1/dstarp1a.root",
36 "http://root.cern.ch/files/h1/dstarp1b.root",
37 "http://root.cern.ch/files/h1/dstarp2.root"};
38
39int mp105_processEntryList()
40{
41
42 // MacOSX may generate connection to WindowServer errors
43 gROOT->SetBatch(kTRUE);
44
45 TStopwatch stp;
46
47#include "mp_H1_lambdas.C"
48
50
51 std::cout << tutname << "creating the entry list \n";
52
53 auto sumElist = pool.Process(files, doH1fillList, "h42");
54
55 // Print the entry list
56 if (sumElist) {
57 sumElist->Print();
58 } else {
59 std::cout << tutname << " ERROR creating the entry list \n";
60 return -1;
61 }
62
63 // Time taken
64 stp.Print();
65 stp.Start();
66
67 // Let's analyse H1 with the list
68 std::cout << tutname << "processing the entry list with a lambda \n";
69
70 // Run the analysis
71 auto hListFun = pool.Process(files, doH1useList, *sumElist, "h42");
72
73 // Check the output
74 if (checkH1(hListFun) < 0)
75 return -1;
76
77 // Do the fit
78 if (doFit(hListFun, logfile.c_str()) < 0)
79 return -1;
80
81 stp.Print();
82 stp.Start();
83
84 // Run the analysis with a selector
85 TString selectorPath = gROOT->GetTutorialDir();
86 selectorPath += "/tree/h1analysisTreeReader.C+";
87 std::cout << tutname << "processing the entry list with selector '" << selectorPath << "'\n";
88 auto sel = TSelector::GetSelector(selectorPath);
89
90 // In a second run we use sel
91 sel->SetOption("useList");
92 gSystem->RedirectOutput(logfile.c_str(), "w", &gRH);
93 auto hListSel = pool.Process(files, *sel, *sumElist, "h42");
94 gSystem->RedirectOutput(0, 0, &gRH);
95
96 // Check the output
97 if (checkH1(hListSel) < 0)
98 return -1;
99
100 // Do the fit
101 if (doFit(hListSel, logfile.c_str()) < 0)
102 return -1;
103
104 stp.Print();
105 stp.Start();
106
107 return 0;
108}
const Bool_t kTRUE
Definition: RtypesCore.h:89
#define gROOT
Definition: TROOT.h:406
R__EXTERN TSystem * gSystem
Definition: TSystem.h:556
This class provides an interface to process a TTree dataset in parallel with multi-process technology...
static TSelector * GetSelector(const char *filename)
The code in filename is loaded (interpreted or compiled, see below), filename must contain a valid cl...
Definition: TSelector.cxx:142
Stopwatch class.
Definition: TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
Basic string class.
Definition: TString.h:131
virtual Int_t RedirectOutput(const char *name, const char *mode="a", RedirectHandle_t *h=nullptr)
Redirect standard output (stdout, stderr) to the specified file.
Definition: TSystem.cxx:1708
<a href="https://nbviewer.jupyter.org/url/root.cern/doc/master/notebooks/mp_H1_lambdas....