CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
L1MuonMatcherAlgoT< Tr > Class Template Reference

#include <L1MuonMatcherAlgo.h>

Public Member Functions

TrajectoryStateOnSurface extrapolate (const reco::Track &tk) const
 Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails. More...
 
TrajectoryStateOnSurface extrapolate (const reco::Candidate &tk) const
 Extrapolate reco::Candidate to the muon station 2, return an invalid TSOS if it fails. More...
 
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 SimVertices to know where to start from. More...
 
TrajectoryStateOnSurface extrapolate (const FreeTrajectoryState &state) const
 Extrapolate a FreeTrajectoryState to the muon station 2, return an invalid TSOS if it fails. More...
 
void init (const edm::EventSetup &iSetup)
 Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators. More...
 
 L1MuonMatcherAlgoT (const edm::ParameterSet &iConfig, edm::ConsumesCollector)
 
bool match (const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
bool match (const reco::Candidate &c, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
bool match (const SimTrack &tk, const edm::SimVertexContainer &vtxs, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
bool match (TrajectoryStateOnSurface &propagated, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi) const
 
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
 
int match (const SimTrack &tk, const edm::SimVertexContainer &vtxs, const std::vector< l1extra::L1MuonParticle > &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
int match (TrajectoryStateOnSurface &propagated, const std::vector< l1extra::L1MuonParticle > &l1, float &deltaR, float &deltaPhi) const
 
int match (const reco::Track &tk, const std::vector< l1t::Muon > &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
int match (const reco::Candidate &c, const std::vector< l1t::Muon > &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
int match (const SimTrack &tk, const edm::SimVertexContainer &vtxs, const std::vector< l1t::Muon > &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
int match (TrajectoryStateOnSurface &propagated, const std::vector< l1t::Muon > &l1, float &deltaR, float &deltaPhi) const
 
template<typename Collection , typename Selector >
int matchGeneric (const reco::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
template<typename Collection , typename Selector >
int matchGeneric (const reco::Candidate &c, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
 
template<typename Collection , typename Selector >
int matchGeneric (TrajectoryStateOnSurface &propagated, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi) const
 
PropagateToMuonpropagatorToMuon ()
 Return the propagator to second muon station (in case it's needed) More...
 
const PropagateToMuonpropagatorToMuon () const
 Return the propagator to second muon station (in case it's needed) More...
 
void setL1PhiOffset (double l1PhiOffset)
 Add this offset to the L1 phi before doing the match, to correct for different scales in L1 vs offline. More...
 
 ~L1MuonMatcherAlgoT ()=default
 

Private Types

typedef StringCutObjectSelector< reco::Candidate, true > L1Selector
 
enum  SortBy {
  SortByDeltaR = 0, SortByDeltaPhi, SortByDeltaEta, SortByPt,
  SortByQual
}
 Sort by deltaPhi or deltaEta instead of deltaR. More...
 

Private Member Functions

template<class T >
int genericQuality (T const &cand) const
 
int genericQuality (l1extra::L1MuonParticle const &cand) const
 
int genericQuality (l1t::Muon const &cand) const
 

Private Attributes

double deltaEta_
 
double deltaPhi_
 
double deltaR2_
 Matching cuts. More...
 
double l1PhiOffset_
 offset to be added to the L1 phi before the match More...
 
L1Selector preselectionCut_
 Preselection cut to apply to L1 candidates before matching. More...
 
PropagateToMuon prop_
 
PropagateToMuonSetupT< Tr > propSetup_
 
SortBy sortBy_
 
bool useStage2L1_
 

Detailed Description

template<edm::Transition Tr>
class L1MuonMatcherAlgoT< Tr >

Definition at line 27 of file L1MuonMatcherAlgo.h.

Member Typedef Documentation

◆ L1Selector

template<edm::Transition Tr>
typedef StringCutObjectSelector<reco::Candidate, true> L1MuonMatcherAlgoT< Tr >::L1Selector
private

Definition at line 249 of file L1MuonMatcherAlgo.h.

Member Enumeration Documentation

◆ SortBy

template<edm::Transition Tr>
enum L1MuonMatcherAlgoT::SortBy
private

Sort by deltaPhi or deltaEta instead of deltaR.

Enumerator
SortByDeltaR 
SortByDeltaPhi 
SortByDeltaEta 
SortByPt 
SortByQual 

Definition at line 257 of file L1MuonMatcherAlgo.h.

Constructor & Destructor Documentation

◆ L1MuonMatcherAlgoT()

template<edm::Transition Tr>
L1MuonMatcherAlgoT< Tr >::L1MuonMatcherAlgoT ( const edm::ParameterSet iConfig,
edm::ConsumesCollector  iCollector 
)
explicit

Definition at line 326 of file L1MuonMatcherAlgo.h.

327  : propSetup_(iConfig, iCollector),
328  useStage2L1_(iConfig.existsAs<bool>("useStage2L1") ? iConfig.getParameter<bool>("useStage2L1") : false),
330  (iConfig.existsAs<std::string>("preselection")) ? iConfig.getParameter<std::string>("preselection") : ""),
331  deltaR2_(std::pow(iConfig.getParameter<double>("maxDeltaR"), 2)),
332  deltaPhi_(iConfig.existsAs<double>("maxDeltaPhi") ? iConfig.getParameter<double>("maxDeltaPhi") : 10),
333  deltaEta_(iConfig.existsAs<double>("maxDeltaEta") ? iConfig.getParameter<double>("maxDeltaEta") : 10),
334  l1PhiOffset_(iConfig.existsAs<double>("l1PhiOffset") ? iConfig.getParameter<double>("l1PhiOffset") : 0) {
335  bool reqQual = iConfig.existsAs<bool>("sortByQual") && iConfig.getParameter<bool>("sortByQual");
336  bool reqPhi = iConfig.existsAs<bool>("sortByDeltaPhi") && iConfig.getParameter<bool>("sortByDeltaPhi");
337  bool reqEta = iConfig.existsAs<bool>("sortByDeltaEta") && iConfig.getParameter<bool>("sortByDeltaEta");
338  bool reqPt = iConfig.existsAs<bool>("sortByPt") && iConfig.getParameter<bool>("sortByPt");
339  std::string sortBy = iConfig.existsAs<std::string>("sortBy") ? iConfig.getParameter<std::string>("sortBy") : "";
340  if (reqPhi + reqEta + reqPt + reqQual > 1)
341  throw cms::Exception("Configuration")
342  << "L1MuonMatcherAlgoT: Can't set more than one 'sortBy<XXX>' parameter to True.\n";
343  if (sortBy == "deltaPhi") {
344  if (reqEta || reqPt || reqQual)
345  throw cms::Exception("Configuration") << "L1MuonMatcherAlgoT: Can't set sortBy = 'deltaPhi' and set also another "
346  "'sortBy<XXX>' parameter to True.\n";
347  else
348  reqPhi = true;
349  }
350  if (sortBy == "deltaEta") {
351  if (reqPhi || reqPt || reqQual)
352  throw cms::Exception("Configuration") << "L1MuonMatcherAlgoT: Can't set sortBy = 'deltaEta' and set also another "
353  "'sortBy<XXX>' parameter to True.\n";
354  else
355  reqEta = true;
356  }
357  if (sortBy == "pt") {
358  if (reqPhi || reqEta || reqQual)
359  throw cms::Exception("Configuration")
360  << "L1MuonMatcherAlgoT: Can't set sortBy = 'pt' and set also another 'sortBy<XXX>' parameter to True.\n";
361  else
362  reqPt = true;
363  }
364  if (sortBy == "quality") {
365  if (reqPhi || reqEta || reqPt)
366  throw cms::Exception("Configuration")
367  << "L1MuonMatcherAlgoT: Can't set sortBy = 'quality' and set also another 'sortBy<XXX>' parameter to True.\n";
368  else
369  reqQual = true;
370  }
371  if (sortBy == "deltaR") {
372  if (reqPhi || reqEta || reqPt || reqQual)
373  throw cms::Exception("Configuration")
374  << "L1MuonMatcherAlgoT: Can't set sortBy = 'deltaR' and set also another 'sortBy<XXX>' parameter to True.\n";
375  }
376  // so, if we're here there's no ambiguity in what the user may want. either everything is false, or exactly one req is true.
377  if (reqEta)
379  else if (reqPhi)
381  else if (reqQual)
383  else if (reqPt)
384  sortBy_ = SortByPt;
385  else
387 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
double l1PhiOffset_
offset to be added to the L1 phi before the match
double deltaR2_
Matching cuts.
PropagateToMuonSetupT< Tr > propSetup_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ ~L1MuonMatcherAlgoT()

template<edm::Transition Tr>
L1MuonMatcherAlgoT< Tr >::~L1MuonMatcherAlgoT ( )
default

Member Function Documentation

◆ extrapolate() [1/4]

template<edm::Transition Tr>
TrajectoryStateOnSurface L1MuonMatcherAlgoT< Tr >::extrapolate ( const reco::Track tk) const
inline

Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.

Definition at line 36 of file L1MuonMatcherAlgo.h.

Referenced by L1MuonMatcherAlgoT< edm::Transition::BeginRun >::match(), and L1MuonMatcherAlgoT< edm::Transition::BeginRun >::matchGeneric().

36 { return prop_.extrapolate(tk); }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
PropagateToMuon prop_

◆ extrapolate() [2/4]

template<edm::Transition Tr>
TrajectoryStateOnSurface L1MuonMatcherAlgoT< Tr >::extrapolate ( const reco::Candidate tk) const
inline

Extrapolate reco::Candidate to the muon station 2, return an invalid TSOS if it fails.

Definition at line 39 of file L1MuonMatcherAlgo.h.

39 { return prop_.extrapolate(tk); }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
PropagateToMuon prop_

◆ extrapolate() [3/4]

template<edm::Transition Tr>
TrajectoryStateOnSurface L1MuonMatcherAlgoT< Tr >::extrapolate ( const SimTrack tk,
const edm::SimVertexContainer vtx 
) const
inline

Extrapolate a SimTrack to the muon station 2, return an invalid TSOS if it fails. Requires SimVertices to know where to start from.

Definition at line 42 of file L1MuonMatcherAlgo.h.

42  {
43  return prop_.extrapolate(tk, vtx);
44  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
PropagateToMuon prop_

◆ extrapolate() [4/4]

template<edm::Transition Tr>
TrajectoryStateOnSurface L1MuonMatcherAlgoT< Tr >::extrapolate ( const FreeTrajectoryState state) const
inline

Extrapolate a FreeTrajectoryState to the muon station 2, return an invalid TSOS if it fails.

Definition at line 47 of file L1MuonMatcherAlgo.h.

47 { return prop_.extrapolate(state); }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
PropagateToMuon prop_

◆ genericQuality() [1/3]

template<edm::Transition Tr>
template<class T >
int L1MuonMatcherAlgoT< Tr >::genericQuality ( T const &  cand) const
inlineprivate

Definition at line 237 of file L1MuonMatcherAlgo.h.

237  {
238  return 0;
239  }

◆ genericQuality() [2/3]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::genericQuality ( l1extra::L1MuonParticle const &  cand) const
inlineprivate

Definition at line 241 of file L1MuonMatcherAlgo.h.

241 { return cand.gmtMuonCand().rank(); }

◆ genericQuality() [3/3]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::genericQuality ( l1t::Muon const &  cand) const
inlineprivate

Definition at line 242 of file L1MuonMatcherAlgo.h.

242 { return cand.hwQual(); }

◆ init()

template<edm::Transition Tr>
void L1MuonMatcherAlgoT< Tr >::init ( const edm::EventSetup iSetup)

Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators.

Definition at line 390 of file L1MuonMatcherAlgo.h.

Referenced by HLTMuonPlotter::beginRun().

390  {
391  prop_ = propSetup_.init(iSetup);
392 }
PropagateToMuonSetupT< Tr > propSetup_
PropagateToMuon prop_

◆ match() [1/12]

template<edm::Transition Tr>
bool L1MuonMatcherAlgoT< Tr >::match ( const reco::Track tk,
const l1extra::L1MuonParticle l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Try to match one track to one L1. Return true if succeeded (and update deltaR, deltaPhi and propagated TSOS accordingly) The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 56 of file L1MuonMatcherAlgo.h.

Referenced by L1MuonMatcherAlgoT< edm::Transition::BeginRun >::match().

60  {
61  propagated = extrapolate(tk);
62  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
63  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [2/12]

template<edm::Transition Tr>
bool L1MuonMatcherAlgoT< Tr >::match ( const reco::Candidate c,
const l1extra::L1MuonParticle l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Try to match one track to one L1. Return true if succeeded (and update deltaR, deltaPhi and propagated TSOS accordingly) The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 67 of file L1MuonMatcherAlgo.h.

71  {
72  propagated = extrapolate(c);
73  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
74  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [3/12]

template<edm::Transition Tr>
bool L1MuonMatcherAlgoT< Tr >::match ( const SimTrack tk,
const edm::SimVertexContainer vtxs,
const l1extra::L1MuonParticle l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Try to match one simtrack to one L1. Return true if succeeded (and update deltaR, deltaPhi and propagated TSOS accordingly) The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 78 of file L1MuonMatcherAlgo.h.

83  {
84  propagated = extrapolate(tk, vtxs);
85  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
86  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [4/12]

template<edm::Transition Tr>
bool L1MuonMatcherAlgoT< Tr >::match ( TrajectoryStateOnSurface propagated,
const l1extra::L1MuonParticle l1,
float &  deltaR,
float &  deltaPhi 
) const

Try to match one track to one L1. Return true if succeeded (and update deltaR, deltaPhi accordingly) The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 395 of file L1MuonMatcherAlgo.h.

398  {
399  if (preselectionCut_(l1)) {
400  GlobalPoint pos = propagated.globalPosition();
401  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi() + l1PhiOffset_);
402  double thisDeltaEta = pos.eta() - l1.eta();
403  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi() + l1PhiOffset_);
404  if ((fabs(thisDeltaPhi) < deltaPhi_) && (fabs(thisDeltaEta) < deltaEta_) && (thisDeltaR2 < deltaR2_)) {
405  deltaR = std::sqrt(thisDeltaR2);
406  deltaPhi = thisDeltaPhi;
407  return true;
408  }
409  }
410  return false;
411 }
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
double l1PhiOffset_
offset to be added to the L1 phi before the match
GlobalPoint globalPosition() const
T sqrt(T t)
Definition: SSEVec.h:19
double deltaR2_
Matching cuts.
double phi() const final
momentum azimuthal angle
double eta() const final
momentum pseudorapidity

◆ match() [5/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( const reco::Track tk,
const std::vector< l1extra::L1MuonParticle > &  l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 100 of file L1MuonMatcherAlgo.h.

104  {
105  propagated = extrapolate(tk);
106  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
107  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [6/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( const reco::Candidate c,
const std::vector< l1extra::L1MuonParticle > &  l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 112 of file L1MuonMatcherAlgo.h.

116  {
117  propagated = extrapolate(c);
118  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
119  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [7/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( const SimTrack tk,
const edm::SimVertexContainer vtxs,
const std::vector< l1extra::L1MuonParticle > &  l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 124 of file L1MuonMatcherAlgo.h.

129  {
130  propagated = extrapolate(tk, vtxs);
131  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
132  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [8/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( TrajectoryStateOnSurface propagated,
const std::vector< l1extra::L1MuonParticle > &  l1,
float &  deltaR,
float &  deltaPhi 
) const

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi accordingly) Returns -1 if the match fails The preselection cut on L1, if specified in the config, is applied before the match

Definition at line 414 of file L1MuonMatcherAlgo.h.

417  {
418  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
419 }
Definition: L1Scalers.h:16
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
int matchGeneric(const reco::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [9/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( const reco::Track tk,
const std::vector< l1t::Muon > &  l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to stage2 L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails The preselection cut on stage2 L1, if specified in the config, is applied before the match

Definition at line 147 of file L1MuonMatcherAlgo.h.

151  {
152  propagated = extrapolate(tk);
153  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
154  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [10/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( const reco::Candidate c,
const std::vector< l1t::Muon > &  l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to stage2 L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails The preselection cut on stage2 L1, if specified in the config, is applied before the match

Definition at line 159 of file L1MuonMatcherAlgo.h.

163  {
164  propagated = extrapolate(c);
165  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
166  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [11/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( const SimTrack tk,
const edm::SimVertexContainer vtxs,
const std::vector< l1t::Muon > &  l1,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to stage2 L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails The preselection cut on stage 2 L1, if specified in the config, is applied before the match

Definition at line 171 of file L1MuonMatcherAlgo.h.

176  {
177  propagated = extrapolate(tk, vtxs);
178  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
179  }
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ match() [12/12]

template<edm::Transition Tr>
int L1MuonMatcherAlgoT< Tr >::match ( TrajectoryStateOnSurface propagated,
const std::vector< l1t::Muon > &  l1,
float &  deltaR,
float &  deltaPhi 
) const

Find the best match to stage 2 L1, and return its index in the vector (and update deltaR, deltaPhi accordingly) Returns -1 if the match fails The preselection cut on stage 2 L1, if specified in the config, is applied before the match

Definition at line 422 of file L1MuonMatcherAlgo.h.

425  {
426  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
427 }
Definition: L1Scalers.h:16
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
int matchGeneric(const reco::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ matchGeneric() [1/3]

template<edm::Transition Tr>
template<typename Collection , typename Selector >
int L1MuonMatcherAlgoT< Tr >::matchGeneric ( const reco::Track tk,
const Collection &  l1,
const Selector sel,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails Only the objects passing the selector will be allowed for the match. If you don't need a selector, just use an AnySelector (CommonTools/Utils) which accepts everything

Definition at line 196 of file L1MuonMatcherAlgo.h.

Referenced by L1MuonMatcherAlgoT< edm::Transition::BeginRun >::matchGeneric().

201  {
202  propagated = extrapolate(tk);
203  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
204  }
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::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ matchGeneric() [2/3]

template<edm::Transition Tr>
template<typename Collection , typename Selector >
int L1MuonMatcherAlgoT< Tr >::matchGeneric ( const reco::Candidate c,
const Collection &  l1,
const Selector sel,
float &  deltaR,
float &  deltaPhi,
TrajectoryStateOnSurface propagated 
) const
inline

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi and propagated TSOS accordingly) Returns -1 if the match fails Only the objects passing the selector will be allowed for the match. If you don't need a selector, just use an AnySelector (CommonTools/Utils) which accepts everything

Definition at line 211 of file L1MuonMatcherAlgo.h.

216  {
217  propagated = extrapolate(c);
218  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
219  }
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::Track &tk, const Collection &l1, const Selector &sel, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const

◆ matchGeneric() [3/3]

template<edm::Transition Tr>
template<typename Collection , typename Selector >
int L1MuonMatcherAlgoT< Tr >::matchGeneric ( TrajectoryStateOnSurface propagated,
const Collection &  l1,
const Selector sel,
float &  deltaR,
float &  deltaPhi 
) const

Find the best match to L1, and return its index in the vector (and update deltaR, deltaPhi accordingly) Returns -1 if the match fails Only the objects passing the selector will be allowed for the match. The selector defaults to an AnySelector (CommonTools/Utils) which just accepts everything

Definition at line 266 of file L1MuonMatcherAlgo.h.

270  {
271  typedef typename Collection::value_type obj;
272 
273  int match = -1;
274  double minDeltaPhi = deltaPhi_;
275  double minDeltaEta = deltaEta_;
276  double minDeltaR2 = deltaR2_;
277  double maxPt = -1.0;
278  int maxQual = 0;
279  GlobalPoint pos = propagated.globalPosition();
280  for (int i = 0, n = l1s.size(); i < n; ++i) {
281  const obj &l1 = l1s[i];
282  if (sel(l1)) {
283  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi() + l1PhiOffset_);
284  double thisDeltaEta = pos.eta() - l1.eta();
285  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi() + l1PhiOffset_);
286  int thisQual = genericQuality<obj>(l1);
287  double thisPt = l1.pt();
288  if ((fabs(thisDeltaPhi) < deltaPhi_) && (fabs(thisDeltaEta) < deltaEta_) &&
289  (thisDeltaR2 < deltaR2_)) { // check all
290  bool betterMatch = (match == -1);
291  switch (sortBy_) {
292  case SortByDeltaR:
293  betterMatch = (thisDeltaR2 < minDeltaR2);
294  break;
295  case SortByDeltaEta:
296  betterMatch = (fabs(thisDeltaEta) < fabs(minDeltaEta));
297  break;
298  case SortByDeltaPhi:
299  betterMatch = (fabs(thisDeltaPhi) < fabs(minDeltaPhi));
300  break;
301  case SortByPt:
302  betterMatch = (thisPt > maxPt);
303  break;
304  // Quality is an int, adding sorting by pT in case of identical qualities
305  case SortByQual:
306  betterMatch = (thisQual > maxQual || ((thisQual == maxQual) && (thisPt > maxPt)));
307  break;
308  }
309  if (betterMatch) { // sort on one
310  match = i;
311  deltaR = std::sqrt(thisDeltaR2);
312  deltaPhi = thisDeltaPhi;
313  minDeltaR2 = thisDeltaR2;
314  minDeltaEta = thisDeltaEta;
315  minDeltaPhi = thisDeltaPhi;
316  maxQual = thisQual;
317  maxPt = thisPt;
318  }
319  }
320  }
321  }
322  return match;
323 }
Definition: L1Scalers.h:16
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
double l1PhiOffset_
offset to be added to the L1 phi before the match
GlobalPoint globalPosition() const
T sqrt(T t)
Definition: SSEVec.h:19
double deltaR2_
Matching cuts.
maxPt
Definition: PV_cfg.py:223

◆ propagatorToMuon() [1/2]

template<edm::Transition Tr>
PropagateToMuon& L1MuonMatcherAlgoT< Tr >::propagatorToMuon ( )
inline

Return the propagator to second muon station (in case it's needed)

Definition at line 50 of file L1MuonMatcherAlgo.h.

50 { return prop_; }
PropagateToMuon prop_

◆ propagatorToMuon() [2/2]

template<edm::Transition Tr>
const PropagateToMuon& L1MuonMatcherAlgoT< Tr >::propagatorToMuon ( ) const
inline

Return the propagator to second muon station (in case it's needed)

Definition at line 52 of file L1MuonMatcherAlgo.h.

52 { return prop_; }
PropagateToMuon prop_

◆ setL1PhiOffset()

template<edm::Transition Tr>
void L1MuonMatcherAlgoT< Tr >::setL1PhiOffset ( double  l1PhiOffset)
inline

Add this offset to the L1 phi before doing the match, to correct for different scales in L1 vs offline.

Definition at line 233 of file L1MuonMatcherAlgo.h.

double l1PhiOffset_
offset to be added to the L1 phi before the match

Member Data Documentation

◆ deltaEta_

template<edm::Transition Tr>
double L1MuonMatcherAlgoT< Tr >::deltaEta_
private

Definition at line 254 of file L1MuonMatcherAlgo.h.

◆ deltaPhi_

template<edm::Transition Tr>
double L1MuonMatcherAlgoT< Tr >::deltaPhi_
private

Definition at line 254 of file L1MuonMatcherAlgo.h.

◆ deltaR2_

template<edm::Transition Tr>
double L1MuonMatcherAlgoT< Tr >::deltaR2_
private

Matching cuts.

Definition at line 254 of file L1MuonMatcherAlgo.h.

◆ l1PhiOffset_

template<edm::Transition Tr>
double L1MuonMatcherAlgoT< Tr >::l1PhiOffset_
private

offset to be added to the L1 phi before the match

Definition at line 261 of file L1MuonMatcherAlgo.h.

Referenced by L1MuonMatcherAlgoT< edm::Transition::BeginRun >::setL1PhiOffset().

◆ preselectionCut_

template<edm::Transition Tr>
L1Selector L1MuonMatcherAlgoT< Tr >::preselectionCut_
private

Preselection cut to apply to L1 candidates before matching.

Definition at line 251 of file L1MuonMatcherAlgo.h.

◆ prop_

template<edm::Transition Tr>
PropagateToMuon L1MuonMatcherAlgoT< Tr >::prop_
private

◆ propSetup_

template<edm::Transition Tr>
PropagateToMuonSetupT<Tr> L1MuonMatcherAlgoT< Tr >::propSetup_
private

Definition at line 244 of file L1MuonMatcherAlgo.h.

◆ sortBy_

template<edm::Transition Tr>
SortBy L1MuonMatcherAlgoT< Tr >::sortBy_
private

◆ useStage2L1_

template<edm::Transition Tr>
bool L1MuonMatcherAlgoT< Tr >::useStage2L1_
private

Definition at line 247 of file L1MuonMatcherAlgo.h.