File indexing completed on 2024-09-07 04:37:24
0001 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0002 #include "DataFormats/Math/interface/deltaR.h"
0003 #include "FWCore/Framework/interface/global/EDProducer.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007
0008 #include <memory>
0009
0010 class ElectronMatchedCandidateProducer : public edm::global::EDProducer<> {
0011 public:
0012 explicit ElectronMatchedCandidateProducer(const edm::ParameterSet &);
0013
0014 private:
0015 void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
0016
0017
0018
0019 edm::EDGetTokenT<edm::View<reco::GsfElectron>> electronCollectionToken_;
0020 edm::EDGetTokenT<edm::View<reco::Candidate>> scCollectionToken_;
0021 double delRMatchingCut_;
0022 };
0023
0024 ElectronMatchedCandidateProducer::ElectronMatchedCandidateProducer(const edm::ParameterSet ¶ms) {
0025 const edm::InputTag allelectrons("gsfElectrons");
0026 electronCollectionToken_ = consumes<edm::View<reco::GsfElectron>>(
0027 params.getUntrackedParameter<edm::InputTag>("ReferenceElectronCollection", allelectrons));
0028 scCollectionToken_ = consumes<edm::View<reco::Candidate>>(params.getParameter<edm::InputTag>("src"));
0029
0030 delRMatchingCut_ = params.getUntrackedParameter<double>("deltaR", 0.30);
0031
0032 produces<edm::PtrVector<reco::Candidate>>();
0033 produces<edm::RefToBaseVector<reco::Candidate>>();
0034 }
0035
0036
0037
0038
0039
0040
0041
0042 void ElectronMatchedCandidateProducer::produce(edm::StreamID, edm::Event &event, const edm::EventSetup &) const {
0043
0044 auto outColRef = std::make_unique<edm::RefToBaseVector<reco::Candidate>>();
0045 auto outColPtr = std::make_unique<edm::PtrVector<reco::Candidate>>();
0046
0047
0048 edm::Handle<edm::View<reco::GsfElectron>> electrons;
0049 event.getByToken(electronCollectionToken_, electrons);
0050
0051
0052 edm::Handle<edm::View<reco::Candidate>> recoCandColl;
0053 event.getByToken(scCollectionToken_, recoCandColl);
0054
0055 unsigned int counter = 0;
0056
0057
0058 for (edm::View<reco::Candidate>::const_iterator scIt = recoCandColl->begin(); scIt != recoCandColl->end();
0059 ++scIt, ++counter) {
0060
0061 for (edm::View<reco::GsfElectron>::const_iterator elec = electrons->begin(); elec != electrons->end(); ++elec) {
0062 reco::SuperClusterRef eSC = elec->superCluster();
0063
0064 double dRval = reco::deltaR((float)eSC->eta(), (float)eSC->phi(), scIt->eta(), scIt->phi());
0065
0066 if (dRval < delRMatchingCut_) {
0067
0068 outColRef->push_back(recoCandColl->refAt(counter));
0069 outColPtr->push_back(recoCandColl->ptrAt(counter));
0070 }
0071 }
0072
0073 }
0074
0075 event.put(std::move(outColRef));
0076 event.put(std::move(outColPtr));
0077 }
0078
0079
0080 #include "FWCore/Framework/interface/MakerMacros.h"
0081 DEFINE_FWK_MODULE(ElectronMatchedCandidateProducer);