CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
L1MuonMatcherAlgo Class Reference

Matcher of reconstructed objects to L1 Muons. More...

#include "MuonAnalysis/MuonAssociators/interface/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...
 
 L1MuonMatcherAlgo (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...
 
 ~L1MuonMatcherAlgo ()
 

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 (const T &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_
 
SortBy sortBy_
 
bool useStage2L1_
 

Detailed Description

Matcher of reconstructed objects to L1 Muons.

"HLTriggerOffline/Muon/interface/L1MuonMatcherAlgo.h"

Author
Giovanni Petrucciani
Giovanni Petrucciani
Version
Id:
L1MuonMatcherAlgo.h,v 1.8 2010/07/01 07:40:52 gpetrucc Exp

Definition at line 29 of file L1MuonMatcherAlgo.h.

Member Typedef Documentation

Definition at line 245 of file L1MuonMatcherAlgo.h.

Member Enumeration Documentation

Sort by deltaPhi or deltaEta instead of deltaR.

Enumerator
SortByDeltaR 
SortByDeltaPhi 
SortByDeltaEta 
SortByPt 
SortByQual 

Definition at line 253 of file L1MuonMatcherAlgo.h.

Constructor & Destructor Documentation

L1MuonMatcherAlgo::L1MuonMatcherAlgo ( const edm::ParameterSet iConfig,
edm::ConsumesCollector  iCollector 
)
explicit

Definition at line 5 of file L1MuonMatcherAlgo.cc.

References Exception, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), sortBy_, SortByDeltaEta, SortByDeltaPhi, SortByDeltaR, SortByPt, SortByQual, and AlCaHLTBitMon_QueryRunRegistry::string.

6  : prop_(iConfig, iCollector),
7  useStage2L1_(iConfig.existsAs<bool>("useStage2L1") ? iConfig.getParameter<bool>("useStage2L1") : false),
9  (iConfig.existsAs<std::string>("preselection")) ? iConfig.getParameter<std::string>("preselection") : ""),
10  deltaR2_(std::pow(iConfig.getParameter<double>("maxDeltaR"), 2)),
11  deltaPhi_(iConfig.existsAs<double>("maxDeltaPhi") ? iConfig.getParameter<double>("maxDeltaPhi") : 10),
12  deltaEta_(iConfig.existsAs<double>("maxDeltaEta") ? iConfig.getParameter<double>("maxDeltaEta") : 10),
13  l1PhiOffset_(iConfig.existsAs<double>("l1PhiOffset") ? iConfig.getParameter<double>("l1PhiOffset") : 0) {
14  bool reqQual = iConfig.existsAs<bool>("sortByQual") && iConfig.getParameter<bool>("sortByQual");
15  bool reqPhi = iConfig.existsAs<bool>("sortByDeltaPhi") && iConfig.getParameter<bool>("sortByDeltaPhi");
16  bool reqEta = iConfig.existsAs<bool>("sortByDeltaEta") && iConfig.getParameter<bool>("sortByDeltaEta");
17  bool reqPt = iConfig.existsAs<bool>("sortByPt") && iConfig.getParameter<bool>("sortByPt");
18  std::string sortBy = iConfig.existsAs<std::string>("sortBy") ? iConfig.getParameter<std::string>("sortBy") : "";
19  if (reqPhi + reqEta + reqPt + reqQual > 1)
20  throw cms::Exception("Configuration")
21  << "L1MuonMatcherAlgo: Can't set more than one 'sortBy<XXX>' parameter to True.\n";
22  if (sortBy == "deltaPhi") {
23  if (reqEta || reqPt || reqQual)
24  throw cms::Exception("Configuration")
25  << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaPhi' and set also another 'sortBy<XXX>' parameter to True.\n";
26  else
27  reqPhi = true;
28  }
29  if (sortBy == "deltaEta") {
30  if (reqPhi || reqPt || reqQual)
31  throw cms::Exception("Configuration")
32  << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaEta' and set also another 'sortBy<XXX>' parameter to True.\n";
33  else
34  reqEta = true;
35  }
36  if (sortBy == "pt") {
37  if (reqPhi || reqEta || reqQual)
38  throw cms::Exception("Configuration")
39  << "L1MuonMatcherAlgo: Can't set sortBy = 'pt' and set also another 'sortBy<XXX>' parameter to True.\n";
40  else
41  reqPt = true;
42  }
43  if (sortBy == "quality") {
44  if (reqPhi || reqEta || reqPt)
45  throw cms::Exception("Configuration")
46  << "L1MuonMatcherAlgo: Can't set sortBy = 'quality' and set also another 'sortBy<XXX>' parameter to True.\n";
47  else
48  reqQual = true;
49  }
50  if (sortBy == "deltaR") {
51  if (reqPhi || reqEta || reqPt || reqQual)
52  throw cms::Exception("Configuration")
53  << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaR' and set also another 'sortBy<XXX>' parameter to True.\n";
54  }
55  // 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.
56  if (reqEta)
58  else if (reqPhi)
60  else if (reqQual)
62  else if (reqPt)
63  sortBy_ = SortByPt;
64  else
66 }
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
PropagateToMuon prop_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double l1PhiOffset_
offset to be added to the L1 phi before the match
double deltaR2_
Matching cuts.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
L1MuonMatcherAlgo::~L1MuonMatcherAlgo ( )

