Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:33

0001 /* \class JetMatcherDR
0002  *
0003  * Producer for association map:
0004  * class to match two collections of jet with one-to-one matching
0005  * All elements of class "matched" are matched to each element of class "source" minimizing DeltaR
0006  *
0007  */
0008 
0009 #include "FWCore/Framework/interface/global/EDProducer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/makeRefToBaseProdFrom.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "DataFormats/Common/interface/AssociationMap.h"
0014 #include "DataFormats/JetReco/interface/JetCollection.h"
0015 #include "DataFormats/Math/interface/deltaR.h"
0016 
0017 class JetMatcherDR : public edm::global::EDProducer<> {
0018 public:
0019   typedef edm::AssociationMap<edm::OneToOne<reco::JetView, reco::JetView> > JetMatchMap;
0020   explicit JetMatcherDR(const edm::ParameterSet& iConfig)
0021       : sourceToken_(consumes<reco::JetView>(iConfig.getParameter<edm::InputTag>("source"))),
0022         matchedToken_(consumes<reco::JetView>(iConfig.getParameter<edm::InputTag>("matched"))) {
0023     produces<JetMatchMap>();
0024   };
0025   ~JetMatcherDR() override{};
0026   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0027 
0028 private:
0029   const edm::EDGetTokenT<reco::JetView> sourceToken_;
0030   const edm::EDGetTokenT<reco::JetView> matchedToken_;
0031 };
0032 
0033 void JetMatcherDR::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0034   const auto& source = iEvent.getHandle(sourceToken_);
0035   const auto& matched = iEvent.getHandle(matchedToken_);
0036   auto matching = std::make_unique<JetMatchMap>(source, matched);
0037   for (size_t i = 0; i < source->size(); i++) {
0038     std::pair<int, float> m(-1, 999.f);
0039     for (size_t j = 0; j < matched->size(); j++) {
0040       const auto dR = deltaR((*source)[i].p4(), (*matched)[j].p4());
0041       if (dR < m.second)
0042         m = {j, dR};
0043     }
0044     matching->insert(source->refAt(i), m.first >= 0 ? matched->refAt(m.first) : reco::JetBaseRef());
0045   }
0046   iEvent.put(std::move(matching));
0047 }
0048 
0049 #include "FWCore/Framework/interface/MakerMacros.h"
0050 DEFINE_FWK_MODULE(JetMatcherDR);