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 
15 #include <cmath>
26 
27 class L1MuonMatcherAlgo {
28  public:
29  explicit L1MuonMatcherAlgo(const edm::ParameterSet & iConfig) ;
31 
33  void init(const edm::EventSetup &iSetup) ;
34 
37 
40 
42  TrajectoryStateOnSurface extrapolate(const SimTrack &tk, const edm::SimVertexContainer &vtx) const { return prop_.extrapolate(tk, vtx); }
43 
46 
50  const PropagateToMuon & propagatorToMuon() const { return prop_; }
51 
54  bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
55  propagated = extrapolate(tk);
56  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
57  }
58 
61  bool match(const reco::Candidate &c, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
62  propagated = extrapolate(c);
63  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
64  }
65 
68  bool match(const SimTrack &tk, const edm::SimVertexContainer &vtxs, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
69  propagated = extrapolate(tk, vtxs);
70  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
71  }
72 
73 
76  bool match(TrajectoryStateOnSurface & propagated, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi) const ;
77 
81  int match(const reco::Track &tk, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
82  propagated = extrapolate(tk);
83  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
84  }
85 
89  int match(const reco::Candidate &c, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
90  propagated = extrapolate(c);
91  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
92  }
93 
97  int match(const SimTrack &tk, const edm::SimVertexContainer &vtxs, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
98  propagated = extrapolate(tk, vtxs);
99  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
100  }
101 
102 
106  int match(TrajectoryStateOnSurface &propagated, const std::vector<l1extra::L1MuonParticle> &l1, float &deltaR, float &deltaPhi) const ;
107 
108 
113  template<typename Collection, typename Selector>
114  int matchGeneric(const reco::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
115  propagated = extrapolate(tk);
116  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
117  }
118 
123  template<typename Collection, typename Selector>
124  int matchGeneric(const reco::Candidate &c, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const {
125  propagated = extrapolate(c);
126  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
127  }
128 
133  template<typename Collection, typename Selector>
134  int matchGeneric(TrajectoryStateOnSurface &propagated, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi) const ;
135 
136 
138  void setL1PhiOffset(double l1PhiOffset) { l1PhiOffset_ = l1PhiOffset; }
139 
140  private:
142 
146 
149 
153 
155  double l1PhiOffset_;
156 };
157 
158 template<typename Collection, typename Selector>
159 int
160 L1MuonMatcherAlgo::matchGeneric(TrajectoryStateOnSurface &propagated, const Collection &l1s, const Selector &sel, float &deltaR, float &deltaPhi) const {
161  typedef typename Collection::value_type obj;
162  int match = -1;
163  double minDeltaPhi = deltaPhi_;
164  double minDeltaEta = deltaEta_;
165  double minDeltaR2 = deltaR2_;
166  double minPt = -1.0;
167  GlobalPoint pos = propagated.globalPosition();
168  for (int i = 0, n = l1s.size(); i < n; ++i) {
169  const obj &l1 = l1s[i];
170  if (sel(l1)) {
171  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi()+l1PhiOffset_);
172  double thisDeltaEta = pos.eta() - l1.eta();
173  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi()+l1PhiOffset_);
174  double thisPt = l1.pt();
175  if ((fabs(thisDeltaPhi) < deltaPhi_) && (fabs(thisDeltaEta) < deltaEta_) && (thisDeltaR2 < deltaR2_)) { // check all
176  bool betterMatch = (match == -1);
177  switch (sortBy_) {
178  case SortByDeltaR: betterMatch = (thisDeltaR2 < minDeltaR2); break;
179  case SortByDeltaEta: betterMatch = (fabs(thisDeltaEta) < fabs(minDeltaEta)); break;
180  case SortByDeltaPhi: betterMatch = (fabs(thisDeltaPhi) < fabs(minDeltaPhi)); break;
181  case SortByPt: betterMatch = (thisPt > minPt); break;
182  }
183  if (betterMatch) { // sort on one
184  match = i;
185  deltaR = std::sqrt(thisDeltaR2);
186  deltaPhi = thisDeltaPhi;
187  minDeltaR2 = thisDeltaR2;
188  minDeltaEta = thisDeltaEta;
189  minDeltaPhi = thisDeltaPhi;
190  minPt = thisPt;
191  }
192  }
193  }
194  }
195  return match;
196 }
197 
198 #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
int match(const SimTrack &tk, const edm::SimVertexContainer &vtxs, const std::vector< l1extra::L1MuonParticle > &l1, 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
TrajectoryStateOnSurface extrapolate(const SimTrack &tk, const edm::SimVertexContainer &vtx) const
Extrapolate a SimTrack to the muon station 2, return an invalid TSOS if it fails. Requires SimVertice...
GlobalPoint globalPosition() const
const PropagateToMuon & propagatorToMuon() const
Return the propagator to second muon station (in case it&#39;s needed)
bool match(const SimTrack &tk, const edm::SimVertexContainer &vtxs, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
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)
SortBy
Sort by deltaPhi or deltaEta instead of deltaR.
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
void setL1PhiOffset(double l1PhiOffset)
Add this offset to the L1 phi before doing the match, to correct for different scales in L1 vs offlin...
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
std::vector< SimVertex > SimVertexContainer
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 l1PhiOffset_
offset to be added to the L1 phi before the match
double deltaR2_
Matching cuts.