Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef PhysicsTools_PFCandProducer_GenJetClosestMatchSelectorDefinition
0002 #define PhysicsTools_PFCandProducer_GenJetClosestMatchSelectorDefinition
0003 
0004 #include "FWCore/Framework/interface/ConsumesCollector.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 
0009 #include "DataFormats/JetReco/interface/GenJet.h"
0010 
0011 #include "DataFormats/Math/interface/deltaR.h"
0012 
0013 #include <iostream>
0014 
0015 struct GenJetClosestMatchSelectorDefinition {
0016   typedef reco::GenJetCollection collection;
0017   typedef edm::Handle<collection> HandleToCollection;
0018   typedef std::vector<reco::GenJet *> container;
0019   typedef container::const_iterator const_iterator;
0020 
0021   GenJetClosestMatchSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC) {
0022     matchTo_ = iC.consumes<edm::View<reco::Candidate>>(cfg.getParameter<edm::InputTag>("MatchTo"));
0023   }
0024 
0025   const_iterator begin() const { return selected_.begin(); }
0026 
0027   const_iterator end() const { return selected_.end(); }
0028 
0029   void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s) {
0030     selected_.clear();
0031 
0032     edm::Handle<edm::View<reco::Candidate>> matchCandidates;
0033     e.getByToken(matchTo_, matchCandidates);
0034 
0035     unsigned key = 0;
0036 
0037     //    std::cout<<"number of candidates
0038     //    "<<matchCandidates->size()<<std::endl;
0039 
0040     typedef edm::View<reco::Candidate>::const_iterator IC;
0041     for (IC ic = matchCandidates->begin(); ic != matchCandidates->end(); ++ic) {
0042       double eta2 = ic->eta();
0043       double phi2 = ic->phi();
0044 
0045       //      std::cout<<"cand "<<eta2<<" "<<phi2<<std::endl;
0046 
0047       // look for the closest gen jet
0048       double deltaR2Min = 9999;
0049       collection::const_iterator closest = hc->end();
0050       for (collection::const_iterator genjet = hc->begin(); genjet != hc->end(); ++genjet, ++key) {
0051         reco::GenJetRef genJetRef(hc, key);
0052 
0053         // is it matched?
0054 
0055         double eta1 = genjet->eta();
0056         double phi1 = genjet->phi();
0057 
0058         double deltaR2 = reco::deltaR2(eta1, phi1, eta2, phi2);
0059 
0060         // std::cout<<"  genjet "<<eta1<<" "<<phi1<<" "<<deltaR2<<std::endl;
0061 
0062         // cut should be a parameter
0063         if (deltaR2 < deltaR2Min) {
0064           deltaR2Min = deltaR2;
0065           closest = genjet;
0066         }
0067       }
0068 
0069       if (deltaR2Min < 0.01) {
0070         // std::cout<<deltaR2Min<<std::endl;
0071         selected_.push_back(new reco::GenJet(*closest));
0072       }
0073     }  // end collection iteration
0074 
0075     // std::cout<<selected_.size()<<std::endl;
0076   }  // end select()
0077 
0078   size_t size() const { return selected_.size(); }
0079 
0080 private:
0081   container selected_;
0082   edm::EDGetTokenT<edm::View<reco::Candidate>> matchTo_;
0083 };
0084 
0085 #endif