File indexing completed on 2024-12-21 03:54:28
0001 #include "CommonTools/PileupAlgos/interface/PuppiAlgo.h"
0002 #include "CommonTools/PileupAlgos/interface/PuppiContainer.h"
0003 #include "DataFormats/Candidate/interface/Candidate.h"
0004 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0005 #include "DataFormats/Common/interface/Association.h"
0006 #include "DataFormats/Common/interface/ValueMap.h"
0007 #include "DataFormats/Common/interface/View.h"
0008 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0009 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0010 #include "DataFormats/Math/interface/LorentzVector.h"
0011 #include "DataFormats/Math/interface/PtEtaPhiMass.h"
0012 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0013 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/Framework/interface/Frameworkfwd.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/Framework/interface/stream/EDProducer.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021
0022 #include <memory>
0023
0024
0025 class PuppiProducer : public edm::stream::EDProducer<> {
0026 public:
0027 explicit PuppiProducer(const edm::ParameterSet&);
0028
0029 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0030 typedef math::XYZTLorentzVector LorentzVector;
0031 typedef std::vector<LorentzVector> LorentzVectorCollection;
0032 typedef reco::VertexCollection VertexCollection;
0033 typedef edm::View<reco::Candidate> CandidateView;
0034 typedef std::vector<reco::PFCandidate> PFInputCollection;
0035 typedef std::vector<reco::PFCandidate> PFOutputCollection;
0036 typedef std::vector<pat::PackedCandidate> PackedOutputCollection;
0037 typedef edm::View<reco::PFCandidate> PFView;
0038 typedef edm::Association<reco::VertexCollection> CandToVertex;
0039
0040 private:
0041 void produce(edm::Event&, const edm::EventSetup&) override;
0042
0043 edm::EDGetTokenT<CandidateView> tokenPFCandidates_;
0044 edm::EDGetTokenT<VertexCollection> tokenVertices_;
0045 edm::EDGetTokenT<CandToVertex> tokenVertexAssociation_;
0046 edm::EDGetTokenT<edm::ValueMap<int>> tokenVertexAssociationQuality_;
0047 edm::EDGetTokenT<PuppiContainer> tokenPuppiContainer_;
0048 edm::EDGetTokenT<PFOutputCollection> tokenPuppiCandidates_;
0049 edm::EDGetTokenT<PackedOutputCollection> tokenPackedPuppiCandidates_;
0050 edm::EDGetTokenT<double> puProxyValueToken_;
0051 edm::EDPutTokenT<edm::ValueMap<float>> ptokenPupOut_;
0052 edm::EDPutTokenT<edm::ValueMap<LorentzVector>> ptokenP4PupOut_;
0053 edm::EDPutTokenT<edm::ValueMap<reco::CandidatePtr>> ptokenValues_;
0054 edm::EDPutTokenT<pat::PackedCandidateCollection> ptokenPackedPuppiCandidates_;
0055 edm::EDPutTokenT<reco::PFCandidateCollection> ptokenPuppiCandidates_;
0056 edm::EDPutTokenT<double> ptokenNalgos_;
0057 edm::EDPutTokenT<std::vector<double>> ptokenRawAlphas_;
0058 edm::EDPutTokenT<std::vector<double>> ptokenAlphas_;
0059 edm::EDPutTokenT<std::vector<double>> ptokenAlphasMed_;
0060 edm::EDPutTokenT<std::vector<double>> ptokenAlphasRms_;
0061 std::string fPuppiName;
0062 std::string fPFName;
0063 std::string fPVName;
0064 bool fUseVertexAssociation;
0065 int vertexAssociationQuality_;
0066 bool fPuppiDiagnostics;
0067 bool fPuppiNoLep;
0068 bool fUseFromPVLooseTight;
0069 bool fUseFromPV2Recovery;
0070 bool fUseDZ;
0071 bool fUseDZforPileup;
0072 double fDZCut;
0073 double fEtaMinUseDZ;
0074 double fPtMaxCharged;
0075 double fEtaMaxCharged;
0076 double fPtMaxPhotons;
0077 double fEtaMaxPhotons;
0078 double fPtMinForFromPV2Recovery;
0079 uint fNumOfPUVtxsForCharged;
0080 double fDZCutForChargedFromPUVtxs;
0081 bool fUseExistingWeights;
0082 bool fApplyPhotonProtectionForExistingWeights;
0083 bool fClonePackedCands;
0084 int fVtxNdofCut;
0085 double fVtxZCut;
0086 bool fUsePUProxyValue;
0087 PuppiContainer fPuppiContainer;
0088 };
0089
0090
0091 PuppiProducer::PuppiProducer(const edm::ParameterSet& iConfig) : fPuppiContainer(iConfig) {
0092 fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
0093 fPuppiNoLep = iConfig.getParameter<bool>("puppiNoLep");
0094 fUseFromPVLooseTight = iConfig.getParameter<bool>("UseFromPVLooseTight");
0095 fUseFromPV2Recovery = iConfig.getParameter<bool>("UseFromPV2Recovery");
0096 fUseDZ = iConfig.getParameter<bool>("UseDeltaZCut");
0097 fUseDZforPileup = iConfig.getParameter<bool>("UseDeltaZCutForPileup");
0098 fDZCut = iConfig.getParameter<double>("DeltaZCut");
0099 fEtaMinUseDZ = iConfig.getParameter<double>("EtaMinUseDeltaZ");
0100 fPtMaxCharged = iConfig.getParameter<double>("PtMaxCharged");
0101 fEtaMaxCharged = iConfig.getParameter<double>("EtaMaxCharged");
0102 fPtMaxPhotons = iConfig.getParameter<double>("PtMaxPhotons");
0103 fEtaMaxPhotons = iConfig.getParameter<double>("EtaMaxPhotons");
0104 fPtMinForFromPV2Recovery = iConfig.getParameter<double>("PtMinForFromPV2Recovery");
0105 fNumOfPUVtxsForCharged = iConfig.getParameter<uint>("NumOfPUVtxsForCharged");
0106 fDZCutForChargedFromPUVtxs = iConfig.getParameter<double>("DeltaZCutForChargedFromPUVtxs");
0107 fUseExistingWeights = iConfig.getParameter<bool>("useExistingWeights");
0108 fApplyPhotonProtectionForExistingWeights = iConfig.getParameter<bool>("applyPhotonProtectionForExistingWeights");
0109 fClonePackedCands = iConfig.getParameter<bool>("clonePackedCands");
0110 fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
0111 fVtxZCut = iConfig.getParameter<double>("vtxZCut");
0112
0113 tokenPFCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
0114 tokenVertices_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
0115 fUseVertexAssociation = iConfig.getParameter<bool>("useVertexAssociation");
0116 vertexAssociationQuality_ = iConfig.getParameter<int>("vertexAssociationQuality");
0117 if (fUseVertexAssociation) {
0118 tokenVertexAssociation_ = consumes<CandToVertex>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
0119 tokenVertexAssociationQuality_ =
0120 consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
0121 }
0122
0123 fUsePUProxyValue = iConfig.getParameter<bool>("usePUProxyValue");
0124
0125 if (fUsePUProxyValue) {
0126 puProxyValueToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("PUProxyValue"));
0127 }
0128
0129 ptokenPupOut_ = produces<edm::ValueMap<float>>();
0130 ptokenP4PupOut_ = produces<edm::ValueMap<LorentzVector>>();
0131 ptokenValues_ = produces<edm::ValueMap<reco::CandidatePtr>>();
0132
0133 if (fUseExistingWeights || fClonePackedCands)
0134 ptokenPackedPuppiCandidates_ = produces<pat::PackedCandidateCollection>();
0135 else {
0136 ptokenPuppiCandidates_ = produces<reco::PFCandidateCollection>();
0137 }
0138
0139 if (fPuppiDiagnostics) {
0140 ptokenNalgos_ = produces<double>("PuppiNAlgos");
0141 ptokenRawAlphas_ = produces<std::vector<double>>("PuppiRawAlphas");
0142 ptokenAlphas_ = produces<std::vector<double>>("PuppiAlphas");
0143 ptokenAlphasMed_ = produces<std::vector<double>>("PuppiAlphasMed");
0144 ptokenAlphasRms_ = produces<std::vector<double>>("PuppiAlphasRms");
0145 }
0146 }
0147
0148 void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0149
0150 edm::Handle<CandidateView> hPFProduct;
0151 iEvent.getByToken(tokenPFCandidates_, hPFProduct);
0152 const CandidateView* pfCol = hPFProduct.product();
0153
0154
0155 edm::Handle<reco::VertexCollection> hVertexProduct;
0156 iEvent.getByToken(tokenVertices_, hVertexProduct);
0157 const reco::VertexCollection* pvCol = hVertexProduct.product();
0158
0159 edm::Association<reco::VertexCollection> associatedPV;
0160 edm::ValueMap<int> associationQuality;
0161 if ((fUseVertexAssociation) && (!fUseExistingWeights)) {
0162 associatedPV = iEvent.get(tokenVertexAssociation_);
0163 associationQuality = iEvent.get(tokenVertexAssociationQuality_);
0164 }
0165
0166 double puProxyValue = 0.;
0167 if (fUsePUProxyValue) {
0168 puProxyValue = iEvent.get(puProxyValueToken_);
0169 } else {
0170 for (auto const& vtx : *pvCol) {
0171 if (!vtx.isFake() && vtx.ndof() >= fVtxNdofCut && std::abs(vtx.z()) <= fVtxZCut)
0172 ++puProxyValue;
0173 }
0174 }
0175
0176 std::vector<RecoObj> recoObjCollection;
0177
0178 PuppiContainer::Weights weightsInfo;
0179 std::vector<double> lWeights;
0180 if (!fUseExistingWeights) {
0181
0182 recoObjCollection.reserve(pfCol->size());
0183 int iCand = 0;
0184 for (auto const& aPF : *pfCol) {
0185 RecoObj pReco;
0186 pReco.pt = aPF.pt();
0187 pReco.eta = aPF.eta();
0188 pReco.phi = aPF.phi();
0189 pReco.m = aPF.mass();
0190 pReco.rapidity = aPF.rapidity();
0191 pReco.charge = aPF.charge();
0192 pReco.pdgId = aPF.pdgId();
0193 const reco::Vertex* closestVtx = nullptr;
0194 double pDZ = -9999;
0195 double pD0 = -9999;
0196 uint pVtxId = 0;
0197 bool isLepton = ((std::abs(pReco.pdgId) == 11) || (std::abs(pReco.pdgId) == 13));
0198 const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
0199
0200 if (fUseVertexAssociation) {
0201 const reco::VertexRef& PVOrig = associatedPV[reco::CandidatePtr(hPFProduct, iCand)];
0202 int quality = associationQuality[reco::CandidatePtr(hPFProduct, iCand)];
0203 if (PVOrig.isNonnull() && (quality >= vertexAssociationQuality_)) {
0204 closestVtx = PVOrig.get();
0205 pVtxId = PVOrig.key();
0206 }
0207 if (std::abs(pReco.charge) == 0)
0208 pReco.id = 0;
0209 else if (fPuppiNoLep && isLepton)
0210 pReco.id = 3;
0211 else if (closestVtx != nullptr && pVtxId == 0)
0212 pReco.id = 1;
0213 else if (closestVtx != nullptr && pVtxId > 0)
0214 pReco.id = 2;
0215 else
0216 pReco.id = 0;
0217 } else if (lPack == nullptr) {
0218 const reco::PFCandidate* pPF = dynamic_cast<const reco::PFCandidate*>(&aPF);
0219 double curdz = 9999;
0220 int closestVtxForUnassociateds = -9999;
0221 const reco::TrackRef aTrackRef = pPF->trackRef();
0222 bool lFirst = true;
0223 for (auto const& aV : *pvCol) {
0224 if (lFirst) {
0225 if (aTrackRef.isNonnull()) {
0226 pDZ = aTrackRef->dz(aV.position());
0227 pD0 = aTrackRef->d0();
0228 } else if (pPF->gsfTrackRef().isNonnull()) {
0229 pDZ = pPF->gsfTrackRef()->dz(aV.position());
0230 pD0 = pPF->gsfTrackRef()->d0();
0231 }
0232 lFirst = false;
0233 if (pDZ > -9999)
0234 pVtxId = 0;
0235 }
0236 if (aTrackRef.isNonnull() && aV.trackWeight(pPF->trackRef()) > 0) {
0237 closestVtx = &aV;
0238 break;
0239 }
0240
0241 double tmpdz = 99999;
0242 if (aTrackRef.isNonnull())
0243 tmpdz = aTrackRef->dz(aV.position());
0244 else if (pPF->gsfTrackRef().isNonnull())
0245 tmpdz = pPF->gsfTrackRef()->dz(aV.position());
0246 if (std::abs(tmpdz) < curdz) {
0247 curdz = std::abs(tmpdz);
0248 closestVtxForUnassociateds = pVtxId;
0249 }
0250 pVtxId++;
0251 }
0252 int tmpFromPV = 0;
0253
0254 if (std::abs(pReco.charge) > 0) {
0255 if (closestVtx != nullptr && pVtxId > 0)
0256 tmpFromPV = 0;
0257 if (closestVtx != nullptr && pVtxId == 0)
0258 tmpFromPV = 3;
0259 if (closestVtx == nullptr && closestVtxForUnassociateds == 0)
0260 tmpFromPV = 2;
0261 if (closestVtx == nullptr && closestVtxForUnassociateds != 0)
0262 tmpFromPV = 1;
0263 }
0264 pReco.dZ = pDZ;
0265 pReco.d0 = pD0;
0266 pReco.id = 0;
0267 if (std::abs(pReco.charge) == 0) {
0268 pReco.id = 0;
0269 } else {
0270 if (fPuppiNoLep && isLepton)
0271 pReco.id = 3;
0272 else if (tmpFromPV == 0) {
0273 pReco.id = 2;
0274 if (fNumOfPUVtxsForCharged > 0 and (pVtxId <= fNumOfPUVtxsForCharged) and
0275 (std::abs(pDZ) < fDZCutForChargedFromPUVtxs))
0276 pReco.id = 1;
0277 } else if (tmpFromPV == 3)
0278 pReco.id = 1;
0279 else if (tmpFromPV == 1 || tmpFromPV == 2) {
0280 pReco.id = 0;
0281 if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
0282 pReco.id = 1;
0283 else if (std::abs(pReco.eta) > fEtaMaxCharged)
0284 pReco.id = 1;
0285 else if ((fUseDZ) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) < fDZCut))
0286 pReco.id = 1;
0287 else if (fUseFromPV2Recovery && tmpFromPV == 2 && (pReco.pt > fPtMinForFromPV2Recovery))
0288 pReco.id = 1;
0289 else if ((fUseDZforPileup) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) >= fDZCut))
0290 pReco.id = 2;
0291 else if (fUseFromPVLooseTight && tmpFromPV == 1)
0292 pReco.id = 2;
0293 else if (fUseFromPVLooseTight && tmpFromPV == 2)
0294 pReco.id = 1;
0295 }
0296 }
0297 } else if (lPack->vertexRef().isNonnull()) {
0298 pDZ = lPack->dz();
0299 pD0 = lPack->dxy();
0300 pReco.dZ = pDZ;
0301 pReco.d0 = pD0;
0302
0303 pReco.id = 0;
0304 if (std::abs(pReco.charge) == 0) {
0305 pReco.id = 0;
0306 }
0307 if (std::abs(pReco.charge) > 0) {
0308 if (fPuppiNoLep && isLepton) {
0309 pReco.id = 3;
0310 } else if (lPack->fromPV() == 0) {
0311 pReco.id = 2;
0312 if ((fNumOfPUVtxsForCharged > 0) and (std::abs(pDZ) < fDZCutForChargedFromPUVtxs)) {
0313 for (size_t puVtx_idx = 1; puVtx_idx <= fNumOfPUVtxsForCharged && puVtx_idx < pvCol->size();
0314 ++puVtx_idx) {
0315 if (lPack->fromPV(puVtx_idx) >= 2) {
0316 pReco.id = 1;
0317 break;
0318 }
0319 }
0320 }
0321 } else if (lPack->fromPV() == (pat::PackedCandidate::PVUsedInFit)) {
0322 pReco.id = 1;
0323 } else if (lPack->fromPV() == (pat::PackedCandidate::PVTight) ||
0324 lPack->fromPV() == (pat::PackedCandidate::PVLoose)) {
0325 pReco.id = 0;
0326 if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
0327 pReco.id = 1;
0328 else if (std::abs(pReco.eta) > fEtaMaxCharged)
0329 pReco.id = 1;
0330 else if ((fUseDZ) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) < fDZCut))
0331 pReco.id = 1;
0332 else if (fUseFromPV2Recovery && lPack->fromPV() == (pat::PackedCandidate::PVTight) &&
0333 (pReco.pt > fPtMinForFromPV2Recovery))
0334 pReco.id = 1;
0335 else if ((fUseDZforPileup) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) >= fDZCut))
0336 pReco.id = 2;
0337 else if (fUseFromPVLooseTight && lPack->fromPV() == (pat::PackedCandidate::PVLoose))
0338 pReco.id = 2;
0339 else if (fUseFromPVLooseTight && lPack->fromPV() == (pat::PackedCandidate::PVTight))
0340 pReco.id = 1;
0341 }
0342 }
0343 }
0344
0345 recoObjCollection.push_back(pReco);
0346 iCand++;
0347 }
0348
0349
0350 weightsInfo = fPuppiContainer.calculatePuppiWeights(recoObjCollection, puProxyValue);
0351 lWeights = std::move(weightsInfo.weights);
0352 } else {
0353
0354 lWeights.reserve(pfCol->size());
0355 for (auto const& aPF : *pfCol) {
0356 const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
0357 float curpupweight = -1.;
0358 if (lPack == nullptr) {
0359
0360 throw edm::Exception(edm::errors::LogicError,
0361 "PuppiProducer: cannot get weights since inputs are not PackedCandidates");
0362 } else {
0363 if (fPuppiNoLep) {
0364 curpupweight = lPack->puppiWeightNoLep();
0365 } else {
0366 curpupweight = lPack->puppiWeight();
0367 }
0368 }
0369
0370 if (fApplyPhotonProtectionForExistingWeights && (fPtMaxPhotons > 0) && (lPack->pdgId() == 22) &&
0371 (std::abs(lPack->eta()) < fEtaMaxPhotons) && (lPack->pt() > fPtMaxPhotons))
0372 curpupweight = 1;
0373 lWeights.push_back(curpupweight);
0374 }
0375 }
0376
0377
0378 edm::ValueMap<float> lPupOut;
0379 edm::ValueMap<float>::Filler lPupFiller(lPupOut);
0380 lPupFiller.insert(hPFProduct, lWeights.begin(), lWeights.end());
0381 lPupFiller.fill();
0382
0383
0384
0385
0386 static const reco::PFCandidate dummySinceTranslateIsNotStatic;
0387
0388
0389
0390
0391
0392 PFOutputCollection fPuppiCandidates;
0393 PackedOutputCollection fPackedPuppiCandidates;
0394
0395 edm::ValueMap<LorentzVector> p4PupOut;
0396 LorentzVectorCollection puppiP4s;
0397 std::vector<reco::CandidatePtr> values(hPFProduct->size());
0398
0399 int iCand = -1;
0400 puppiP4s.reserve(hPFProduct->size());
0401 if (fUseExistingWeights || fClonePackedCands)
0402 fPackedPuppiCandidates.reserve(hPFProduct->size());
0403 else
0404 fPuppiCandidates.reserve(hPFProduct->size());
0405 for (auto const& aCand : *hPFProduct) {
0406 ++iCand;
0407 std::unique_ptr<pat::PackedCandidate> pCand;
0408 std::unique_ptr<reco::PFCandidate> pfCand;
0409
0410 if (fUseExistingWeights || fClonePackedCands) {
0411 const pat::PackedCandidate* cand = dynamic_cast<const pat::PackedCandidate*>(&aCand);
0412 if (!cand)
0413 throw edm::Exception(edm::errors::LogicError, "PuppiProducer: inputs are not PackedCandidates");
0414 pCand = std::make_unique<pat::PackedCandidate>(*cand);
0415 } else {
0416 auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(aCand.pdgId());
0417 const reco::PFCandidate* cand = dynamic_cast<const reco::PFCandidate*>(&aCand);
0418 pfCand = std::make_unique<reco::PFCandidate>(cand ? *cand : reco::PFCandidate(aCand.charge(), aCand.p4(), id));
0419 }
0420
0421
0422 if (fClonePackedCands && (!fUseExistingWeights)) {
0423 if (fPuppiNoLep)
0424 pCand->setPuppiWeight(pCand->puppiWeight(), lWeights[iCand]);
0425 else
0426 pCand->setPuppiWeight(lWeights[iCand], pCand->puppiWeightNoLep());
0427 }
0428
0429 puppiP4s.emplace_back(lWeights[iCand] * aCand.px(),
0430 lWeights[iCand] * aCand.py(),
0431 lWeights[iCand] * aCand.pz(),
0432 lWeights[iCand] * aCand.energy());
0433
0434
0435
0436
0437 if (fUseExistingWeights || fClonePackedCands) {
0438 pCand->setP4(puppiP4s.back());
0439 pCand->setSourceCandidatePtr(aCand.sourceCandidatePtr(0));
0440 fPackedPuppiCandidates.push_back(*pCand);
0441 } else {
0442 pfCand->setP4(puppiP4s.back());
0443 pfCand->setSourceCandidatePtr(aCand.sourceCandidatePtr(0));
0444 fPuppiCandidates.push_back(*pfCand);
0445 }
0446 }
0447
0448
0449 edm::ValueMap<LorentzVector>::Filler p4PupFiller(p4PupOut);
0450 p4PupFiller.insert(hPFProduct, puppiP4s.begin(), puppiP4s.end());
0451 p4PupFiller.fill();
0452
0453 iEvent.emplace(ptokenPupOut_, std::move(lPupOut));
0454 iEvent.emplace(ptokenP4PupOut_, std::move(p4PupOut));
0455 if (fUseExistingWeights || fClonePackedCands) {
0456 edm::OrphanHandle<pat::PackedCandidateCollection> oh =
0457 iEvent.emplace(ptokenPackedPuppiCandidates_, fPackedPuppiCandidates);
0458 for (unsigned int ic = 0, nc = oh->size(); ic < nc; ++ic) {
0459 reco::CandidatePtr pkref(oh, ic);
0460 values[ic] = pkref;
0461 }
0462 } else {
0463 edm::OrphanHandle<reco::PFCandidateCollection> oh = iEvent.emplace(ptokenPuppiCandidates_, fPuppiCandidates);
0464 for (unsigned int ic = 0, nc = oh->size(); ic < nc; ++ic) {
0465 reco::CandidatePtr pkref(oh, ic);
0466 values[ic] = pkref;
0467 }
0468 }
0469 edm::ValueMap<reco::CandidatePtr> pfMap_p;
0470 edm::ValueMap<reco::CandidatePtr>::Filler filler(pfMap_p);
0471 filler.insert(hPFProduct, values.begin(), values.end());
0472 filler.fill();
0473 iEvent.emplace(ptokenValues_, std::move(pfMap_p));
0474
0475
0476 if (fPuppiDiagnostics && !fUseExistingWeights) {
0477
0478
0479 double nalgos(fPuppiContainer.puppiNAlgos());
0480
0481 iEvent.emplace(ptokenRawAlphas_, std::move(weightsInfo.puppiRawAlphas));
0482 iEvent.emplace(ptokenNalgos_, nalgos);
0483 iEvent.emplace(ptokenAlphas_, std::move(weightsInfo.puppiAlphas));
0484 iEvent.emplace(ptokenAlphasMed_, std::move(weightsInfo.puppiAlphasMed));
0485 iEvent.emplace(ptokenAlphasRms_, std::move(weightsInfo.puppiAlphasRMS));
0486 }
0487 }
0488
0489 void PuppiProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0490 edm::ParameterSetDescription desc;
0491 desc.add<bool>("puppiDiagnostics", false);
0492 desc.add<bool>("puppiNoLep", false);
0493 desc.add<bool>("UseFromPVLooseTight", false);
0494 desc.add<bool>("UseFromPV2Recovery", false);
0495 desc.add<bool>("UseDeltaZCut", true);
0496 desc.add<bool>("UseDeltaZCutForPileup", true);
0497 desc.add<double>("DeltaZCut", 0.3);
0498 desc.add<double>("EtaMinUseDeltaZ", 0.);
0499 desc.add<double>("PtMaxCharged", -1.);
0500 desc.add<double>("EtaMaxCharged", 99999.);
0501 desc.add<double>("PtMaxPhotons", -1.);
0502 desc.add<double>("EtaMaxPhotons", 2.5);
0503 desc.add<double>("PtMaxNeutrals", 200.);
0504 desc.add<double>("PtMaxNeutralsStartSlope", 0.);
0505 desc.add<double>("PtMinForFromPV2Recovery", 0.);
0506 desc.add<uint>("NumOfPUVtxsForCharged", 0);
0507 desc.add<double>("DeltaZCutForChargedFromPUVtxs", 0.2);
0508 desc.add<bool>("useExistingWeights", false);
0509 desc.add<bool>("applyPhotonProtectionForExistingWeights", false);
0510 desc.add<bool>("clonePackedCands", false);
0511 desc.add<int>("vtxNdofCut", 4);
0512 desc.add<double>("vtxZCut", 24);
0513 desc.add<edm::InputTag>("candName", edm::InputTag("particleFlow"));
0514 desc.add<edm::InputTag>("vertexName", edm::InputTag("offlinePrimaryVertices"));
0515 desc.add<bool>("useVertexAssociation", false);
0516 desc.add<int>("vertexAssociationQuality", 0);
0517 desc.add<edm::InputTag>("vertexAssociation", edm::InputTag(""));
0518 desc.add<bool>("applyCHS", true);
0519 desc.add<bool>("invertPuppi", false);
0520 desc.add<bool>("useExp", false);
0521 desc.add<double>("MinPuppiWeight", .01);
0522 desc.add<bool>("usePUProxyValue", false);
0523 desc.add<edm::InputTag>("PUProxyValue", edm::InputTag(""));
0524
0525 PuppiAlgo::fillDescriptionsPuppiAlgo(desc);
0526
0527 descriptions.add("PuppiProducer", desc);
0528 }
0529
0530 DEFINE_FWK_MODULE(PuppiProducer);