Definition at line 68 of file L1MuonMatcherAlgo.cc.

68 {}

Member Function Documentation

TrajectoryStateOnSurface L1MuonMatcherAlgo::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 38 of file L1MuonMatcherAlgo.h.

References PropagateToMuon::extrapolate(), and prop_.

Referenced by match(), and matchGeneric().

38 { 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_
TrajectoryStateOnSurface L1MuonMatcherAlgo::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 41 of file L1MuonMatcherAlgo.h.

References PropagateToMuon::extrapolate(), and prop_.

41 { 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_
TrajectoryStateOnSurface L1MuonMatcherAlgo::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 44 of file L1MuonMatcherAlgo.h.

References PropagateToMuon::extrapolate(), and prop_.

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

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

Definition at line 49 of file L1MuonMatcherAlgo.h.

References PropagateToMuon::extrapolate(), and prop_.

49 { 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_
template<class T >
int L1MuonMatcherAlgo::genericQuality ( const T cand) const
inlineprivate

Definition at line 271 of file L1MuonMatcherAlgo.h.

271  {
272  return 0;
273 }
void L1MuonMatcherAlgo::init ( const edm::EventSetup iSetup)

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

Definition at line 70 of file L1MuonMatcherAlgo.cc.

References PropagateToMuon::init(), and prop_.

Referenced by pat::HLTL1MuonMatcher::produce().

70 { prop_.init(iSetup); }
PropagateToMuon prop_
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each event, to initialize geometry, magnetic field and propagato...
bool L1MuonMatcherAlgo::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 58 of file L1MuonMatcherAlgo.h.

References extrapolate(), and TrajectoryStateOnSurface::isValid().

Referenced by match(), and matchGeneric().

62  {
63  propagated = extrapolate(tk);
64  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
65  }
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
bool L1MuonMatcherAlgo::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 69 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

73  {
74  propagated = extrapolate(c);
75  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
76  }
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
bool L1MuonMatcherAlgo::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 80 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

85  {
86  propagated = extrapolate(tk, vtxs);
87  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : false;
88  }
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
bool L1MuonMatcherAlgo::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 72 of file L1MuonMatcherAlgo.cc.

References deltaEta_, srCondWrite_cfg::deltaPhi, deltaPhi_, reco::deltaR2(), deltaR2_, PV3DBase< T, PVType, FrameType >::eta(), reco::LeafCandidate::eta(), TrajectoryStateOnSurface::globalPosition(), l1PhiOffset_, PV3DBase< T, PVType, FrameType >::phi(), reco::LeafCandidate::phi(), preselectionCut_, and mathSSE::sqrt().

75  {
76  if (preselectionCut_(l1)) {
77  GlobalPoint pos = propagated.globalPosition();
78  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi() + l1PhiOffset_);
79  double thisDeltaEta = pos.eta() - l1.eta();
80  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi() + l1PhiOffset_);
81  if ((fabs(thisDeltaPhi) < deltaPhi_) && (fabs(thisDeltaEta) < deltaEta_) && (thisDeltaR2 < deltaR2_)) {
82  deltaR = std::sqrt(thisDeltaR2);
83  deltaPhi = thisDeltaPhi;
84  return true;
85  }
86  }
87  return false;
88 }
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
GlobalPoint globalPosition() const
T sqrt(T t)
Definition: SSEVec.h:19
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
T eta() const
Definition: PV3DBase.h:73
double l1PhiOffset_
offset to be added to the L1 phi before the match
double phi() const final
momentum azimuthal angle
double deltaR2_
Matching cuts.
double eta() const final
momentum pseudorapidity
int L1MuonMatcherAlgo::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 102 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

106  {
107  propagated = extrapolate(tk);
108  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
109  }
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
int L1MuonMatcherAlgo::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 114 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

118  {
119  propagated = extrapolate(c);
120  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
121  }
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
int L1MuonMatcherAlgo::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 126 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

131  {
132  propagated = extrapolate(tk, vtxs);
133  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
134  }
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
int L1MuonMatcherAlgo::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 90 of file L1MuonMatcherAlgo.cc.

References matchGeneric(), and preselectionCut_.

93  {
94  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
95 }
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
int L1MuonMatcherAlgo::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 149 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

153  {
154  propagated = extrapolate(tk);
155  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
156  }
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
int L1MuonMatcherAlgo::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 161 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

165  {
166  propagated = extrapolate(c);
167  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
168  }
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
int L1MuonMatcherAlgo::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 173 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and match().

178  {
179  propagated = extrapolate(tk, vtxs);
180  return propagated.isValid() ? match(propagated, l1, deltaR, deltaPhi) : -1;
181  }
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
int L1MuonMatcherAlgo::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 97 of file L1MuonMatcherAlgo.cc.

References matchGeneric(), and preselectionCut_.

100  {
101  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
102 }
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
template<typename Collection , typename Selector >
int L1MuonMatcherAlgo::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 198 of file L1MuonMatcherAlgo.h.

References extrapolate(), and TrajectoryStateOnSurface::isValid().

Referenced by match(), matchGeneric(), and pat::HLTL1MuonMatcher::produce().

203  {
204  propagated = extrapolate(tk);
205  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
206  }
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
template<typename Collection , typename Selector >
int L1MuonMatcherAlgo::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 213 of file L1MuonMatcherAlgo.h.

References extrapolate(), TrajectoryStateOnSurface::isValid(), and matchGeneric().

218  {
219  propagated = extrapolate(c);
220  return propagated.isValid() ? matchGeneric(propagated, l1, sel, deltaR, deltaPhi) : -1;
221  }
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
template<typename Collection , typename Selector >
int L1MuonMatcherAlgo::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 276 of file L1MuonMatcherAlgo.h.

References deltaEta_, srCondWrite_cfg::deltaPhi, deltaPhi_, reco::deltaR2(), deltaR2_, PV3DBase< T, PVType, FrameType >::eta(), TrajectoryStateOnSurface::globalPosition(), mps_fire::i, l1PhiOffset_, match(), HLT_FULL_cff::minDeltaEta, dqmiodumpmetadata::n, getGTfromDQMFile::obj, PV3DBase< T, PVType, FrameType >::phi(), EgammaValidation_Wenu_cff::sel, sortBy_, SortByDeltaEta, SortByDeltaPhi, SortByDeltaR, SortByQual, and mathSSE::sqrt().

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

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

Definition at line 52 of file L1MuonMatcherAlgo.h.

References prop_.

52 { return prop_; }
PropagateToMuon prop_
const PropagateToMuon& L1MuonMatcherAlgo::propagatorToMuon ( ) const
inline

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

Definition at line 54 of file L1MuonMatcherAlgo.h.

References prop_.

54 { return prop_; }
PropagateToMuon prop_
void L1MuonMatcherAlgo::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 235 of file L1MuonMatcherAlgo.h.

References l1PhiOffset_.

235 { l1PhiOffset_ = l1PhiOffset; }
double l1PhiOffset_
offset to be added to the L1 phi before the match

Member Data Documentation

double L1MuonMatcherAlgo::deltaEta_
private

Definition at line 250 of file L1MuonMatcherAlgo.h.

Referenced by match(), and matchGeneric().

double L1MuonMatcherAlgo::deltaPhi_
private

Definition at line 250 of file L1MuonMatcherAlgo.h.

Referenced by match(), and matchGeneric().

double L1MuonMatcherAlgo::deltaR2_
private

Matching cuts.

Definition at line 250 of file L1MuonMatcherAlgo.h.

Referenced by match(), and matchGeneric().

double L1MuonMatcherAlgo::l1PhiOffset_
private

offset to be added to the L1 phi before the match

Definition at line 257 of file L1MuonMatcherAlgo.h.

Referenced by match(), matchGeneric(), and setL1PhiOffset().

L1Selector L1MuonMatcherAlgo::preselectionCut_
private

Preselection cut to apply to L1 candidates before matching.

Definition at line 247 of file L1MuonMatcherAlgo.h.

Referenced by match().

PropagateToMuon L1MuonMatcherAlgo::prop_
private

Definition at line 241 of file L1MuonMatcherAlgo.h.

Referenced by extrapolate(), init(), and propagatorToMuon().

SortBy L1MuonMatcherAlgo::sortBy_
private

Definition at line 254 of file L1MuonMatcherAlgo.h.

Referenced by L1MuonMatcherAlgo(), and matchGeneric().

bool L1MuonMatcherAlgo::useStage2L1_
private

Definition at line 243 of file L1MuonMatcherAlgo.h.