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  sortByDeltaPhi_(iConfig.existsAs<bool>("sortByDeltaPhi") ? iConfig.getParameter<bool>("sortByDeltaPhi") : false)
11 {
12 }
13 
15 
16 void
18  prop_.init(iSetup);
19 }
20 
21 bool
23  if (preselectionCut_(l1)) {
24  GlobalPoint pos = propagated.globalPosition();
25  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi());
26  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi());
27  if ((fabs(thisDeltaPhi) < deltaPhi_) && (thisDeltaR2 < deltaR2_)) {
28  deltaR = std::sqrt(thisDeltaR2);
29  deltaPhi = thisDeltaPhi;
30  return true;
31  }
32  }
33  return false;
34 }
35 
36 int
37 L1MuonMatcherAlgo::match(TrajectoryStateOnSurface & propagated, const std::vector<l1extra::L1MuonParticle> &l1s, float &deltaR, float &deltaPhi) const {
38  return matchGeneric(propagated, l1s, preselectionCut_, deltaR, deltaPhi);
39 /*
40  int match = -1;
41  double minDeltaPhi = deltaPhi_;
42  double minDeltaR2 = deltaR2_;
43  GlobalPoint pos = propagated.globalPosition();
44  for (int i = 0, n = l1s.size(); i < n; ++i) {
45  const l1extra::L1MuonParticle &l1 = l1s[i];
46  if (preselectionCut_(l1)) {
47  double thisDeltaPhi = ::deltaPhi(double(pos.phi()), l1.phi());
48  double thisDeltaR2 = ::deltaR2(double(pos.eta()), double(pos.phi()), l1.eta(), l1.phi());
49  if ((fabs(thisDeltaPhi) < deltaPhi_) && (thisDeltaR2 < deltaR2_)) { // check both
50  if (sortByDeltaPhi_ ? (fabs(thisDeltaPhi) < fabs(minDeltaPhi)) : (thisDeltaR2 < minDeltaR2)) { // sort on one
51  match = i;
52  deltaR = std::sqrt(thisDeltaR2);
53  deltaPhi = thisDeltaPhi;
54  if (sortByDeltaPhi_) minDeltaPhi = thisDeltaPhi; else minDeltaR2 = thisDeltaR2;
55  }
56  }
57  }
58  }
59  return match;
60 */
61 }
62 
63 
L1Selector preselectionCut_
Preselection cut to apply to L1 candidates before matching.
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators...
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)
volatile std::atomic< bool > shutdown_flag false
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