CMS 3D CMS Logo

L1MuonMatcherAlgo.cc
Go to the documentation of this file.
2 
4 
6  prop_(iConfig),
7  useStage2L1_(iConfig.existsAs<bool>("useStage2L1") ? iConfig.getParameter<bool>("useStage2L1") : false),
8  preselectionCut_((iConfig.existsAs<std::string>("preselection")) ? iConfig.getParameter<std::string>("preselection") : ""),
9  deltaR2_(std::pow(iConfig.getParameter<double>("maxDeltaR"),2)),
10  deltaPhi_(iConfig.existsAs<double>("maxDeltaPhi") ? iConfig.getParameter<double>("maxDeltaPhi") : 10),
11  deltaEta_(iConfig.existsAs<double>("maxDeltaEta") ? iConfig.getParameter<double>("maxDeltaEta") : 10),
12  l1PhiOffset_(iConfig.existsAs<double>("l1PhiOffset") ? iConfig.getParameter<double>("l1PhiOffset") : 0)
13 {
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) throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set more than one 'sortBy<XXX>' parameter to True.\n";
20  if (sortBy == "deltaPhi") {
21  if (reqEta || reqPt || reqQual)
22  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaPhi' and set also another 'sortBy<XXX>' parameter to True.\n";
23  else reqPhi = true;
24  }
25  if (sortBy == "deltaEta") {
26  if(reqPhi || reqPt || reqQual)
27  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaEta' and set also another 'sortBy<XXX>' parameter to True.\n";
28  else reqEta = true;
29  }
30  if (sortBy == "pt") {
31  if(reqPhi || reqEta || reqQual)
32  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'pt' and set also another 'sortBy<XXX>' parameter to True.\n";
33  else reqPt = true;
34  }
35  if (sortBy == "quality") {
36  if(reqPhi || reqEta || reqPt)
37  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'quality' and set also another 'sortBy<XXX>' parameter to True.\n";
38  else reqQual = true;
39  }
40  if (sortBy == "deltaR") {
41  if(reqPhi || reqEta || reqPt || reqQual)
42  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaR' and set also another 'sortBy<XXX>' parameter to True.\n";
43  }
44  // 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.
45  if (reqEta) sortBy_ = SortByDeltaEta;
46  else if (reqPhi) sortBy_ = SortByDeltaPhi;
47  else if (reqQual) sortBy_ = SortByQual;
48  else if (reqPt) sortBy_ = SortByPt;
49  else sortBy_ = SortByDeltaR;
50 }
51 
52 
54 
55 void
57  prop_.init(iSetup);
58 }
59 
60 bool
61 L1MuonMatcherAlgo::match(TrajectoryStateOnSurface & propagated, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi) const {
62  if (preselectionCut_(l1)) {
63  GlobalPoint pos = propagated.globalPosition();
64  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi()+l1PhiOffset_);
65  double thisDeltaEta = pos.eta() - l1.eta();
66  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi()+l1PhiOffset_);
67  if ((fabs(thisDeltaPhi) < deltaPhi_) && (fabs(thisDeltaEta) < deltaEta_) && (thisDeltaR2 < deltaR2_)) {
68  deltaR = std::sqrt(thisDeltaR2);
69  deltaPhi = thisDeltaPhi;
70  return true;
71  }
72  }
73  return false;
74 }
75 
76 int
77 L1MuonMatcherAlgo::match(TrajectoryStateOnSurface & propagated, const std::vector<l1extra::L1MuonParticle> &l1s, float &deltaR, float &deltaPhi) const {
78  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
79 }
80 
81 int
82 L1MuonMatcherAlgo::match(TrajectoryStateOnSurface & propagated, const std::vector<l1t::Muon> &l1s, float &deltaR, float &deltaPhi) const {
83  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
84 }
85 
86 
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
T getParameter(std::string const &) const
Definition: L1Scalers.h:17
void init(const edm::EventSetup &iSetup)
double eta() const final
momentum pseudorapidity
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:161
PropagateToMuon prop_
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
GlobalPoint globalPosition() const
T sqrt(T t)
Definition: SSEVec.h:18
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
void init(const edm::EventSetup &iSetup)
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
L1MuonMatcherAlgo(const edm::ParameterSet &iConfig)
double l1PhiOffset_
offset to be added to the L1 phi before the match
double phi() const final
momentum azimuthal angle
double deltaR2_
Matching cuts.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40