TopQuarkAnalysis
Examples
bin
TopJetFWLiteAnalyzer.cc
Go to the documentation of this file.
1
#include <memory>
2
#include <string>
3
#include <cstdlib>
4
#include <sstream>
5
#include <fstream>
6
#include <iostream>
7
8
#include "
FWCore/FWLite/interface/FWLiteEnabler.h
"
9
#include "
DataFormats/PatCandidates/interface/Jet.h
"
10
11
#include "
TopQuarkAnalysis/Examples/bin/NiceStyle.cc
"
12
#include "
TopQuarkAnalysis/Examples/interface/RootSystem.h
"
13
#include "
TopQuarkAnalysis/Examples/interface/RootHistograms.h
"
14
#include "
TopQuarkAnalysis/Examples/interface/RootPostScript.h
"
15
16
int
main
(
int
argc
,
char
*
argv
[]) {
17
if
(
argc
< 3) {
18
// -------------------------------------------------
19
std::cerr
<<
"ERROR:: "
20
<<
"Wrong number of arguments! Please specify:"
<< std::endl
21
<<
" * filepath"
<< std::endl
22
<<
" * process name"
<< std::endl;
23
// -------------------------------------------------
24
return
-1;
25
}
26
27
// load framework libraries
28
gSystem->Load(
"libFWCoreFWLite"
);
29
FWLiteEnabler::enable
();
30
31
// set nice style for histograms
32
setNiceStyle
();
33
34
// define some histograms
35
TH1I* noJets =
new
TH1I(
"noJets"
,
"N_{Jets}"
, 10, 0, 10);
36
TH1F* ptJets =
new
TH1F(
"ptJets"
,
"pt_{Jets}"
, 100, 0., 300.);
37
TH1F* enJets =
new
TH1F(
"enJets"
,
"energy_{Jets}"
, 100, 0., 300.);
38
TH1F*
etaJets
=
new
TH1F(
"etaJets"
,
"eta_{Jets}"
, 100, -3., 3.);
39
TH1F* phiJets =
new
TH1F(
"phiJets"
,
"phi_{Jets}"
, 100, -5., 5.);
40
41
// -------------------------------------------------
42
std::cout
<<
"open file: "
<<
argv
[1] << std::endl;
43
// -------------------------------------------------
44
TFile* inFile = TFile::Open(
argv
[1]);
45
TTree* events_ =
nullptr
;
46
if
(inFile)
47
inFile->GetObject(
"Events"
, events_);
48
if
(events_ ==
nullptr
) {
49
// -------------------------------------------------
50
std::cerr
<<
"ERROR:: "
51
<<
"Unable to retrieve TTree Events!"
<< std::endl
52
<<
" Eighter wrong file name or the the tree doesn't exists"
<< std::endl;
53
// -------------------------------------------------
54
return
-1;
55
}
56
57
// acess branch of elecs
58
char
jetsName[50];
59
sprintf(jetsName,
"patJets_selectedPatJets__%s.obj"
,
argv
[2]);
60
TBranch* jets_ = events_->GetBranch(jetsName);
61
assert
(jets_ !=
nullptr
);
62
63
// loop over events and fill histograms
64
std::vector<pat::Jet>
jets
;
65
int
nevt
= events_->GetEntries();
66
67
// -------------------------------------------------
68
std::cout
<<
"start looping "
<<
nevt
<<
" events..."
<< std::endl;
69
// -------------------------------------------------
70
for
(
int
evt = 0; evt <
nevt
; ++evt) {
71
// set branch address
72
jets_->SetAddress(&
jets
);
73
// get event
74
jets_->GetEntry(evt);
75
events_->GetEntry(evt, 0);
76
77
// -------------------------------------------------
78
if
(evt > 0 && !(evt % 10))
79
std::cout
<<
" processing event: "
<< evt << std::endl;
80
// -------------------------------------------------
81
82
// fill histograms
83
noJets->Fill(
jets
.size());
84
for
(
unsigned
idx
= 0;
idx
<
jets
.size(); ++
idx
) {
85
// fill histograms
86
ptJets->Fill(
jets
[
idx
].
pt
());
87
enJets->Fill(
jets
[
idx
].
energy
());
88
etaJets
->Fill(
jets
[
idx
].
eta
());
89
phiJets->Fill(
jets
[
idx
].
phi
());
90
}
91
}
92
// -------------------------------------------------
93
std::cout
<<
"close file"
<< std::endl;
94
// -------------------------------------------------
95
inFile->Close();
96
97
// save histograms to file
98
TFile
outFile
(
"analyzeJets.root"
,
"recreate"
);
99
outFile
.mkdir(
"analyzeJet"
);
100
outFile
.cd(
"analyzeJet"
);
101
noJets->Write();
102
ptJets->Write();
103
enJets->Write();
104
etaJets
->Write();
105
phiJets->Write();
106
outFile
.Close();
107
108
// free allocated space
109
delete
noJets;
110
delete
ptJets;
111
delete
enJets;
112
delete
etaJets
;
113
delete
phiJets;
114
115
return
0;
116
}
RootPostScript.h
RootHistograms.h
cmsBatch.argv
argv
Definition:
cmsBatch.py:279
NiceStyle.cc
dir2webdir.argc
argc
Definition:
dir2webdir.py:39
DiDispStaMuonMonitor_cfi.pt
pt
Definition:
DiDispStaMuonMonitor_cfi.py:39
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
cms::cuda::assert
assert(be >=bs)
singleTopDQM_cfi.jets
jets
Definition:
singleTopDQM_cfi.py:42
heavyIonCSV_trainingSettings.idx
idx
Definition:
heavyIonCSV_trainingSettings.py:5
setNiceStyle
void setNiceStyle()
Definition:
NiceStyle.cc:3
PVValHelper::eta
Definition:
PVValidationHelpers.h:69
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition:
HCALHighEnergyHPDFilter_cfi.py:5
L1TdeCSCTF_cfi.outFile
outFile
Definition:
L1TdeCSCTF_cfi.py:5
FWLiteEnabler::enable
static void enable()
enable automatic library loading
Definition:
FWLiteEnabler.cc:46
PVValHelper::phi
Definition:
PVValidationHelpers.h:68
Jet.h
b2gHadronicHLTEventValidation_cfi.etaJets
etaJets
Definition:
b2gHadronicHLTEventValidation_cfi.py:11
RootSystem.h
nevt
int nevt
Definition:
ReggeGribovPartonMCHadronizer.h:66
main
int main(int argc, char *argv[])
Definition:
TopJetFWLiteAnalyzer.cc:16
EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.cerr
cerr
Definition:
EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.py:8
FWLiteEnabler.h
Generated for CMSSW Reference Manual by
1.8.16