Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:49

0001 #include "DataFormats/Common/interface/ValueMap.h"
0002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0003 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0004 #include "FWCore/Framework/interface/global/EDProducer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 
0010 class GEDGsfElectronValueMapProducer : public edm::global::EDProducer<> {
0011 public:
0012   explicit GEDGsfElectronValueMapProducer(const edm::ParameterSet&);
0013   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0014   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0015 
0016 private:
0017   const edm::EDGetTokenT<reco::GsfElectronCollection> electronsToken_;
0018   const edm::EDGetTokenT<reco::PFCandidateCollection> pfCandsToken_;
0019   const edm::EDPutTokenT<edm::ValueMap<reco::GsfElectronRef>> putToken_;
0020 };
0021 
0022 GEDGsfElectronValueMapProducer::GEDGsfElectronValueMapProducer(const edm::ParameterSet& cfg)
0023     : electronsToken_(consumes<reco::GsfElectronCollection>(cfg.getParameter<edm::InputTag>("gedGsfElectrons"))),
0024       pfCandsToken_(consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("egmPFCandidatesTag"))),
0025       putToken_{produces<edm::ValueMap<reco::GsfElectronRef>>()} {}
0026 
0027 void GEDGsfElectronValueMapProducer::produce(edm::StreamID, edm::Event& event, const edm::EventSetup&) const {
0028   auto electrons = event.getHandle(electronsToken_);
0029   auto pfCandidatesHandle = event.getHandle(pfCandsToken_);
0030 
0031   // ValueMap
0032   edm::ValueMap<reco::GsfElectronRef> valMap{};
0033   edm::ValueMap<reco::GsfElectronRef>::Filler valMapFiller(valMap);
0034 
0035   //Loop over the collection of PFFCandidates
0036   std::vector<reco::GsfElectronRef> values;
0037 
0038   for (auto const& pfCandidate : *pfCandidatesHandle) {
0039     reco::GsfElectronRef myRef;
0040     // First check that the GsfTrack is non null
0041     if (pfCandidate.gsfTrackRef().isNonnull()) {
0042       // now look for the corresponding GsfElectron
0043       const auto itcheck = std::find_if(electrons->begin(), electrons->end(), [&pfCandidate](const auto& ele) {
0044         return (ele.gsfTrack() == pfCandidate.gsfTrackRef());
0045       });
0046       if (itcheck != electrons->end()) {
0047         // Build the Ref from the handle and the index
0048         myRef = reco::GsfElectronRef(electrons, itcheck - electrons->begin());
0049       }
0050     }
0051     values.push_back(myRef);
0052   }
0053   valMapFiller.insert(pfCandidatesHandle, values.begin(), values.end());
0054 
0055   valMapFiller.fill();
0056   event.emplace(putToken_, valMap);
0057 }
0058 
0059 void GEDGsfElectronValueMapProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0060   edm::ParameterSetDescription desc;
0061   desc.add<edm::InputTag>("gedGsfElectrons", {"gedGsfElectronsTmp"});
0062   desc.add<edm::InputTag>("egmPFCandidatesTag", {"particleFlowEGamma"});
0063   descriptions.add("gedGsfElectronValueMapsTmp", desc);
0064 }
0065 
0066 #include "FWCore/Framework/interface/MakerMacros.h"
0067 DEFINE_FWK_MODULE(GEDGsfElectronValueMapProducer);