CMS 3D CMS Logo

Functions
PFB Namespace Reference

Functions

template<typename C , typename M >
void match (const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, bool matchCharge=false, float dRMax=-1)
 
template<typename C , typename M >
void match (const C &candCollection, const M &matchedCandCollection, std::vector< int > &matchIndices, const edm::ParameterSet &parameterSet, edm::View< reco::Muon > muonMatchedCandCollection, bool matchCharge=false, float dRMax=-1)
 

Function Documentation

template<typename C , typename M >
void PFB::match ( const C &  candCollection,
const M &  matchedCandCollection,
std::vector< int > &  matchIndices,
bool  matchCharge = false,
float  dRMax = -1 
)

Definition at line 18 of file Matchers.h.

References ALCARECOTkAlJpsiMuMu_cff::charge, reco::deltaR2(), dR2Max, allElectronIsolations_cfi::dRMax, and mps_fire::i.

Referenced by PFCandidateMonitor::fill(), PFJetMonitor::fill(), and PFCandidateManager::fill().

22  {
23 
24  // compute distance to each candidate in the matchedCandCollection.
25 
26  float dR2Max = 0;
27  if(dRMax>0) dR2Max = dRMax*dRMax;
28 
29  matchIndices.clear();
30  matchIndices.resize( candCollection.size(), -1);
31 
32  for( unsigned i=0; i<candCollection.size(); ++i) {
33 
34  static const double bigNumber = 1e14;
35  double dR2min = bigNumber;
36  int jMin = -1;
37  for( unsigned jm=0; jm<matchedCandCollection.size(); ++jm) {
38 
39  if( matchCharge &&
40  candCollection[i].charge() != matchedCandCollection[jm].charge() )
41  continue;
42 
43  double dR2 = reco::deltaR2( candCollection[i],
44  matchedCandCollection[jm] );
45  if( dR2<dR2min ) {
46  dR2min = dR2;
47  jMin = jm;
48  }
49  }
50 
51  if( (dR2Max>0 && dR2min < dR2Max) || dRMax<=0 ) {
52  matchIndices[i] = jMin;
53  /* std::cout<<"match "<<dR2min<<std::endl; */
54  }
55  // store the closest match, no cut on deltaR.
56  }
57  }
const double dR2Max
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
template<typename C , typename M >
void PFB::match ( const C &  candCollection,
const M &  matchedCandCollection,
std::vector< int > &  matchIndices,
const edm::ParameterSet parameterSet,
edm::View< reco::Muon muonMatchedCandCollection,
bool  matchCharge = false,
float  dRMax = -1 
)

Definition at line 62 of file Matchers.h.

References ALCARECOTkAlJpsiMuMu_cff::charge, reco::deltaR2(), dR2Max, allElectronIsolations_cfi::dRMax, edm::ParameterSet::getParameter(), mps_fire::i, muon::isGoodMuon(), EnergyCorrector::pt, and muon::RPCMuLoose.

69  {
70 
71  // compute distance to each candidate in the matchedCandCollection.
72 
73  float dR2Max = 0;
74  if(dRMax>0) dR2Max = dRMax*dRMax;
75 
76  matchIndices.clear();
77  matchIndices.resize( candCollection.size(), -1);
78 
79  for( unsigned i=0; i<candCollection.size(); ++i) {
80 
81  static const double bigNumber = 1e14;
82  double dR2min = bigNumber;
83  int jMin = -1;
84  for( unsigned jm=0; jm<matchedCandCollection.size(); ++jm) {
85 
86  if ( parameterSet.getParameter<bool>("slimmedLikeSelection") ) {
87  if ( !( muonMatchedCandCollection[jm].pt() > parameterSet.getParameter<double>("ptBase") ||
88  muonMatchedCandCollection[jm].isPFMuon() ||
89  ( muonMatchedCandCollection[jm].pt() > parameterSet.getParameter<double>("ptNotPF") &&
90  (muonMatchedCandCollection[jm].isGlobalMuon() || muonMatchedCandCollection[jm].isStandAloneMuon() || muonMatchedCandCollection[jm].numberOfMatches() > 0 || muon::isGoodMuon(muonMatchedCandCollection[jm], muon::RPCMuLoose) ) )
91  ) )
92  continue ;
93  }
94 
95  if( matchCharge &&
96  candCollection[i].charge() != matchedCandCollection[jm].charge() )
97  continue ;
98 
99  double dR2 = reco::deltaR2( candCollection[i],
100  matchedCandCollection[jm] );
101  if( dR2<dR2min ) {
102  dR2min = dR2;
103  jMin = jm;
104  }
105  }
106 
107  if( (dR2Max>0 && dR2min < dR2Max) || dRMax<=0 ) {
108  matchIndices[i] = jMin;
109  /* std::cout<<"match "<<dR2min<<std::endl; */
110  }
111  // store the closest match, no cut on deltaR.
112  }
113  }
T getParameter(std::string const &) const
const double dR2Max
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36