CMS 3D CMS Logo

GenJetClosestMatchSelectorDefinition.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PFCandProducer_GenJetClosestMatchSelectorDefinition
2 #define PhysicsTools_PFCandProducer_GenJetClosestMatchSelectorDefinition
3 
8 
10 
12 
13 #include <iostream>
14 
16 
17 
20  typedef std::vector< reco::GenJet *> container;
21  typedef container::const_iterator const_iterator;
22 
24 
25  matchTo_ = iC.consumes< edm::View<reco::Candidate> >(cfg.getParameter< edm::InputTag >( "MatchTo" ));
26  }
27 
28  const_iterator begin() const { return selected_.begin(); }
29 
30  const_iterator end() const { return selected_.end(); }
31 
32  void select( const HandleToCollection & hc,
33  const edm::Event & e,
34  const edm::EventSetup& s)
35  {
36 
37  selected_.clear();
38 
40  e.getByToken( matchTo_, matchCandidates);
41 
42 
43  unsigned key=0;
44 
45  // std::cout<<"number of candidates "<<matchCandidates->size()<<std::endl;
46 
48  for( IC ic = matchCandidates->begin();
49  ic!= matchCandidates->end(); ++ic ) {
50 
51  double eta2 = ic->eta();
52  double phi2 = ic->phi();
53 
54  // std::cout<<"cand "<<eta2<<" "<<phi2<<std::endl;
55 
56 
57  // look for the closest gen jet
58  double deltaR2Min = 9999;
59  collection::const_iterator closest = hc->end();
60  for( collection::const_iterator genjet = hc->begin();
61  genjet != hc->end();
62  ++genjet, ++key) {
63 
64  reco::GenJetRef genJetRef(hc, key);
65 
66  // is it matched?
67 
68  double eta1 = genjet->eta();
69  double phi1 = genjet->phi();
70 
71 
72  double deltaR2 = reco::deltaR2(eta1, phi1, eta2, phi2);
73 
74  // std::cout<<" genjet "<<eta1<<" "<<phi1<<" "<<deltaR2<<std::endl;
75 
76  // cut should be a parameter
77  if( deltaR2<deltaR2Min ) {
78  deltaR2Min = deltaR2;
79  closest = genjet;
80  }
81  }
82 
83  if(deltaR2Min<0.01 ) {
84  // std::cout<<deltaR2Min<<std::endl;
85  selected_.push_back( new reco::GenJet(*closest) );
86  }
87  } // end collection iteration
88 
89  // std::cout<<selected_.size()<<std::endl;
90  } // end select()
91 
92  size_t size() const { return selected_.size(); }
93 
94 private:
95  container selected_;
97 };
98 
99 #endif
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
std::vector< GenJet > GenJetCollection
collection of GenJet objects
edm::EDGetTokenT< edm::View< reco::Candidate > > matchTo_
Jets made from MC generator particles.
Definition: GenJet.h:24
GenJetClosestMatchSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
susybsm::HSCParticleCollection hc
Definition: classes.h:25