CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1MuonMatcherAlgo.h
Go to the documentation of this file.
1 #ifndef MuonAnalysis_MuonAssociators_interface_L1MuonMatcherAlgo_h
2 #define MuonAnalysis_MuonAssociators_interface_L1MuonMatcherAlgo_h
3 //
4 //
5 
14 #include <cmath>
25 
27  public:
28  explicit L1MuonMatcherAlgo(const edm::ParameterSet & iConfig) ;
30 
32  void init(const edm::EventSetup &iSetup) ;
33 
36 
39 
42 
46  const PropagateToMuon & propagatorToMuon() const { return prop_; }
47 
50  bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
51  propagated = extrapolate(tk);
52  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
53  }
54 
57  bool match(const reco::Candidate &c, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
58  propagated = extrapolate(c);
59  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
60  }
61 
64  bool match(TrajectoryStateOnSurface & propagated, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi) const ;
65 
69  int match(const reco::Track &tk, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
70  propagated = extrapolate(tk);
71  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
72  }
73 
77  int match(const reco::Candidate &c, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
78  propagated = extrapolate(c);
79  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
80  }
81 
85  int match(TrajectoryStateOnSurface &propagated, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi) const ;
86 
87 
92  template<typename Collection, typename Selector>
93  int matchGeneric(const reco::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
94  propagated = extrapolate(tk);
95  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
96  }
97 
102  template<typename Collection, typename Selector>
103  int matchGeneric(const reco::Candidate &c, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
104  propagated = extrapolate(c);
105  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
106  }
107 
112  template<typename Collection, typename Selector>
113  int matchGeneric(TrajectoryStateOnSurface &propagated, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi) const ;
114 
115 
116  private:
118 
122 
125 
128 };
129 
130 template<typename Collection, typename Selector>
131 int
132 L1MuonMatcherAlgo::matchGeneric(TrajectoryStateOnSurface &propagated, const Collection &l1s, const Selector &sel, float &deltaR, float &deltaPhi) const {
133  typedef typename Collection::value_type obj;
134  int match = -1;
135  double minDeltaPhi = deltaPhi_;
136  double minDeltaR2 = deltaR2_;
137  GlobalPoint pos = propagated.globalPosition();
138  for (int i = 0, n = l1s.size(); i < n; ++i) {
139  const obj &l1 = l1s[i];
140  if (sel(l1)) {
141  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi());
142  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi());
143  if ((fabs(thisDeltaPhi) < deltaPhi_) && (thisDeltaR2 < deltaR2_)) { // check both
144  if (sortByDeltaPhi_ ? (fabs(thisDeltaPhi) < fabs(minDeltaPhi)) : (thisDeltaR2 < minDeltaR2)) { // sort on one
145  match = i;
146  deltaR = std::sqrt(thisDeltaR2);
147  deltaPhi = thisDeltaPhi;
148  if (sortByDeltaPhi_) minDeltaPhi = thisDeltaPhi; else minDeltaR2 = thisDeltaR2;
149  }
150  }
151  }
152  }
153  return match;
154 }
155 
156 #endif
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
int i
Definition: DBlmapReader.cc:9
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &state) const
Extrapolate a FreeTrajectoryState to the muon station 2, return an invalid TSOS if it fails...
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators...
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
int matchGeneric(const reco::Candidate &c, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
PropagateToMuon prop_
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
GlobalPoint globalPosition() const
bool sortByDeltaPhi_
Sort by deltaPhi instead of deltaR.
const PropagateToMuon & propagatorToMuon() const
Return the propagator to second muon station (in case it&#39;s needed)
Propagate an object (usually a track) to the second muon station. Support for other muon stations wil...
int match(const reco::Track &tk, const std::vector< l1extra::L1MuonParticle > &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
int match(const reco::Candidate &c, const std::vector< l1extra::L1MuonParticle > &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
T sqrt(T t)
Definition: SSEVec.h:48
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
Container::value_type value_type
Functor that operates on &lt;T&gt;
Definition: Selector.h:24
PropagateToMuon & propagatorToMuon()
Return the propagator to second muon station (in case it&#39;s needed)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
bool match(const reco::Candidate &c, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
int matchGeneric(const reco::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
T eta() const
Definition: PV3DBase.h:76
StringCutObjectSelector< l1extra::L1MuonParticle > L1Selector
L1MuonMatcherAlgo(const edm::ParameterSet &iConfig)
Matcher of reconstructed objects to L1 Muons.
TrajectoryStateOnSurface extrapolate(const reco::Candidate &tk) const
Extrapolate reco::Candidate to the muon station 2, return an invalid TSOS if it fails.
double deltaR2_
Matching cuts.