Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:15:54

0001 /* \class CandOneToManyDeltaRMatcher
0002  *
0003  * Producer for simple match map:
0004  * class to match two collections of candidate
0005  * with one-to-Many matching
0006  * All elements of class "matched" are matched to each element
0007  * of class "source" orderd in DeltaR
0008  *
0009  */
0010 
0011 #include "FWCore/Framework/interface/global/EDProducer.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
0013 
0014 #include "DataFormats/Candidate/interface/Candidate.h"
0015 
0016 #include <vector>
0017 #include <iostream>
0018 
0019 class CandOneToManyDeltaRMatcher : public edm::global::EDProducer<> {
0020 public:
0021   CandOneToManyDeltaRMatcher(const edm::ParameterSet&);
0022   ~CandOneToManyDeltaRMatcher() override;
0023 
0024 private:
0025   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0026 
0027   edm::EDGetTokenT<reco::CandidateCollection> sourceToken_;
0028   edm::EDGetTokenT<reco::CandidateCollection> matchedToken_;
0029   bool printdebug_;
0030 };
0031 
0032 #include "FWCore/Framework/interface/ESHandle.h"
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/Framework/interface/EventSetup.h"
0035 #include "FWCore/Utilities/interface/InputTag.h"
0036 #include "FWCore/Utilities/interface/EDMException.h"
0037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0039 
0040 #include "DataFormats/Common/interface/Handle.h"
0041 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0042 #include "DataFormats/Candidate/interface/CandMatchMap.h"
0043 #include "DataFormats/Candidate/interface/CandMatchMapMany.h"
0044 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0045 
0046 #include <Math/VectorUtil.h>
0047 #include <TMath.h>
0048 
0049 using namespace edm;
0050 using namespace std;
0051 using namespace reco;
0052 using namespace ROOT::Math::VectorUtil;
0053 
0054 namespace reco {
0055   namespace helper {
0056     typedef pair<size_t, double> MatchPair;
0057 
0058     struct SortBySecond {
0059       bool operator()(const MatchPair& p1, const MatchPair& p2) const { return p1.second < p2.second; }
0060     };
0061   }  // namespace helper
0062 }  // namespace reco
0063 
0064 CandOneToManyDeltaRMatcher::CandOneToManyDeltaRMatcher(const ParameterSet& cfg)
0065     : sourceToken_(consumes<CandidateCollection>(cfg.getParameter<InputTag>("src"))),
0066       matchedToken_(consumes<CandidateCollection>(cfg.getParameter<InputTag>("matched"))),
0067       printdebug_(cfg.getUntrackedParameter<bool>("printDebug", false)) {
0068   produces<CandMatchMapMany>();
0069 }
0070 
0071 CandOneToManyDeltaRMatcher::~CandOneToManyDeltaRMatcher() {}
0072 
0073 void CandOneToManyDeltaRMatcher::produce(edm::StreamID, Event& evt, const EventSetup& es) const {
0074   Handle<CandidateCollection> source;
0075   Handle<CandidateCollection> matched;
0076   evt.getByToken(sourceToken_, source);
0077   evt.getByToken(matchedToken_, matched);
0078 
0079   if (printdebug_) {
0080     for (CandidateCollection::const_iterator c = source->begin(); c != source->end(); ++c) {
0081       cout << "[CandOneToManyDeltaRMatcher] Et source  " << c->et() << endl;
0082     }
0083     for (CandidateCollection::const_iterator c = matched->begin(); c != matched->end(); ++c) {
0084       cout << "[CandOneToManyDeltaRMatcher] Et matched " << c->et() << endl;
0085     }
0086   }
0087 
0088   auto matchMap = std::make_unique<CandMatchMapMany>(
0089       CandMatchMapMany::ref_type(CandidateRefProd(source), CandidateRefProd(matched)));
0090   for (size_t c = 0; c != source->size(); ++c) {
0091     const Candidate& src = (*source)[c];
0092     if (printdebug_)
0093       cout << "[CandOneToManyDeltaRMatcher] source (Et,Eta,Phi) =(" << src.et() << "," << src.eta() << "," << src.phi()
0094            << ")" << endl;
0095     vector<reco::helper::MatchPair> v;
0096     for (size_t m = 0; m != matched->size(); ++m) {
0097       const Candidate& match = (*matched)[m];
0098       double dist = DeltaR(src.p4(), match.p4());
0099       v.push_back(make_pair(m, dist));
0100     }
0101     if (!v.empty()) {
0102       sort(v.begin(), v.end(), reco::helper::SortBySecond());
0103       for (size_t m = 0; m != v.size(); ++m) {
0104         if (printdebug_)
0105           cout << "[CandOneToManyDeltaRMatcher]       match (Et,Eta,Phi) =(" << (*matched)[v[m].first].et() << ","
0106                << (*matched)[v[m].first].eta() << "," << (*matched)[v[m].first].phi() << ") DeltaR=" << v[m].second
0107                << endl;
0108         matchMap->insert(CandidateRef(source, c), make_pair(CandidateRef(matched, v[m].first), v[m].second));
0109       }
0110     }
0111   }
0112 
0113   evt.put(std::move(matchMap));
0114 }
0115 
0116 #include "FWCore/PluginManager/interface/ModuleDef.h"
0117 #include "FWCore/Framework/interface/MakerMacros.h"
0118 
0119 DEFINE_FWK_MODULE(CandOneToManyDeltaRMatcher);