-
Notifications
You must be signed in to change notification settings - Fork 4.4k
/
Copy pathRun3ScoutingParticleToRecoPFCandidateProducer.cc
286 lines (243 loc) · 12 KB
/
Run3ScoutingParticleToRecoPFCandidateProducer.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
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
// system include files
#include <memory>
// user include files
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/stream/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
#include "DataFormats/Scouting/interface/Run3ScoutingParticle.h"
#include "DataFormats/Common/interface/OrphanHandle.h"
#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
#include "fastjet/contrib/SoftKiller.hh"
class Run3ScoutingParticleToRecoPFCandidateProducer : public edm::stream::EDProducer<> {
public:
explicit Run3ScoutingParticleToRecoPFCandidateProducer(const edm::ParameterSet &);
~Run3ScoutingParticleToRecoPFCandidateProducer() override;
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
void beginStream(edm::StreamID) override {}
void produce(edm::Event &iEvent, edm::EventSetup const &setup) override;
void endStream() override {}
void createPFCandidates(edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
std::unique_ptr<reco::PFCandidateCollection> &pfcands);
void createPFCandidatesSK(edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
std::unique_ptr<reco::PFCandidateCollection> &pfcands);
reco::PFCandidate createPFCand(Run3ScoutingParticle scoutingparticle);
void clearVars();
private:
const edm::EDGetTokenT<std::vector<Run3ScoutingParticle>> input_scoutingparticle_token_;
const edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> particletable_token_;
bool use_softKiller_;
bool use_CHS_;
const HepPDT::ParticleDataTable *pdTable_;
std::vector<int8_t> vertexIndex_;
std::vector<float> normchi2_;
std::vector<float> dz_;
std::vector<float> dxy_;
std::vector<float> dzsig_;
std::vector<float> dxysig_;
std::vector<int> lostInnerHits_;
std::vector<int> quality_;
std::vector<float> trkPt_;
std::vector<float> trkEta_;
std::vector<float> trkPhi_;
};
//
// constructors and destructor
//
Run3ScoutingParticleToRecoPFCandidateProducer::Run3ScoutingParticleToRecoPFCandidateProducer(
const edm::ParameterSet &iConfig)
: input_scoutingparticle_token_(consumes(iConfig.getParameter<edm::InputTag>("scoutingparticle"))),
particletable_token_(esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>()),
use_softKiller_(iConfig.getParameter<bool>("softKiller")),
use_CHS_(iConfig.getParameter<bool>("CHS")) {
//register products
produces<reco::PFCandidateCollection>();
produces<edm::ValueMap<int>>("vertexIndex");
produces<edm::ValueMap<float>>("normchi2");
produces<edm::ValueMap<float>>("dz");
produces<edm::ValueMap<float>>("dxy");
produces<edm::ValueMap<float>>("dzsig");
produces<edm::ValueMap<float>>("dxysig");
produces<edm::ValueMap<int>>("lostInnerHits");
produces<edm::ValueMap<int>>("quality");
produces<edm::ValueMap<float>>("trkPt");
produces<edm::ValueMap<float>>("trkEta");
produces<edm::ValueMap<float>>("trkPhi");
}
Run3ScoutingParticleToRecoPFCandidateProducer::~Run3ScoutingParticleToRecoPFCandidateProducer() = default;
reco::PFCandidate Run3ScoutingParticleToRecoPFCandidateProducer::createPFCand(Run3ScoutingParticle scoutingparticle) {
auto m = pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId())) != nullptr
? pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId()))->mass()
: -99.f;
auto q = pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId())) != nullptr
? pdTable_->particle(HepPDT::ParticleID(scoutingparticle.pdgId()))->charge()
: -99.f;
if (m < -90 or q < -90) {
LogDebug("createPFCand") << "<Run3ScoutingParticleToRecoPFCandidateProducer::createPFCand>:" << std::endl
<< "Unrecognisable pdgId - skipping particle" << std::endl;
return reco::PFCandidate();
}
float px = scoutingparticle.pt() * cos(scoutingparticle.phi());
float py = scoutingparticle.pt() * sin(scoutingparticle.phi());
float pz = scoutingparticle.pt() * sinh(scoutingparticle.eta());
float p = scoutingparticle.pt() * cosh(scoutingparticle.eta());
float energy = std::sqrt(p * p + m * m);
reco::Particle::LorentzVector p4(px, py, pz, energy);
static const reco::PFCandidate dummy;
auto pfcand = reco::PFCandidate(q, p4, dummy.translatePdgIdToType(scoutingparticle.pdgId()));
bool relativeTrackVars = scoutingparticle.relative_trk_vars();
vertexIndex_.push_back(scoutingparticle.vertex());
normchi2_.push_back(scoutingparticle.normchi2());
dz_.push_back(scoutingparticle.dz());
dxy_.push_back(scoutingparticle.dxy());
dzsig_.push_back(scoutingparticle.dzsig());
dxysig_.push_back(scoutingparticle.dxysig());
lostInnerHits_.push_back(scoutingparticle.lostInnerHits());
quality_.push_back(scoutingparticle.quality());
trkPt_.push_back(relativeTrackVars ? scoutingparticle.trk_pt() + scoutingparticle.pt() : scoutingparticle.trk_pt());
trkEta_.push_back(relativeTrackVars ? scoutingparticle.trk_eta() + scoutingparticle.eta()
: scoutingparticle.trk_eta());
trkPhi_.push_back(relativeTrackVars ? scoutingparticle.trk_phi() + scoutingparticle.phi()
: scoutingparticle.trk_phi());
return pfcand;
}
void Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidates(
edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
std::unique_ptr<reco::PFCandidateCollection> &pfcands) {
for (unsigned int icand = 0; icand < scoutingparticleHandle->size(); ++icand) {
auto &scoutingparticle = (*scoutingparticleHandle)[icand];
if (use_CHS_ and scoutingparticle.vertex() > 0)
continue;
auto pfcand = createPFCand(scoutingparticle);
if (pfcand.energy() != 0)
pfcands->push_back(pfcand);
}
}
void Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidatesSK(
edm::Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle,
std::unique_ptr<reco::PFCandidateCollection> &pfcands) {
std::vector<fastjet::PseudoJet> fj;
for (auto iter = scoutingparticleHandle->begin(),
ibegin = scoutingparticleHandle->begin(),
iend = scoutingparticleHandle->end();
iter != iend;
++iter) {
auto m = pdTable_->particle(HepPDT::ParticleID(iter->pdgId())) != nullptr
? pdTable_->particle(HepPDT::ParticleID(iter->pdgId()))->mass()
: -99.f;
if (m < -90) {
LogDebug("createPFCandidatesSK") << "<Run3ScoutingParticleToRecoPFCandidateProducer::createPFCandidatesSK>:"
<< std::endl
<< "Unrecognisable pdgId - skipping particle" << std::endl;
continue;
}
math::PtEtaPhiMLorentzVector p4(iter->pt(), iter->eta(), iter->phi(), m);
fj.push_back(fastjet::PseudoJet(p4.px(), p4.py(), p4.pz(), p4.energy()));
fj.back().set_user_index(iter - ibegin);
}
fastjet::contrib::SoftKiller soft_killer(5, 0.4);
std::vector<fastjet::PseudoJet> soft_killed_particles = soft_killer(fj);
for (auto &particle : soft_killed_particles) {
const Run3ScoutingParticle scoutingparticle = scoutingparticleHandle->at(particle.user_index());
auto pfcand = createPFCand(scoutingparticle);
if (pfcand.energy() != 0)
pfcands->push_back(pfcand);
}
}
// ------------ method called to produce the data ------------
void Run3ScoutingParticleToRecoPFCandidateProducer::produce(edm::Event &iEvent, edm::EventSetup const &setup) {
using namespace edm;
auto pdt = setup.getHandle(particletable_token_);
pdTable_ = pdt.product();
Handle<std::vector<Run3ScoutingParticle>> scoutingparticleHandle;
iEvent.getByToken(input_scoutingparticle_token_, scoutingparticleHandle);
auto pfcands = std::make_unique<reco::PFCandidateCollection>();
if (use_softKiller_) {
createPFCandidatesSK(scoutingparticleHandle, pfcands);
} else {
createPFCandidates(scoutingparticleHandle, pfcands);
}
edm::OrphanHandle<reco::PFCandidateCollection> oh = iEvent.put(std::move(pfcands));
std::unique_ptr<edm::ValueMap<int>> vertexIndex_VM(new edm::ValueMap<int>());
edm::ValueMap<int>::Filler filler_vertexIndex(*vertexIndex_VM);
filler_vertexIndex.insert(oh, vertexIndex_.begin(), vertexIndex_.end());
filler_vertexIndex.fill();
iEvent.put(std::move(vertexIndex_VM), "vertexIndex");
std::unique_ptr<edm::ValueMap<float>> normchi2_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_normchi2(*normchi2_VM);
filler_normchi2.insert(oh, normchi2_.begin(), normchi2_.end());
filler_normchi2.fill();
iEvent.put(std::move(normchi2_VM), "normchi2");
std::unique_ptr<edm::ValueMap<float>> dz_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_dz(*dz_VM);
filler_dz.insert(oh, dz_.begin(), dz_.end());
filler_dz.fill();
iEvent.put(std::move(dz_VM), "dz");
std::unique_ptr<edm::ValueMap<float>> dxy_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_dxy(*dxy_VM);
filler_dxy.insert(oh, dxy_.begin(), dxy_.end());
filler_dxy.fill();
iEvent.put(std::move(dxy_VM), "dxy");
std::unique_ptr<edm::ValueMap<float>> dzsig_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_dzsig(*dzsig_VM);
filler_dzsig.insert(oh, dzsig_.begin(), dzsig_.end());
filler_dzsig.fill();
iEvent.put(std::move(dzsig_VM), "dzsig");
std::unique_ptr<edm::ValueMap<float>> dxysig_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_dxysig(*dxysig_VM);
filler_dxysig.insert(oh, dxysig_.begin(), dxysig_.end());
filler_dxysig.fill();
iEvent.put(std::move(dxysig_VM), "dxysig");
std::unique_ptr<edm::ValueMap<int>> lostInnerHits_VM(new edm::ValueMap<int>());
edm::ValueMap<int>::Filler filler_lostInnerHits(*lostInnerHits_VM);
filler_lostInnerHits.insert(oh, lostInnerHits_.begin(), lostInnerHits_.end());
filler_lostInnerHits.fill();
iEvent.put(std::move(lostInnerHits_VM), "lostInnerHits");
std::unique_ptr<edm::ValueMap<int>> quality_VM(new edm::ValueMap<int>());
edm::ValueMap<int>::Filler filler_quality(*quality_VM);
filler_quality.insert(oh, quality_.begin(), quality_.end());
filler_quality.fill();
iEvent.put(std::move(quality_VM), "quality");
std::unique_ptr<edm::ValueMap<float>> trkPt_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_trkPt(*trkPt_VM);
filler_trkPt.insert(oh, trkPt_.begin(), trkPt_.end());
filler_trkPt.fill();
iEvent.put(std::move(trkPt_VM), "trkPt");
std::unique_ptr<edm::ValueMap<float>> trkEta_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_trkEta(*trkEta_VM);
filler_trkEta.insert(oh, trkEta_.begin(), trkEta_.end());
filler_trkEta.fill();
iEvent.put(std::move(trkEta_VM), "trkEta");
std::unique_ptr<edm::ValueMap<float>> trkPhi_VM(new edm::ValueMap<float>());
edm::ValueMap<float>::Filler filler_trkPhi(*trkPhi_VM);
filler_trkPhi.insert(oh, trkPhi_.begin(), trkPhi_.end());
filler_trkPhi.fill();
iEvent.put(std::move(trkPhi_VM), "trkPhi");
clearVars();
}
void Run3ScoutingParticleToRecoPFCandidateProducer::clearVars() {
vertexIndex_.clear();
normchi2_.clear();
dz_.clear();
dxy_.clear();
dzsig_.clear();
dxysig_.clear();
lostInnerHits_.clear();
quality_.clear();
trkPt_.clear();
trkEta_.clear();
trkPhi_.clear();
}
// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void Run3ScoutingParticleToRecoPFCandidateProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
edm::ParameterSetDescription desc;
desc.add<edm::InputTag>("scoutingparticle", edm::InputTag("hltScoutingPFPacker"));
desc.add<bool>("softKiller", false);
desc.add<bool>("CHS", false);
descriptions.addWithDefaultLabel(desc);
}
// declare this class as a framework plugin
DEFINE_FWK_MODULE(Run3ScoutingParticleToRecoPFCandidateProducer);