Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-05 23:51:54

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