CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1MuonMatcherAlgo.cc
Go to the documentation of this file.
2 
4 
6  prop_(iConfig),
7  preselectionCut_(iConfig.existsAs<std::string>("preselection") ? iConfig.getParameter<std::string>("preselection") : ""),
8  deltaR2_(std::pow(iConfig.getParameter<double>("maxDeltaR"),2)),
9  deltaPhi_(iConfig.existsAs<double>("maxDeltaPhi") ? iConfig.getParameter<double>("maxDeltaPhi") : 10),
10  deltaEta_(iConfig.existsAs<double>("maxDeltaEta") ? iConfig.getParameter<double>("maxDeltaEta") : 10),
11  l1PhiOffset_(iConfig.existsAs<double>("l1PhiOffset") ? iConfig.getParameter<double>("l1PhiOffset") : 0)
12 {
13  bool reqPhi = iConfig.existsAs<bool>("sortByDeltaPhi") && iConfig.getParameter<bool>("sortByDeltaPhi");
14  bool reqEta = iConfig.existsAs<bool>("sortByDeltaEta") && iConfig.getParameter<bool>("sortByDeltaEta");
15  bool reqPt = iConfig.existsAs<bool>("sortByPt") && iConfig.getParameter<bool>("sortByPt");
16  std::string sortBy = iConfig.existsAs<std::string>("sortBy") ? iConfig.getParameter<std::string>("sortBy") : "";
17  if (reqPhi + reqEta + reqPt > 1) throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set more than one 'sortBy<XXX>' parameter to True.\n";
18  if (sortBy == "deltaPhi") {
19  if (reqEta || reqPt)
20  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaPhi' and set also another 'sortBy<XXX>' parameter to True.\n";
21  else reqPhi = true;
22  }
23  if (sortBy == "deltaEta") {
24  if(reqPhi || reqPt)
25  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaEta' and set also another 'sortBy<XXX>' parameter to True.\n";
26  else reqEta = true;
27  }
28  if (sortBy == "pt") {
29  if(reqPhi || reqEta)
30  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'pt' and set also another 'sortBy<XXX>' parameter to True.\n";
31  else reqPt = true;
32  }
33  if (sortBy == "deltaR") {
34  if(reqPhi || reqEta || reqPt)
35  throw cms::Exception("Configuration") << "L1MuonMatcherAlgo: Can't set sortBy = 'deltaR' and set also another 'sortBy<XXX>' parameter to True.\n";
36  }
37  // 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.
38  if (reqEta) sortBy_ = SortByDeltaEta;
39  else if (reqPhi) sortBy_ = SortByDeltaPhi;
40  else if (reqPt) sortBy_ = SortByPt;
41  else sortBy_ = SortByDeltaR;
42 }
43 
44 
46 
47 void
49  prop_.init(iSetup);
50 }
51 
52 bool
53 L1MuonMatcherAlgo::match(TrajectoryStateOnSurface & propagated, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi) const {
54  if (preselectionCut_(l1)) {
55  GlobalPoint pos = propagated.globalPosition();
56  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi()+l1PhiOffset_);
57  double thisDeltaEta = pos.eta() - l1.eta();
58  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi()+l1PhiOffset_);
59  if ((fabs(thisDeltaPhi) < deltaPhi_) && (fabs(thisDeltaEta) < deltaEta_) && (thisDeltaR2 < deltaR2_)) {
60  deltaR = std::sqrt(thisDeltaR2);
61  deltaPhi = thisDeltaPhi;
62  return true;
63  }
64  }
65  return false;
66 }
67 
68 int
69 L1MuonMatcherAlgo::match(TrajectoryStateOnSurface & propagated, const std::vector<l1extra::L1MuonParticle> &l1s, float &deltaR, float &deltaPhi) const {
70  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
71 }
72 
73 
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
T getParameter(std::string const &) const
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators...
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
PropagateToMuon prop_
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
GlobalPoint globalPosition() const
virtual double eta() const
momentum pseudorapidity
T sqrt(T t)
Definition: SSEVec.h:48
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
bool match(const reco::Track &tk, const l1extra::L1MuonParticle &l1, float &deltaR, float &deltaPhi, TrajectoryStateOnSurface &propagated) const
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators...
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
virtual double phi() const
momentum azimuthal angle
double deltaR2_
Matching cuts.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40