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
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
|
#include "CommonTools/PileupAlgos/interface/PuppiAlgo.h"
#include "CommonTools/PileupAlgos/interface/PuppiContainer.h"
#include "DataFormats/Candidate/interface/Candidate.h"
#include "DataFormats/Candidate/interface/CandidateFwd.h"
#include "DataFormats/Common/interface/Association.h"
#include "DataFormats/Common/interface/ValueMap.h"
#include "DataFormats/Common/interface/View.h"
#include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
#include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
#include "DataFormats/Math/interface/LorentzVector.h"
#include "DataFormats/Math/interface/PtEtaPhiMass.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
#include "DataFormats/PatCandidates/interface/PackedCandidate.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/stream/EDProducer.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include <memory>
// ------------------------------------------------------------------------------------------
class PuppiProducer : public edm::stream::EDProducer<> {
public:
explicit PuppiProducer(const edm::ParameterSet&);
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
typedef math::XYZTLorentzVector LorentzVector;
typedef std::vector<LorentzVector> LorentzVectorCollection;
typedef reco::VertexCollection VertexCollection;
typedef edm::View<reco::Candidate> CandidateView;
typedef std::vector<reco::PFCandidate> PFInputCollection;
typedef std::vector<reco::PFCandidate> PFOutputCollection;
typedef std::vector<pat::PackedCandidate> PackedOutputCollection;
typedef edm::View<reco::PFCandidate> PFView;
typedef edm::Association<reco::VertexCollection> CandToVertex;
private:
void produce(edm::Event&, const edm::EventSetup&) override;
edm::EDGetTokenT<CandidateView> tokenPFCandidates_;
edm::EDGetTokenT<VertexCollection> tokenVertices_;
edm::EDGetTokenT<CandToVertex> tokenVertexAssociation_;
edm::EDGetTokenT<edm::ValueMap<int>> tokenVertexAssociationQuality_;
edm::EDGetTokenT<PuppiContainer> tokenPuppiContainer_;
edm::EDGetTokenT<PFOutputCollection> tokenPuppiCandidates_;
edm::EDGetTokenT<PackedOutputCollection> tokenPackedPuppiCandidates_;
edm::EDGetTokenT<double> puProxyValueToken_;
edm::EDPutTokenT<edm::ValueMap<float>> ptokenPupOut_;
edm::EDPutTokenT<edm::ValueMap<LorentzVector>> ptokenP4PupOut_;
edm::EDPutTokenT<edm::ValueMap<reco::CandidatePtr>> ptokenValues_;
edm::EDPutTokenT<pat::PackedCandidateCollection> ptokenPackedPuppiCandidates_;
edm::EDPutTokenT<reco::PFCandidateCollection> ptokenPuppiCandidates_;
edm::EDPutTokenT<double> ptokenNalgos_;
edm::EDPutTokenT<std::vector<double>> ptokenRawAlphas_;
edm::EDPutTokenT<std::vector<double>> ptokenAlphas_;
edm::EDPutTokenT<std::vector<double>> ptokenAlphasMed_;
edm::EDPutTokenT<std::vector<double>> ptokenAlphasRms_;
std::string fPuppiName;
std::string fPFName;
std::string fPVName;
bool fUseVertexAssociation;
int vertexAssociationQuality_;
bool fPuppiDiagnostics;
bool fPuppiNoLep;
bool fUseFromPVLooseTight;
bool fUseFromPV2Recovery;
bool fUseDZ;
bool fUseDZforPileup;
double fDZCut;
double fEtaMinUseDZ;
double fPtMaxCharged;
double fEtaMaxCharged;
double fPtMaxPhotons;
double fEtaMaxPhotons;
double fPtMinForFromPV2Recovery;
uint fNumOfPUVtxsForCharged;
double fDZCutForChargedFromPUVtxs;
bool fUseExistingWeights;
bool fApplyPhotonProtectionForExistingWeights;
bool fClonePackedCands;
int fVtxNdofCut;
double fVtxZCut;
bool fUsePUProxyValue;
PuppiContainer fPuppiContainer;
};
// ------------------------------------------------------------------------------------------
PuppiProducer::PuppiProducer(const edm::ParameterSet& iConfig) : fPuppiContainer(iConfig) {
fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
fPuppiNoLep = iConfig.getParameter<bool>("puppiNoLep");
fUseFromPVLooseTight = iConfig.getParameter<bool>("UseFromPVLooseTight");
fUseFromPV2Recovery = iConfig.getParameter<bool>("UseFromPV2Recovery");
fUseDZ = iConfig.getParameter<bool>("UseDeltaZCut");
fUseDZforPileup = iConfig.getParameter<bool>("UseDeltaZCutForPileup");
fDZCut = iConfig.getParameter<double>("DeltaZCut");
fEtaMinUseDZ = iConfig.getParameter<double>("EtaMinUseDeltaZ");
fPtMaxCharged = iConfig.getParameter<double>("PtMaxCharged");
fEtaMaxCharged = iConfig.getParameter<double>("EtaMaxCharged");
fPtMaxPhotons = iConfig.getParameter<double>("PtMaxPhotons");
fEtaMaxPhotons = iConfig.getParameter<double>("EtaMaxPhotons");
fPtMinForFromPV2Recovery = iConfig.getParameter<double>("PtMinForFromPV2Recovery");
fNumOfPUVtxsForCharged = iConfig.getParameter<uint>("NumOfPUVtxsForCharged");
fDZCutForChargedFromPUVtxs = iConfig.getParameter<double>("DeltaZCutForChargedFromPUVtxs");
fUseExistingWeights = iConfig.getParameter<bool>("useExistingWeights");
fApplyPhotonProtectionForExistingWeights = iConfig.getParameter<bool>("applyPhotonProtectionForExistingWeights");
fClonePackedCands = iConfig.getParameter<bool>("clonePackedCands");
fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
fVtxZCut = iConfig.getParameter<double>("vtxZCut");
tokenPFCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
tokenVertices_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
fUseVertexAssociation = iConfig.getParameter<bool>("useVertexAssociation");
vertexAssociationQuality_ = iConfig.getParameter<int>("vertexAssociationQuality");
if (fUseVertexAssociation) {
tokenVertexAssociation_ = consumes<CandToVertex>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
tokenVertexAssociationQuality_ =
consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
}
fUsePUProxyValue = iConfig.getParameter<bool>("usePUProxyValue");
if (fUsePUProxyValue) {
puProxyValueToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("PUProxyValue"));
}
ptokenPupOut_ = produces<edm::ValueMap<float>>();
ptokenP4PupOut_ = produces<edm::ValueMap<LorentzVector>>();
ptokenValues_ = produces<edm::ValueMap<reco::CandidatePtr>>();
if (fUseExistingWeights || fClonePackedCands)
ptokenPackedPuppiCandidates_ = produces<pat::PackedCandidateCollection>();
else {
ptokenPuppiCandidates_ = produces<reco::PFCandidateCollection>();
}
if (fPuppiDiagnostics) {
ptokenNalgos_ = produces<double>("PuppiNAlgos");
ptokenRawAlphas_ = produces<std::vector<double>>("PuppiRawAlphas");
ptokenAlphas_ = produces<std::vector<double>>("PuppiAlphas");
ptokenAlphasMed_ = produces<std::vector<double>>("PuppiAlphasMed");
ptokenAlphasRms_ = produces<std::vector<double>>("PuppiAlphasRms");
}
}
// ------------------------------------------------------------------------------------------
void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
// Get PFCandidate Collection
edm::Handle<CandidateView> hPFProduct;
iEvent.getByToken(tokenPFCandidates_, hPFProduct);
const CandidateView* pfCol = hPFProduct.product();
// Get vertex collection w/PV as the first entry?
edm::Handle<reco::VertexCollection> hVertexProduct;
iEvent.getByToken(tokenVertices_, hVertexProduct);
const reco::VertexCollection* pvCol = hVertexProduct.product();
edm::Association<reco::VertexCollection> associatedPV;
edm::ValueMap<int> associationQuality;
if ((fUseVertexAssociation) && (!fUseExistingWeights)) {
associatedPV = iEvent.get(tokenVertexAssociation_);
associationQuality = iEvent.get(tokenVertexAssociationQuality_);
}
double puProxyValue = 0.;
if (fUsePUProxyValue) {
puProxyValue = iEvent.get(puProxyValueToken_);
} else {
for (auto const& vtx : *pvCol) {
if (!vtx.isFake() && vtx.ndof() >= fVtxNdofCut && std::abs(vtx.z()) <= fVtxZCut)
++puProxyValue;
}
}
std::vector<RecoObj> recoObjCollection;
PuppiContainer::Weights weightsInfo;
std::vector<double> lWeights;
if (!fUseExistingWeights) {
//Fill the reco objects
recoObjCollection.reserve(pfCol->size());
int iCand = 0;
for (auto const& aPF : *pfCol) {
RecoObj pReco;
pReco.pt = aPF.pt();
pReco.eta = aPF.eta();
pReco.phi = aPF.phi();
pReco.m = aPF.mass();
pReco.rapidity = aPF.rapidity();
pReco.charge = aPF.charge();
pReco.pdgId = aPF.pdgId();
const reco::Vertex* closestVtx = nullptr;
double pDZ = -9999;
double pD0 = -9999;
uint pVtxId = 0;
bool isLepton = ((std::abs(pReco.pdgId) == 11) || (std::abs(pReco.pdgId) == 13));
const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
if (fUseVertexAssociation) {
const reco::VertexRef& PVOrig = associatedPV[reco::CandidatePtr(hPFProduct, iCand)];
int quality = associationQuality[reco::CandidatePtr(hPFProduct, iCand)];
if (PVOrig.isNonnull() && (quality >= vertexAssociationQuality_)) {
closestVtx = PVOrig.get();
pVtxId = PVOrig.key();
}
if (std::abs(pReco.charge) == 0)
pReco.id = 0;
else if (fPuppiNoLep && isLepton)
pReco.id = 3;
else if (closestVtx != nullptr && pVtxId == 0)
pReco.id = 1; // Associated to main vertex
else if (closestVtx != nullptr && pVtxId > 0)
pReco.id = 2; // Associated to PU
else
pReco.id = 0; // Unassociated
} else if (lPack == nullptr) {
const reco::PFCandidate* pPF = dynamic_cast<const reco::PFCandidate*>(&aPF);
double curdz = 9999;
int closestVtxForUnassociateds = -9999;
const reco::TrackRef aTrackRef = pPF->trackRef();
bool lFirst = true;
for (auto const& aV : *pvCol) {
if (lFirst) {
if (aTrackRef.isNonnull()) {
pDZ = aTrackRef->dz(aV.position());
pD0 = aTrackRef->d0();
} else if (pPF->gsfTrackRef().isNonnull()) {
pDZ = pPF->gsfTrackRef()->dz(aV.position());
pD0 = pPF->gsfTrackRef()->d0();
}
lFirst = false;
if (pDZ > -9999)
pVtxId = 0;
}
if (aTrackRef.isNonnull() && aV.trackWeight(pPF->trackRef()) > 0) {
closestVtx = &aV;
break;
}
// in case it's unassocciated, keep more info
double tmpdz = 99999;
if (aTrackRef.isNonnull())
tmpdz = aTrackRef->dz(aV.position());
else if (pPF->gsfTrackRef().isNonnull())
tmpdz = pPF->gsfTrackRef()->dz(aV.position());
if (std::abs(tmpdz) < curdz) {
curdz = std::abs(tmpdz);
closestVtxForUnassociateds = pVtxId;
}
pVtxId++;
}
int tmpFromPV = 0;
// mocking the miniAOD definitions
if (std::abs(pReco.charge) > 0) {
if (closestVtx != nullptr && pVtxId > 0)
tmpFromPV = 0;
if (closestVtx != nullptr && pVtxId == 0)
tmpFromPV = 3;
if (closestVtx == nullptr && closestVtxForUnassociateds == 0)
tmpFromPV = 2;
if (closestVtx == nullptr && closestVtxForUnassociateds != 0)
tmpFromPV = 1;
}
pReco.dZ = pDZ;
pReco.d0 = pD0;
pReco.id = 0;
if (std::abs(pReco.charge) == 0) {
pReco.id = 0;
} else {
if (fPuppiNoLep && isLepton)
pReco.id = 3;
else if (tmpFromPV == 0) {
pReco.id = 2;
if (fNumOfPUVtxsForCharged > 0 and (pVtxId <= fNumOfPUVtxsForCharged) and
(std::abs(pDZ) < fDZCutForChargedFromPUVtxs))
pReco.id = 1;
} else if (tmpFromPV == 3)
pReco.id = 1;
else if (tmpFromPV == 1 || tmpFromPV == 2) {
pReco.id = 0;
if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
pReco.id = 1;
else if (std::abs(pReco.eta) > fEtaMaxCharged)
pReco.id = 1;
else if ((fUseDZ) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) < fDZCut))
pReco.id = 1;
else if (fUseFromPV2Recovery && tmpFromPV == 2 && (pReco.pt > fPtMinForFromPV2Recovery))
pReco.id = 1;
else if ((fUseDZforPileup) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) >= fDZCut))
pReco.id = 2;
else if (fUseFromPVLooseTight && tmpFromPV == 1)
pReco.id = 2;
else if (fUseFromPVLooseTight && tmpFromPV == 2)
pReco.id = 1;
}
}
} else if (lPack->vertexRef().isNonnull()) {
pDZ = lPack->dz();
pD0 = lPack->dxy();
pReco.dZ = pDZ;
pReco.d0 = pD0;
pReco.id = 0;
if (std::abs(pReco.charge) == 0) {
pReco.id = 0;
}
if (std::abs(pReco.charge) > 0) {
if (fPuppiNoLep && isLepton) {
pReco.id = 3;
} else if (lPack->fromPV() == 0) {
pReco.id = 2;
if ((fNumOfPUVtxsForCharged > 0) and (std::abs(pDZ) < fDZCutForChargedFromPUVtxs)) {
for (size_t puVtx_idx = 1; puVtx_idx <= fNumOfPUVtxsForCharged && puVtx_idx < pvCol->size();
++puVtx_idx) {
if (lPack->fromPV(puVtx_idx) >= 2) {
pReco.id = 1;
break;
}
}
}
} else if (lPack->fromPV() == (pat::PackedCandidate::PVUsedInFit)) {
pReco.id = 1;
} else if (lPack->fromPV() == (pat::PackedCandidate::PVTight) ||
lPack->fromPV() == (pat::PackedCandidate::PVLoose)) {
pReco.id = 0;
if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
pReco.id = 1;
else if (std::abs(pReco.eta) > fEtaMaxCharged)
pReco.id = 1;
else if ((fUseDZ) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) < fDZCut))
pReco.id = 1;
else if (fUseFromPV2Recovery && lPack->fromPV() == (pat::PackedCandidate::PVTight) &&
(pReco.pt > fPtMinForFromPV2Recovery))
pReco.id = 1;
else if ((fUseDZforPileup) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) >= fDZCut))
pReco.id = 2;
else if (fUseFromPVLooseTight && lPack->fromPV() == (pat::PackedCandidate::PVLoose))
pReco.id = 2;
else if (fUseFromPVLooseTight && lPack->fromPV() == (pat::PackedCandidate::PVTight))
pReco.id = 1;
}
}
}
recoObjCollection.push_back(pReco);
iCand++;
}
//Compute the weights and get the particles
weightsInfo = fPuppiContainer.calculatePuppiWeights(recoObjCollection, puProxyValue);
lWeights = std::move(weightsInfo.weights);
} else {
//Use the existing weights
lWeights.reserve(pfCol->size());
for (auto const& aPF : *pfCol) {
const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
float curpupweight = -1.;
if (lPack == nullptr) {
// throw error
throw edm::Exception(edm::errors::LogicError,
"PuppiProducer: cannot get weights since inputs are not PackedCandidates");
} else {
if (fPuppiNoLep) {
curpupweight = lPack->puppiWeightNoLep();
} else {
curpupweight = lPack->puppiWeight();
}
}
// Optional: Protect high pT photons (important for gamma to hadronic recoil balance) for existing weights.
if (fApplyPhotonProtectionForExistingWeights && (fPtMaxPhotons > 0) && (lPack->pdgId() == 22) &&
(std::abs(lPack->eta()) < fEtaMaxPhotons) && (lPack->pt() > fPtMaxPhotons))
curpupweight = 1;
lWeights.push_back(curpupweight);
}
}
//Fill it into the event
edm::ValueMap<float> lPupOut;
edm::ValueMap<float>::Filler lPupFiller(lPupOut);
lPupFiller.insert(hPFProduct, lWeights.begin(), lWeights.end());
lPupFiller.fill();
// This is a dummy to access the "translate" method which is a
// non-static member function even though it doesn't need to be.
// Will fix in the future.
static const reco::PFCandidate dummySinceTranslateIsNotStatic;
// Fill a new PF/Packed Candidate Collection and write out the ValueMap of the new p4s.
// Since the size of the ValueMap must be equal to the input collection, we need
// to search the "puppi" particles to find a match for each input. If none is found,
// the input is set to have a four-vector of 0,0,0,0
PFOutputCollection fPuppiCandidates;
PackedOutputCollection fPackedPuppiCandidates;
edm::ValueMap<LorentzVector> p4PupOut;
LorentzVectorCollection puppiP4s;
std::vector<reco::CandidatePtr> values(hPFProduct->size());
int iCand = -1;
puppiP4s.reserve(hPFProduct->size());
if (fUseExistingWeights || fClonePackedCands)
fPackedPuppiCandidates.reserve(hPFProduct->size());
else
fPuppiCandidates.reserve(hPFProduct->size());
for (auto const& aCand : *hPFProduct) {
++iCand;
std::unique_ptr<pat::PackedCandidate> pCand;
std::unique_ptr<reco::PFCandidate> pfCand;
if (fUseExistingWeights || fClonePackedCands) {
const pat::PackedCandidate* cand = dynamic_cast<const pat::PackedCandidate*>(&aCand);
if (!cand)
throw edm::Exception(edm::errors::LogicError, "PuppiProducer: inputs are not PackedCandidates");
pCand = std::make_unique<pat::PackedCandidate>(*cand);
} else {
auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(aCand.pdgId());
const reco::PFCandidate* cand = dynamic_cast<const reco::PFCandidate*>(&aCand);
pfCand = std::make_unique<reco::PFCandidate>(cand ? *cand : reco::PFCandidate(aCand.charge(), aCand.p4(), id));
}
// Here, we are using new weights computed and putting them in the packed candidates.
if (fClonePackedCands && (!fUseExistingWeights)) {
if (fPuppiNoLep)
pCand->setPuppiWeight(pCand->puppiWeight(), lWeights[iCand]);
else
pCand->setPuppiWeight(lWeights[iCand], pCand->puppiWeightNoLep());
}
puppiP4s.emplace_back(lWeights[iCand] * aCand.px(),
lWeights[iCand] * aCand.py(),
lWeights[iCand] * aCand.pz(),
lWeights[iCand] * aCand.energy());
// Here, we are using existing weights, or we're using packed candidates.
// That is, whether or not we recomputed the weights, we store the
// source candidate appropriately, and set the p4 of the packed candidate.
if (fUseExistingWeights || fClonePackedCands) {
pCand->setP4(puppiP4s.back());
pCand->setSourceCandidatePtr(aCand.sourceCandidatePtr(0));
fPackedPuppiCandidates.push_back(*pCand);
} else {
pfCand->setP4(puppiP4s.back());
pfCand->setSourceCandidatePtr(aCand.sourceCandidatePtr(0));
fPuppiCandidates.push_back(*pfCand);
}
}
//Compute the modified p4s
edm::ValueMap<LorentzVector>::Filler p4PupFiller(p4PupOut);
p4PupFiller.insert(hPFProduct, puppiP4s.begin(), puppiP4s.end());
p4PupFiller.fill();
iEvent.emplace(ptokenPupOut_, std::move(lPupOut));
iEvent.emplace(ptokenP4PupOut_, std::move(p4PupOut));
if (fUseExistingWeights || fClonePackedCands) {
edm::OrphanHandle<pat::PackedCandidateCollection> oh =
iEvent.emplace(ptokenPackedPuppiCandidates_, fPackedPuppiCandidates);
for (unsigned int ic = 0, nc = oh->size(); ic < nc; ++ic) {
reco::CandidatePtr pkref(oh, ic);
values[ic] = pkref;
}
} else {
edm::OrphanHandle<reco::PFCandidateCollection> oh = iEvent.emplace(ptokenPuppiCandidates_, fPuppiCandidates);
for (unsigned int ic = 0, nc = oh->size(); ic < nc; ++ic) {
reco::CandidatePtr pkref(oh, ic);
values[ic] = pkref;
}
}
edm::ValueMap<reco::CandidatePtr> pfMap_p;
edm::ValueMap<reco::CandidatePtr>::Filler filler(pfMap_p);
filler.insert(hPFProduct, values.begin(), values.end());
filler.fill();
iEvent.emplace(ptokenValues_, std::move(pfMap_p));
//////////////////////////////////////////////
if (fPuppiDiagnostics && !fUseExistingWeights) {
// all the different alphas per particle
// THE alpha per particle
double nalgos(fPuppiContainer.puppiNAlgos());
iEvent.emplace(ptokenRawAlphas_, std::move(weightsInfo.puppiRawAlphas));
iEvent.emplace(ptokenNalgos_, nalgos);
iEvent.emplace(ptokenAlphas_, std::move(weightsInfo.puppiAlphas));
iEvent.emplace(ptokenAlphasMed_, std::move(weightsInfo.puppiAlphasMed));
iEvent.emplace(ptokenAlphasRms_, std::move(weightsInfo.puppiAlphasRMS));
}
}
void PuppiProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<bool>("puppiDiagnostics", false);
desc.add<bool>("puppiNoLep", false);
desc.add<bool>("UseFromPVLooseTight", false);
desc.add<bool>("UseFromPV2Recovery", false);
desc.add<bool>("UseDeltaZCut", true);
desc.add<bool>("UseDeltaZCutForPileup", true);
desc.add<double>("DeltaZCut", 0.3);
desc.add<double>("EtaMinUseDeltaZ", 0.);
desc.add<double>("PtMaxCharged", -1.);
desc.add<double>("EtaMaxCharged", 99999.);
desc.add<double>("PtMaxPhotons", -1.);
desc.add<double>("EtaMaxPhotons", 2.5);
desc.add<double>("PtMaxNeutrals", 200.);
desc.add<double>("PtMaxNeutralsStartSlope", 0.);
desc.add<double>("PtMinForFromPV2Recovery", 0.);
desc.add<uint>("NumOfPUVtxsForCharged", 0);
desc.add<double>("DeltaZCutForChargedFromPUVtxs", 0.2);
desc.add<bool>("useExistingWeights", false);
desc.add<bool>("applyPhotonProtectionForExistingWeights", false);
desc.add<bool>("clonePackedCands", false);
desc.add<int>("vtxNdofCut", 4);
desc.add<double>("vtxZCut", 24);
desc.add<edm::InputTag>("candName", edm::InputTag("particleFlow"));
desc.add<edm::InputTag>("vertexName", edm::InputTag("offlinePrimaryVertices"));
desc.add<bool>("useVertexAssociation", false);
desc.add<int>("vertexAssociationQuality", 0);
desc.add<edm::InputTag>("vertexAssociation", edm::InputTag(""));
desc.add<bool>("applyCHS", true);
desc.add<bool>("invertPuppi", false);
desc.add<bool>("useExp", false);
desc.add<double>("MinPuppiWeight", .01);
desc.add<bool>("usePUProxyValue", false);
desc.add<edm::InputTag>("PUProxyValue", edm::InputTag(""));
PuppiAlgo::fillDescriptionsPuppiAlgo(desc);
descriptions.add("PuppiProducer", desc);
}
//define this as a plug-in
DEFINE_FWK_MODULE(PuppiProducer);
|