Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/Framework/interface/MakerMacros.h"
0003 #include "FWCore/Framework/interface/global/EDProducer.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/Utilities/interface/InputTag.h"
0006 #include "DataFormats/Common/interface/PtrVector.h"
0007 #include "DataFormats/Common/interface/RefToBaseVector.h"
0008 #include "DataFormats/Candidate/interface/Candidate.h"
0009 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0010 #include "DataFormats/Math/interface/deltaR.h"
0011 
0012 class CandPtrProjector : public edm::global::EDProducer<> {
0013 public:
0014   explicit CandPtrProjector(edm::ParameterSet const& iConfig);
0015   void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
0016   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0017 
0018 private:
0019   edm::EDGetTokenT<edm::View<reco::Candidate>> candSrcToken_;
0020   edm::EDGetTokenT<edm::View<reco::Candidate>> vetoSrcToken_;
0021   bool useDeltaRforFootprint_;
0022   bool extendVetoBySingleSourcePtr_;
0023 };
0024 
0025 CandPtrProjector::CandPtrProjector(edm::ParameterSet const& iConfig)
0026     : candSrcToken_{consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("src"))},
0027       vetoSrcToken_{consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("veto"))},
0028       useDeltaRforFootprint_(iConfig.getParameter<bool>("useDeltaRforFootprint")),
0029       extendVetoBySingleSourcePtr_(iConfig.getParameter<bool>("extendVetoBySingleSourcePtr")) {
0030   produces<edm::PtrVector<reco::Candidate>>();
0031 }
0032 
0033 void CandPtrProjector::produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const {
0034   using namespace edm;
0035   Handle<View<reco::Candidate>> vetoes;
0036   iEvent.getByToken(vetoSrcToken_, vetoes);
0037 
0038   auto result = std::make_unique<PtrVector<reco::Candidate>>();
0039   std::set<reco::CandidatePtr> vetoedPtrs;
0040   for (auto const& veto : *vetoes) {
0041     auto const n = veto.numberOfSourceCandidatePtrs();
0042     for (size_t j{}; j < n; ++j) {
0043       vetoedPtrs.insert(veto.sourceCandidatePtr(j));
0044     }
0045   }
0046 
0047   Handle<View<reco::Candidate>> cands;
0048   iEvent.getByToken(candSrcToken_, cands);
0049   for (size_t i{}; i < cands->size(); ++i) {
0050     auto const c = cands->ptrAt(i);
0051     if (vetoedPtrs.find(c) == vetoedPtrs.cend()) {
0052       bool addcand = true;
0053       if ((extendVetoBySingleSourcePtr_) && (c->numberOfSourceCandidatePtrs() == 1))
0054         addcand = (vetoedPtrs.find(c->sourceCandidatePtr(0)) == vetoedPtrs.cend());
0055       if (useDeltaRforFootprint_)
0056         for (const auto& it : vetoedPtrs)
0057           if (it.isNonnull() && it.isAvailable() && reco::deltaR2(*it, *c) < 0.00000025) {
0058             addcand = false;
0059             break;
0060           }
0061       if (addcand)
0062         result->push_back(c);
0063     }
0064   }
0065   iEvent.put(std::move(result));
0066 }
0067 
0068 void CandPtrProjector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0069   edm::ParameterSetDescription desc;
0070   desc.add<edm::InputTag>("src");
0071   desc.add<edm::InputTag>("veto");
0072   desc.add<bool>("useDeltaRforFootprint", false);
0073   desc.add<bool>("extendVetoBySingleSourcePtr", true);
0074   descriptions.addWithDefaultLabel(desc);
0075 }
0076 
0077 DEFINE_FWK_MODULE(CandPtrProjector);