CMS 3D CMS Logo

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