Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
mp104_processH1.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_legacy
3/// Illustrate the usage of the multiproc to process the H1 analysis
4/// example.
5///
6/// \macro_code
7///
8/// \author Gerardo Ganis
9
10#include "TString.h"
11#include "TROOT.h"
12#include "TTree.h"
13#include "TH1F.h"
14#include "TH2F.h"
15#include "TEntryList.h"
16#include "TTreeReader.h"
17#include "TTreeReaderArray.h"
18#include "TTreeReaderValue.h"
19#include "TSystem.h"
20#include "TMath.h"
21#include "TCanvas.h"
22#include "TStyle.h"
23#include "TF1.h"
24#include "TLine.h"
25#include "TPaveStats.h"
26#include "TStopwatch.h"
28
29static std::string tutname = "mp104_processH1: ";
30static std::string logfile = "mp104_processH1.log";
32
33std::vector<std::string> files {"root://eospublic.cern.ch//eos/root-eos/h1/dstarmb.root",
34 "root://eospublic.cern.ch//eos/root-eos/h1/dstarp1a.root",
35 "root://eospublic.cern.ch//eos/root-eos/h1/dstarp1b.root",
36 "root://eospublic.cern.ch//eos/root-eos/h1/dstarp2.root"};
37
39{
40
41 // MacOSX may generate connection to WindowServer errors
42 gROOT->SetBatch(kTRUE);
43
45
46 // Check and fit lambdas
47 #include "mp_H1_lambdas.C"
48
50
51 std::cout << tutname << "processing the H1 dataset with a lambda \n";
52
53 auto hListFun = pool.Process(files, doH1, "h42");
54
55 // Check the output
56 if (checkH1(hListFun) < 0)
57 return -1;
58
59 // Do the fit
60 if (doFit(hListFun, logfile.c_str()) < 0)
61 return -1;
62
63 stp.Print();
64 stp.Start();
65
66 // Run the analysis with a selector
67
68 TString selectorPath = gROOT->GetTutorialDir();
69 selectorPath += "/analysis/tree/h1analysisTreeReader.C+";
70 std::cout << tutname << "processing the H1 dataset with selector '" << selectorPath << "'\n";
72
73 // In a second run we use sel
74 gSystem->RedirectOutput(logfile.c_str(), "w", &gRH);
75 auto hListSel = pool.Process(files, *sel, "h42");
76 gSystem->RedirectOutput(nullptr, nullptr, &gRH);
77
78 // Check the output
79 if (checkH1(hListSel) < 0)
80 return -1;
81
82 // Do the fit
83 if (doFit(hListSel, logfile.c_str()) < 0)
84 return -1;
85
86 stp.Print();
87 stp.Start();
88
89 return 0;
90}
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t sel
#define gROOT
Definition TROOT.h:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
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...
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:139
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:1715
Lambdas used to check and fit the result of the H1 analysis.