Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:06

0001 /**
0002   Take as input:
0003      - the electron and photon collections
0004      - the electron and photon pf maps (edm::ValueMap<std::vector<reco::PFCandidateRef>>) pointing to old PF candidates
0005      - one ValueMap<reco::PFCandidateRef> that maps old PF candidates into new PF candidates
0006   Produce as output:
0007      - the new electron and photon pf maps
0008 */
0009 
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/Framework/interface/global/EDProducer.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "DataFormats/Common/interface/ValueMap.h"
0015 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0016 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0017 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0018 
0019 class PFEGammaToCandidateRemapper : public edm::global::EDProducer<> {
0020 public:
0021   explicit PFEGammaToCandidateRemapper(const edm::ParameterSet &iConfig);
0022   ~PFEGammaToCandidateRemapper() override {}
0023 
0024   void produce(edm::StreamID iID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override;
0025 
0026 private:
0027   edm::EDGetTokenT<std::vector<reco::Photon>> photons_;
0028   edm::EDGetTokenT<std::vector<reco::GsfElectron>> electrons_;
0029   edm::EDGetTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> photon2pf_, electron2pf_;
0030   edm::EDGetTokenT<edm::ValueMap<reco::PFCandidateRef>> pf2pf_;
0031 
0032   template <typename T>
0033   void run(edm::Event &iEvent,
0034            const edm::EDGetTokenT<std::vector<T>> &colltoken,
0035            const edm::EDGetTokenT<edm::ValueMap<std::vector<reco::PFCandidateRef>>> &oldmaptoken,
0036            const edm::ValueMap<reco::PFCandidateRef> &pf2pf,
0037            const std::string &name) const {
0038     edm::Handle<std::vector<T>> handle;
0039     iEvent.getByToken(colltoken, handle);
0040 
0041     edm::Handle<edm::ValueMap<std::vector<reco::PFCandidateRef>>> oldmap;
0042     iEvent.getByToken(oldmaptoken, oldmap);
0043 
0044     std::vector<std::vector<reco::PFCandidateRef>> refs(handle->size());
0045     for (unsigned int i = 0, n = handle->size(); i < n; ++i) {
0046       edm::Ref<std::vector<T>> egRef(handle, i);
0047       for (const reco::PFCandidateRef &pfRef : (*oldmap)[egRef]) {
0048         refs[i].push_back(pf2pf[pfRef]);
0049       }
0050     }
0051 
0052     std::unique_ptr<edm::ValueMap<std::vector<reco::PFCandidateRef>>> out(
0053         new edm::ValueMap<std::vector<reco::PFCandidateRef>>());
0054     edm::ValueMap<std::vector<reco::PFCandidateRef>>::Filler filler(*out);
0055     filler.insert(handle, refs.begin(), refs.end());
0056     filler.fill();
0057     iEvent.put(std::move(out), name);
0058   }
0059 };
0060 
0061 PFEGammaToCandidateRemapper::PFEGammaToCandidateRemapper(const edm::ParameterSet &iConfig)
0062     : photons_(consumes<std::vector<reco::Photon>>(iConfig.getParameter<edm::InputTag>("photons"))),
0063       electrons_(consumes<std::vector<reco::GsfElectron>>(iConfig.getParameter<edm::InputTag>("electrons"))),
0064       photon2pf_(
0065           consumes<edm::ValueMap<std::vector<reco::PFCandidateRef>>>(iConfig.getParameter<edm::InputTag>("photon2pf"))),
0066       electron2pf_(consumes<edm::ValueMap<std::vector<reco::PFCandidateRef>>>(
0067           iConfig.getParameter<edm::InputTag>("electron2pf"))),
0068       pf2pf_(consumes<edm::ValueMap<reco::PFCandidateRef>>(iConfig.getParameter<edm::InputTag>("pf2pf"))) {
0069   produces<edm::ValueMap<std::vector<reco::PFCandidateRef>>>("photons");
0070   produces<edm::ValueMap<std::vector<reco::PFCandidateRef>>>("electrons");
0071 }
0072 
0073 void PFEGammaToCandidateRemapper::produce(edm::StreamID iID, edm::Event &iEvent, const edm::EventSetup &iSetup) const {
0074   edm::Handle<edm::ValueMap<reco::PFCandidateRef>> pf2pf;
0075   iEvent.getByToken(pf2pf_, pf2pf);
0076   run<reco::Photon>(iEvent, photons_, photon2pf_, *pf2pf, "photons");
0077   run<reco::GsfElectron>(iEvent, electrons_, electron2pf_, *pf2pf, "electrons");
0078 }
0079 
0080 #include "FWCore/Framework/interface/MakerMacros.h"
0081 DEFINE_FWK_MODULE(PFEGammaToCandidateRemapper);