CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenJetClosestMatchSelectorDefinition.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_PFCandProducer_GenJetClosestMatchSelectorDefinition
2 #define PhysicsTools_PFCandProducer_GenJetClosestMatchSelectorDefinition
3 
5 
7 
8 #include <iostream>
9 
11 
12 
15  typedef std::vector< reco::GenJet *> container;
16  typedef container::const_iterator const_iterator;
17 
19 
20  matchTo_ = cfg.getParameter< edm::InputTag >( "MatchTo" );
21  }
22 
23  const_iterator begin() const { return selected_.begin(); }
24 
25  const_iterator end() const { return selected_.end(); }
26 
27  void select( const HandleToCollection & hc,
28  const edm::Event & e,
29  const edm::EventSetup& s)
30  {
31 
32  selected_.clear();
33 
34  edm::Handle< edm::View<reco::Candidate> > matchCandidates;
35  e.getByLabel( matchTo_, matchCandidates);
36 
37 
38  unsigned key=0;
39 
40  // std::cout<<"number of candidates "<<matchCandidates->size()<<std::endl;
41 
43  for( IC ic = matchCandidates->begin();
44  ic!= matchCandidates->end(); ++ic ) {
45 
46  double eta2 = ic->eta();
47  double phi2 = ic->phi();
48 
49  // std::cout<<"cand "<<eta2<<" "<<phi2<<std::endl;
50 
51 
52  // look for the closest gen jet
53  double deltaR2Min = 9999;
54  collection::const_iterator closest = hc->end();
55  for( collection::const_iterator genjet = hc->begin();
56  genjet != hc->end();
57  ++genjet, ++key) {
58 
59  reco::GenJetRef genJetRef(hc, key);
60 
61  // is it matched?
62 
63  double eta1 = genjet->eta();
64  double phi1 = genjet->phi();
65 
66 
67  double deltaR2 = reco::deltaR2(eta1, phi1, eta2, phi2);
68 
69  // std::cout<<" genjet "<<eta1<<" "<<phi1<<" "<<deltaR2<<std::endl;
70 
71  // cut should be a parameter
72  if( deltaR2<deltaR2Min ) {
73  deltaR2Min = deltaR2;
74  closest = genjet;
75  }
76  }
77 
78  if(deltaR2Min<0.01 ) {
79  // std::cout<<deltaR2Min<<std::endl;
80  selected_.push_back( new reco::GenJet(*closest) );
81  }
82  } // end collection iteration
83 
84  // std::cout<<selected_.size()<<std::endl;
85  } // end select()
86 
87  size_t size() const { return selected_.size(); }
88 
89 private:
92 };
93 
94 #endif
T getParameter(std::string const &) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::vector< GenJet > GenJetCollection
collection of GenJet objects
GenJetClosestMatchSelectorDefinition(const edm::ParameterSet &cfg)
Jets made from MC generator particles.
Definition: GenJet.h:25
void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:58
list key
Definition: combine.py:13