Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 //****************************************************************
0002 //
0003 // A simple class to combine two candidate collections into a single
0004 // collection of ptrs to the candidates
0005 // note: it is a std::vector as the candidates are from different
0006 // collections
0007 //
0008 // collection 1 is added fully while collection 2 is deltaR cross
0009 // cleaned against collection 2
0010 //
0011 // usecase: getting a unified list of e/gammas + jets for jet core
0012 //          regional tracking
0013 //
0014 //****************************************************************
0015 
0016 #include "FWCore/Framework/interface/Frameworkfwd.h"
0017 #include "FWCore/Framework/interface/global/EDProducer.h"
0018 
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/Framework/interface/MakerMacros.h"
0021 
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "FWCore/Utilities/interface/StreamID.h"
0024 
0025 #include "DataFormats/Candidate/interface/Candidate.h"
0026 #include "DataFormats/Common/interface/View.h"
0027 #include "DataFormats/Common/interface/Ptr.h"
0028 #include "DataFormats/Math/interface/deltaR.h"
0029 
0030 class CandMergerCleanOthersByDR : public edm::global::EDProducer<> {
0031 public:
0032   explicit CandMergerCleanOthersByDR(const edm::ParameterSet&);
0033   ~CandMergerCleanOthersByDR() override {}
0034 
0035   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0036 
0037 private:
0038   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0039 
0040 private:
0041   const edm::EDGetTokenT<edm::View<reco::Candidate>> coll1Token_;
0042   const edm::EDGetTokenT<edm::View<reco::Candidate>> coll2Token_;
0043   const float maxDR2ToClean_;
0044 };
0045 
0046 namespace {
0047   double pow2(double val) { return val * val; }
0048 }  // namespace
0049 
0050 CandMergerCleanOthersByDR::CandMergerCleanOthersByDR(const edm::ParameterSet& iConfig)
0051     : coll1Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("coll1"))),
0052       coll2Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("coll2"))),
0053       maxDR2ToClean_(pow2(iConfig.getParameter<double>("maxDRToClean"))) {
0054   produces<std::vector<edm::Ptr<reco::Candidate>>>();
0055 }
0056 
0057 namespace {
0058   template <typename T>
0059   edm::Handle<T> getHandle(const edm::Event& iEvent, const edm::EDGetTokenT<T>& token) {
0060     edm::Handle<T> handle;
0061     iEvent.getByToken(token, handle);
0062     return handle;
0063   }
0064 
0065   bool hasDRMatch(const reco::Candidate& cand,
0066                   const std::vector<std::pair<float, float>>& etaPhisToMatch,
0067                   const float maxDR2) {
0068     const float candEta = cand.eta();
0069     const float candPhi = cand.phi();
0070     for (const auto& etaPhi : etaPhisToMatch) {
0071       if (reco::deltaR2(candEta, candPhi, etaPhi.first, etaPhi.second) <= maxDR2) {
0072         return true;
0073       }
0074     }
0075     return false;
0076   }
0077 }  // namespace
0078 
0079 // ------------ method called to produce the data  ------------
0080 void CandMergerCleanOthersByDR::produce(edm::StreamID streamID,
0081                                         edm::Event& iEvent,
0082                                         const edm::EventSetup& iSetup) const {
0083   auto outColl = std::make_unique<std::vector<edm::Ptr<reco::Candidate>>>();
0084 
0085   auto coll1Handle = getHandle(iEvent, coll1Token_);
0086   auto coll2Handle = getHandle(iEvent, coll2Token_);
0087 
0088   std::vector<std::pair<float, float>> coll1EtaPhis;
0089   for (size_t objNr = 0; objNr < coll1Handle->size(); objNr++) {
0090     edm::Ptr<reco::Candidate> objPtr(coll1Handle, objNr);
0091     coll1EtaPhis.push_back({objPtr->eta(), objPtr->phi()});  //just to speed up the DR match
0092     outColl->push_back(objPtr);
0093   }
0094   for (size_t objNr = 0; objNr < coll2Handle->size(); objNr++) {
0095     edm::Ptr<reco::Candidate> objPtr(coll2Handle, objNr);
0096     if (!hasDRMatch(*objPtr, coll1EtaPhis, maxDR2ToClean_)) {
0097       outColl->push_back(objPtr);
0098     }
0099   }
0100 
0101   iEvent.put(std::move(outColl));
0102 }
0103 
0104 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0105 void CandMergerCleanOthersByDR::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0106   edm::ParameterSetDescription desc;
0107   desc.add<edm::InputTag>("coll1", edm::InputTag("egammasForCoreTracking"));
0108   desc.add<edm::InputTag>("coll2", edm::InputTag("jetsForCoreTracking"));
0109   desc.add<double>("maxDRToClean", 0.05);
0110   descriptions.add("candMergerCleanOthersByDR", desc);
0111 }
0112 
0113 //define this as a plug-in
0114 DEFINE_FWK_MODULE(CandMergerCleanOthersByDR);