Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
Validation
RecoParticleFlow
plugins
PFJetBenchmarkAnalyzer.cc
Go to the documentation of this file.
1
// -*- C++ -*-
2
//
3
// Package:
4
// Class: PFJetBenchmarkAnalyzer.cc
5
//
14
//
15
// Original Author: Michel Della Negra
16
// Created: Wed Jan 23 10:11:13 CET 2008
17
// $Id: PFJetBenchmarkAnalyzer.cc,v 1.3 2010/02/20 21:02:43 wmtan Exp $
18
// Extensions by Joanna Weng
19
//
20
21
22
// system include files
23
#include <memory>
24
25
// user include files
26
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
27
#include "
FWCore/Framework/interface/EDAnalyzer.h
"
28
29
#include "
FWCore/Framework/interface/Event.h
"
30
#include "
FWCore/Framework/interface/MakerMacros.h
"
31
32
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
33
#include "
DataFormats/JetReco/interface/Jet.h
"
34
#include "
DataFormats/JetReco/interface/CaloJetCollection.h
"
35
#include "
DataFormats/JetReco/interface/BasicJetCollection.h
"
36
#include "
DataFormats/JetReco/interface/PFJet.h
"
37
#include "
DataFormats/JetReco/interface/GenJet.h
"
38
#include "
DataFormats/Candidate/interface/Candidate.h
"
39
#include "
DataFormats/Candidate/interface/CandidateFwd.h
"
40
#include "
RecoParticleFlow/Benchmark/interface/PFJetBenchmark.h
"
41
#include "
FWCore/ServiceRegistry/interface/Service.h
"
42
#include "
FWCore/Utilities/interface/InputTag.h
"
43
using namespace
edm;
44
using namespace
reco
;
45
using namespace
std;
46
47
//
48
// class decleration
49
50
51
class
PFJetBenchmarkAnalyzer
:
public
edm::EDAnalyzer
{
52
public
:
53
explicit
PFJetBenchmarkAnalyzer
(
const
edm::ParameterSet
&);
54
~
PFJetBenchmarkAnalyzer
();
55
56
57
private
:
58
virtual
void
beginJob
() ;
59
virtual
void
analyze
(
const
edm::Event
&,
const
edm::EventSetup
&);
60
virtual
void
endJob() ;
61
// ----------member data ---------------------------
62
63
};
65
66
//neuhaus - comment
67
PFJetBenchmark
PFJetBenchmark_
;
68
InputTag
sGenJetAlgo
;
69
InputTag
sJetAlgo
;
70
string
outjetfilename
;
71
bool
pfjBenchmarkDebug
;
72
bool
plotAgainstReco
;
73
bool
onlyTwoJets
;
74
double
deltaRMax
=0.1;
75
string
benchmarkLabel_
;
76
double
recPt
;
77
double
maxEta
;
78
DQMStore
*
dbe_
;
79
//
80
// constants, enums and typedefs
81
//
82
83
//
84
// static data member definitions
85
//
86
87
//
88
// constructors and destructor
89
//
90
PFJetBenchmarkAnalyzer::PFJetBenchmarkAnalyzer
(
const
edm::ParameterSet
& iConfig)
91
92
{
93
//now do what ever initialization is needed
94
sGenJetAlgo
=
95
iConfig.
getParameter
<
InputTag
>(
"InputTruthLabel"
);
96
sJetAlgo
=
97
iConfig.
getParameter
<
InputTag
>(
"InputRecoLabel"
);
98
outjetfilename
=
99
iConfig.
getUntrackedParameter
<
string
>(
"OutputFile"
);
100
pfjBenchmarkDebug
=
101
iConfig.
getParameter
<
bool
>(
"pfjBenchmarkDebug"
);
102
plotAgainstReco
=
103
iConfig.
getParameter
<
bool
>(
"PlotAgainstRecoQuantities"
);
104
onlyTwoJets
=
105
iConfig.
getParameter
<
bool
>(
"OnlyTwoJets"
);
106
deltaRMax
=
107
iConfig.
getParameter
<
double
>(
"deltaRMax"
);
108
benchmarkLabel_
=
109
iConfig.
getParameter
<
string
>(
"BenchmarkLabel"
);
110
recPt
=
111
iConfig.
getParameter
<
double
>(
"recPt"
);
112
maxEta
=
113
iConfig.
getParameter
<
double
>(
"maxEta"
);
114
115
dbe_
=
edm::Service<DQMStore>
().
operator
->();
116
117
PFJetBenchmark_
.
setup
(
118
outjetfilename,
119
pfjBenchmarkDebug,
120
plotAgainstReco,
121
onlyTwoJets,
122
deltaRMax,
123
benchmarkLabel_,
124
recPt,
125
maxEta,
126
dbe_);
127
}
128
129
130
PFJetBenchmarkAnalyzer::~PFJetBenchmarkAnalyzer
()
131
{
132
// do anything here that needs to be done at desctruction time
133
// (e.g. close files, deallocate resources etc.)
134
}
135
136
137
//
138
// member functions
139
//
140
141
// ------------ method called to for each event ------------
142
void
143
PFJetBenchmarkAnalyzer::analyze
(
const
edm::Event
&
iEvent
,
144
const
edm::EventSetup
& iSetup)
145
{
146
// get gen jet collection
147
Handle<GenJetCollection>
genjets;
148
bool
isGen = iEvent.
getByLabel
(
sGenJetAlgo
, genjets);
149
if
(!isGen) {
150
std::cout
<<
"Warning : no Gen jets in input !"
<< std::endl;
151
return
;
152
}
153
154
// get rec PFJet collection
155
Handle<PFJetCollection>
pfjets;
156
bool
isReco = iEvent.
getByLabel
(
sJetAlgo
, pfjets);
157
if
(!isReco) {
158
std::cout
<<
"Warning : no PF jets in input !"
<< std::endl;
159
return
;
160
}
161
// Analyse (no "z" in "analyse" : we are in Europe, dammit!)
162
PFJetBenchmark_
.
process
(*pfjets, *genjets);
163
}
164
165
166
// ------------ method called once each job just before starting event loop ------------
167
void
168
PFJetBenchmarkAnalyzer::beginJob
()
169
{
170
171
}
172
173
// ------------ method called once each job just after ending the event loop ------------
174
void
175
PFJetBenchmarkAnalyzer::endJob
() {
176
// PFJetBenchmark_.save();
177
PFJetBenchmark_
.
write
();
178
}
179
180
//define this as a plug-in
181
DEFINE_FWK_MODULE
(
PFJetBenchmarkAnalyzer
);
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
BasicJetCollection.h
edm::Service< DQMStore >
onlyTwoJets
bool onlyTwoJets
Definition:
PFJetBenchmarkAnalyzer.cc:73
Vispa.Share.Profiling.analyze
def analyze
Definition:
Profiling.py:11
deltaRMax
double deltaRMax
Definition:
PFJetBenchmarkAnalyzer.cc:74
cppFunctionSkipper.operator
string operator
Definition:
cppFunctionSkipper.py:10
PFJetBenchmark::write
void write()
Definition:
PFJetBenchmark.cc:46
Event.h
CaloJetCollection.h
MakerMacros.h
PFJetBenchmarkAnalyzer
Definition:
PFJetBenchmarkAnalyzer.cc:51
PFJetBenchmark_
PFJetBenchmark PFJetBenchmark_
PFJet Benchmark.
Definition:
PFJetBenchmarkAnalyzer.cc:67
maxEta
double maxEta
Definition:
PFJetBenchmarkAnalyzer.cc:77
edm::DEFINE_FWK_MODULE
DEFINE_FWK_MODULE(HiMixingModule)
edm::Handle
Definition:
AssociativeIterator.h:48
outjetfilename
string outjetfilename
Definition:
PFJetBenchmarkAnalyzer.cc:70
Frameworkfwd.h
bk::beginJob
void beginJob()
Definition:
Breakpoints.cc:15
ParameterSet.h
PFJetBenchmarkAnalyzer::beginJob
virtual void beginJob()
Definition:
PFJetBenchmarkAnalyzer.cc:168
Candidate.h
DQMStore
Definition:
DQMStore.h:67
iEvent
int iEvent
Definition:
GenABIO.cc:243
PFJetBenchmark::setup
void setup(std::string Filename, bool debug, bool plotAgainstReco=0, bool onlyTwoJets=1, double deltaRMax=0.1, std::string benchmarkLabel_="ParticleFlow", double recPt=-1, double maxEta=-1, DQMStore *dbe_store=NULL)
Definition:
PFJetBenchmark.cc:63
L1Trigger_dataformats.reco
dictionary reco
Definition:
L1Trigger_dataformats.py:9
Service.h
edm::EventSetup
Definition:
EventSetup.h:44
Jet.h
benchmarkLabel_
string benchmarkLabel_
Definition:
PFJetBenchmarkAnalyzer.cc:75
PFJetBenchmarkAnalyzer::analyze
virtual void analyze(const edm::Event &, const edm::EventSetup &)
Definition:
PFJetBenchmarkAnalyzer.cc:143
edm::Event::getByLabel
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition:
Event.h:361
edm::EDAnalyzer
Definition:
EDAnalyzer.h:15
EDAnalyzer.h
dbe_
DQMStore * dbe_
Definition:
PFJetBenchmarkAnalyzer.cc:78
PFJetBenchmark.h
PFJetBenchmark::process
void process(const reco::PFJetCollection &, const reco::GenJetCollection &)
Definition:
PFJetBenchmark.cc:230
PFJetBenchmarkAnalyzer::PFJetBenchmarkAnalyzer
PFJetBenchmarkAnalyzer(const edm::ParameterSet &)
Definition:
PFJetBenchmarkAnalyzer.cc:90
plotAgainstReco
bool plotAgainstReco
Definition:
PFJetBenchmarkAnalyzer.cc:72
sJetAlgo
InputTag sJetAlgo
Definition:
PFJetBenchmarkAnalyzer.cc:69
pfjBenchmarkDebug
bool pfjBenchmarkDebug
Definition:
PFJetBenchmarkAnalyzer.cc:71
sGenJetAlgo
InputTag sGenJetAlgo
Definition:
PFJetBenchmarkAnalyzer.cc:68
PFJetBenchmarkAnalyzer::~PFJetBenchmarkAnalyzer
~PFJetBenchmarkAnalyzer()
Definition:
PFJetBenchmarkAnalyzer.cc:130
edm::InputTag
Definition:
InputTag.h:17
InputTag.h
edm::ParameterSet
Definition:
ParameterSet.h:35
gather_cfg.cout
tuple cout
Definition:
gather_cfg.py:121
CandidateFwd.h
recPt
double recPt
Definition:
PFJetBenchmarkAnalyzer.cc:76
edm::Event
Definition:
Event.h:56
PFJetBenchmarkAnalyzer::endJob
virtual void endJob()
Definition:
PFJetBenchmarkAnalyzer.cc:175
PFJet.h
GenJet.h
PFJetBenchmark
Definition:
PFJetBenchmark.h:38
Generated for CMSSW Reference Manual by
1.8.5