forked from cms-sw/cmssw
-
Notifications
You must be signed in to change notification settings - Fork 0
/
HLTScoutingPFProducer.cc
263 lines (238 loc) · 10.4 KB
/
HLTScoutingPFProducer.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
// -*- C++ -*-
//
// Package: HLTrigger/JetMET
// Class: HLTScoutingPFProducer
//
/**\class HLTScoutingPFProducer HLTScoutingPFProducer.cc HLTrigger/JetMET/plugins/HLTScoutingPFProducer.cc
Description: Producer for ScoutingPFJets from reco::PFJet objects, ScoutingVertexs from reco::Vertexs and ScoutingParticles from reco::PFCandidates
*/
//
// Original Author: Dustin James Anderson
// Created: Fri, 12 Jun 2015 15:49:20 GMT
//
//
// system include files
#include <memory>
// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/JetReco/interface/PFJet.h"
#include "DataFormats/METReco/interface/MET.h"
#include "DataFormats/METReco/interface/METCollection.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "DataFormats/BTauReco/interface/JetTag.h"
#include "DataFormats/Scouting/interface/ScoutingPFJet.h"
#include "DataFormats/Scouting/interface/ScoutingParticle.h"
#include "DataFormats/Scouting/interface/ScoutingVertex.h"
#include "DataFormats/Math/interface/deltaR.h"
class HLTScoutingPFProducer : public edm::global::EDProducer<> {
public:
explicit HLTScoutingPFProducer(const edm::ParameterSet&);
~HLTScoutingPFProducer();
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
private:
virtual void produce(edm::StreamID sid, edm::Event & iEvent, edm::EventSetup const & setup) const override final;
const edm::EDGetTokenT<reco::PFJetCollection> pfJetCollection_;
const edm::EDGetTokenT<reco::JetTagCollection> pfJetTagCollection_;
const edm::EDGetTokenT<reco::PFCandidateCollection> pfCandidateCollection_;
const edm::EDGetTokenT<reco::VertexCollection> vertexCollection_;
const edm::EDGetTokenT<reco::METCollection> metCollection_;
const edm::EDGetTokenT<double> rho_;
const double pfJetPtCut;
const double pfJetEtaCut;
const double pfCandidatePtCut;
const bool doJetTags;
const bool doCandidates;
const bool doMet;
};
//
// constructors and destructor
//
HLTScoutingPFProducer::HLTScoutingPFProducer(const edm::ParameterSet& iConfig):
pfJetCollection_(consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("pfJetCollection"))),
pfJetTagCollection_(consumes<reco::JetTagCollection>(iConfig.getParameter<edm::InputTag>("pfJetTagCollection"))),
pfCandidateCollection_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfCandidateCollection"))),
vertexCollection_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"))),
metCollection_(consumes<reco::METCollection>(iConfig.getParameter<edm::InputTag>("metCollection"))),
rho_(consumes<double>(iConfig.getParameter<edm::InputTag>("rho"))),
pfJetPtCut(iConfig.getParameter<double>("pfJetPtCut")),
pfJetEtaCut(iConfig.getParameter<double>("pfJetEtaCut")),
pfCandidatePtCut(iConfig.getParameter<double>("pfCandidatePtCut")),
doJetTags(iConfig.getParameter<bool>("doJetTags")),
doCandidates(iConfig.getParameter<bool>("doCandidates")),
doMet(iConfig.getParameter<bool>("doMet"))
{
//register products
produces<ScoutingPFJetCollection>();
produces<ScoutingParticleCollection>();
produces<ScoutingVertexCollection>();
produces<double>("rho");
produces<double>("pfMetPt");
produces<double>("pfMetPhi");
}
HLTScoutingPFProducer::~HLTScoutingPFProducer()
{ }
// ------------ method called to produce the data ------------
void HLTScoutingPFProducer::produce(edm::StreamID sid, edm::Event & iEvent, edm::EventSetup const & setup) const
{
using namespace edm;
//get PF jets
Handle<reco::PFJetCollection> pfJetCollection;
if(!iEvent.getByToken(pfJetCollection_, pfJetCollection)){
edm::LogError ("HLTScoutingPFProducer") << "invalid collection: pfJetCollection" << "\n";
return;
}
//get PF jet tags
Handle<reco::JetTagCollection> pfJetTagCollection;
if(doJetTags && !iEvent.getByToken(pfJetTagCollection_, pfJetTagCollection)){
edm::LogError ("HLTScoutingPFProducer") << "invalid collection: pfJetTagCollection" << "\n";
return;
}
//get PF candidates
Handle<reco::PFCandidateCollection> pfCandidateCollection;
if(doCandidates && !iEvent.getByToken(pfCandidateCollection_, pfCandidateCollection)){
edm::LogError ("HLTScoutingPFProducer") << "invalid collection: pfCandidateCollection" << "\n";
return;
}
//get vertices
Handle<reco::VertexCollection> vertexCollection;
if(!iEvent.getByToken(vertexCollection_, vertexCollection)){
edm::LogError ("HLTScoutingPFProducer") << "invalid collection: vertexCollection" << "\n";
return;
}
//get rho
Handle<double> rho;
if(!iEvent.getByToken(rho_, rho)){
edm::LogError ("HLTScoutingPFProducer") << "invalid collection: rho" << "\n";
return;
}
std::auto_ptr<double> outRho(new double(*rho));
//get MET
Handle<reco::METCollection> metCollection;
if(doMet && !iEvent.getByToken(metCollection_, metCollection)){
edm::LogError ("HLTScoutingPFProducer") << "invalid collection: metCollection" << "\n";
return;
}
//produce vertices
std::auto_ptr<ScoutingVertexCollection> outVertices(new ScoutingVertexCollection());
for(auto &vtx : *vertexCollection){
outVertices->emplace_back(
vtx.x(), vtx.y(), vtx.z(), vtx.zError()
);
}
//produce PF candidates
std::auto_ptr<ScoutingParticleCollection> outPFCandidates(new ScoutingParticleCollection());
if(doCandidates){
for(auto &cand : *pfCandidateCollection){
if(cand.pt() > pfCandidatePtCut){
int vertex_index = -1;
int index_counter = 0;
double dr2 = 0.0001;
for (auto &vtx: *outVertices) {
double tmp_dr2 = pow(vtx.x() - cand.vx(), 2) + pow(vtx.y() - cand.vy(), 2)
+ pow(vtx.z() - cand.vz(), 2);
if (tmp_dr2 < dr2) {
dr2 = tmp_dr2;
vertex_index = index_counter;
}
if (dr2 == 0.0)
break;
++index_counter;
}
outPFCandidates->emplace_back(
cand.pt(), cand.eta(), cand.phi(), cand.mass(), cand.pdgId(), vertex_index
);
}
}
}
//produce PF jets
std::auto_ptr<ScoutingPFJetCollection> outPFJets(new ScoutingPFJetCollection());
for(auto &jet : *pfJetCollection){
if(jet.pt() < pfJetPtCut || fabs(jet.eta()) > pfJetEtaCut) continue;
//find the jet tag corresponding to the jet
float tagValue = -20;
float minDR2 = 0.01;
if(doJetTags){
for(auto &tag : *pfJetTagCollection){
float dR2 = reco::deltaR2(jet, *(tag.first));
if(dR2 < minDR2){
minDR2 = dR2;
tagValue = tag.second;
}
}
}
//get the PF constituents of the jet
std::vector<int> candIndices;
if(doCandidates){
for(auto &cand : jet.getPFConstituents()){
if(cand->pt() > pfCandidatePtCut){
//search for the candidate in the collection
float minDR2 = 0.0001;
int matchIndex = -1;
int outIndex = 0;
for(auto &outCand : *outPFCandidates){
float dR2 = pow(cand->eta() - outCand.eta(), 2) + pow(cand->phi() - outCand.phi(), 2);
if(dR2 < minDR2){
minDR2 = dR2;
matchIndex = outIndex;
}
if(minDR2 == 0){
break;
}
outIndex++;
}
candIndices.push_back(matchIndex);
}
}
}
outPFJets->emplace_back(
jet.pt(), jet.eta(), jet.phi(), jet.mass(), jet.jetArea(),
jet.chargedHadronEnergy(), jet.neutralHadronEnergy(), jet.photonEnergy(),
jet.electronEnergy(), jet.muonEnergy(), jet.HFHadronEnergy(), jet.HFEMEnergy(),
jet.chargedHadronMultiplicity(), jet.neutralHadronMultiplicity(), jet.photonMultiplicity(),
jet.electronMultiplicity(), jet.muonMultiplicity(),
jet.HFHadronMultiplicity(), jet.HFEMMultiplicity(),
jet.hoEnergy(), tagValue, 0.0, candIndices
);
}
double metPt = -999;
double metPhi = -999;
if(doMet){
metPt = metCollection->front().pt();
metPhi = metCollection->front().phi();
}
std::auto_ptr<double> outMetPt(new double(metPt));
std::auto_ptr<double> outMetPhi(new double(metPhi));
//put output
iEvent.put(outVertices);
iEvent.put(outPFCandidates);
iEvent.put(outPFJets);
iEvent.put(outRho, "rho");
iEvent.put(outMetPt, "pfMetPt");
iEvent.put(outMetPhi, "pfMetPhi");
}
// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void HLTScoutingPFProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<edm::InputTag>("pfJetCollection",edm::InputTag("hltAK4PFJets"));
desc.add<edm::InputTag>("pfJetTagCollection",edm::InputTag("hltCombinedSecondaryVertexBJetTagsPF"));
desc.add<edm::InputTag>("pfCandidateCollection", edm::InputTag("hltParticleFlow"));
desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
desc.add<edm::InputTag>("metCollection", edm::InputTag("hltPFMETProducer"));
desc.add<edm::InputTag>("rho", edm::InputTag("hltFixedGridRhoFastjetAll"));
desc.add<double>("pfJetPtCut", 20.0);
desc.add<double>("pfJetEtaCut", 3.0);
desc.add<double>("pfCandidatePtCut", 0.6);
desc.add<bool>("doJetTags", true);
desc.add<bool>("doCandidates", true);
desc.add<bool>("doMet", true);
descriptions.add("hltScoutingPFProducer", desc);
}
//define this as a plug-in
DEFINE_FWK_MODULE(HLTScoutingPFProducer);