CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Matchers.h
Go to the documentation of this file.
1 #ifndef __RecoParticleFlow_Benchmark_Matchers__
2 #define __RecoParticleFlow_Benchmark_Matchers__
3 
5 
9 
10 #include <vector>
11 #include <iostream>
12 
13 
14 namespace PFB {
15 
16  template< typename C, typename M >
17  void match( const C& candCollection,
18  const M& matchedCandCollection,
19  std::vector<int>& matchIndices,
20  bool matchCharge = false,
21  float dRMax=-1 ) {
22 
23  // compute distance to each candidate in the matchedCandCollection.
24 
25  float dR2Max = 0;
26  if(dRMax>0) dR2Max = dRMax*dRMax;
27 
28  matchIndices.clear();
29  matchIndices.resize( candCollection.size(), -1);
30 
31  for( unsigned i=0; i<candCollection.size(); ++i) {
32 
33  static const double bigNumber = 1e14;
34  double dR2min = bigNumber;
35  int jMin = -1;
36  for( unsigned jm=0; jm<matchedCandCollection.size(); ++jm) {
37 
38  if( matchCharge &&
39  candCollection[i].charge() != matchedCandCollection[jm].charge() )
40  continue;
41 
42  double dR2 = reco::deltaR2( candCollection[i],
43  matchedCandCollection[jm] );
44  if( dR2<dR2min ) {
45  dR2min = dR2;
46  jMin = jm;
47  }
48  }
49 
50  if( (dR2Max>0 && dR2min < dR2Max) || dRMax<=0 ) {
51  matchIndices[i] = jMin;
52  /* std::cout<<"match "<<dR2min<<std::endl; */
53  }
54  // store the closest match, no cut on deltaR.
55  }
56  }
57 
58  // needed by muon matching
59  //template< typename C, typename M, typename MM>
60  template< typename C, typename M >
61  void match( const C& candCollection,
62  const M& matchedCandCollection,
63  std::vector<int>& matchIndices,
65  //const MM& muonMatchedCandCollection,
66  edm::View<reco::Muon> muonMatchedCandCollection,
67  bool matchCharge = false,
68  float dRMax=-1) {
69 
70  // compute distance to each candidate in the matchedCandCollection.
71 
72  float dR2Max = 0;
73  if(dRMax>0) dR2Max = dRMax*dRMax;
74 
75  matchIndices.clear();
76  matchIndices.resize( candCollection.size(), -1);
77 
78  for( unsigned i=0; i<candCollection.size(); ++i) {
79 
80  static const double bigNumber = 1e14;
81  double dR2min = bigNumber;
82  int jMin = -1;
83  for( unsigned jm=0; jm<matchedCandCollection.size(); ++jm) {
84 
85  if ( parameterSet.getParameter<bool>("slimmedLikeSelection") ) {
86  if ( !( muonMatchedCandCollection[jm].pt() > parameterSet.getParameter<double>("ptBase") ||
87  muonMatchedCandCollection[jm].isPFMuon() ||
88  ( muonMatchedCandCollection[jm].pt() > parameterSet.getParameter<double>("ptNotPF") &&
89  (muonMatchedCandCollection[jm].isGlobalMuon() || muonMatchedCandCollection[jm].isStandAloneMuon() || muonMatchedCandCollection[jm].numberOfMatches() > 0 || muon::isGoodMuon(muonMatchedCandCollection[jm], muon::RPCMuLoose) ) )
90  ) )
91  continue ;
92  }
93 
94  if( matchCharge &&
95  candCollection[i].charge() != matchedCandCollection[jm].charge() )
96  continue ;
97 
98  double dR2 = reco::deltaR2( candCollection[i],
99  matchedCandCollection[jm] );
100  if( dR2<dR2min ) {
101  dR2min = dR2;
102  jMin = jm;
103  }
104  }
105 
106  if( (dR2Max>0 && dR2min < dR2Max) || dRMax<=0 ) {
107  matchIndices[i] = jMin;
108  /* std::cout<<"match "<<dR2min<<std::endl; */
109  }
110  // store the closest match, no cut on deltaR.
111  }
112  }
113 
114 }
115 
116 #endif
T getParameter(std::string const &) const
const double dR2Max
int i
Definition: DBlmapReader.cc:9
void match(const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
Definition: Matchers.h:17
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11