CMS 3D CMS Logo

L1MuonSeedsMerger.cc
Go to the documentation of this file.
2 
5 
7  theDeltaEtaCut = cfg.getParameter<double>("deltaEtaCut");
8  theDiffRelPtCut = cfg.getParameter<double>("diffRelPtCut");
9 }
10 
11 void L1MuonSeedsMerger::resolve(std::vector<TrackAndHits>& tracks) const {
12  sort(tracks.begin(), tracks.end(), Less());
13  typedef std::vector<TrackAndHits>::iterator Tracks_Itr;
14  Tracks_Itr it1 = tracks.begin();
15  while (it1 != tracks.end()) {
16  for (Tracks_Itr it2 = it1 + 1; it1->first && it2 < tracks.end(); it2++) {
17  if (!it2->first)
18  continue;
19  if (it2->first->eta() - it1->first->eta() > theDeltaEtaCut)
20  break;
21  switch (compare(&(*it1), &(*it2))) {
22  case killFirst:
23  delete it1->first;
24  it1->first = nullptr;
25  break;
26  case killSecond:
27  delete it2->first;
28  it2->first = nullptr;
29  break;
30  case mergeTwo:
31  *it2 = *(merge(&(*it1), &(*it2)));
32  it1->first = nullptr;
33  break;
34  case goAhead:
35  default:
36  break;
37  }
38  }
39  if (nullptr == it1->first)
40  tracks.erase(it1);
41  else
42  it1++;
43  }
44 }
45 
47  return (a.first->eta() < b.first->eta());
48 }
49 
51  // temporary algorith, takes track with bigger pt, to be reimplemented
52  if (a->first->pt() > b->first->pt()) {
53  delete b->first;
54  return a;
55  } else {
56  delete a->first;
57  return b;
58  }
59 }
60 
62  int nshared = 0;
63  const SeedingHitSet& hitsA = a->second;
64  const SeedingHitSet& hitsB = b->second;
65  for (unsigned int iHitA = 0, nHitsA = hitsA.size(); iHitA < nHitsA; ++iHitA) {
66  const TrackingRecHit* trhA = hitsA[iHitA]->hit();
67  for (unsigned int iHitB = 0, nHitsB = hitsB.size(); iHitB < nHitsB; ++iHitB) {
68  const TrackingRecHit* trhB = hitsB[iHitB]->hit();
69  if (trhA == trhB)
70  nshared++;
71  }
72  }
73 
74  if (nshared >= 2) {
75  if (hitsA.size() >= 3 && hitsB.size() >= 3)
76  return (a->first->chi2() > b->first->chi2()) ? killFirst : killSecond;
77  else if (hitsB.size() >= 3)
78  return killFirst;
79  else
80  return killSecond;
81  } else if (nshared >= 1) {
82  if (hitsA.size() != hitsB.size())
83  return (hitsA.size() < hitsB.size()) ? killFirst : killSecond;
84  else if (hitsA.size() >= 3 && a->first->charge() == b->first->charge() &&
85  fabs(a->first->pt() - b->first->pt()) / b->first->pt() < theDiffRelPtCut)
86  return (a->first->chi2() > b->first->chi2()) ? killFirst : killSecond;
87  else if (hitsA.size() == 2)
88  return (a->first->pt() < b->first->pt()) ? killFirst : killSecond;
89  else
90  return goAhead;
91  } else
92  return goAhead;
93 }
L1MuonSeedsMerger::mergeTwo
Definition: L1MuonSeedsMerger.h:21
SeedingHitSet
Definition: SeedingHitSet.h:6
L1MuonSeedsMerger.h
L1MuonSeedsMerger::theDiffRelPtCut
float theDiffRelPtCut
Definition: L1MuonSeedsMerger.h:30
L1MuonSeedsMerger::killFirst
Definition: L1MuonSeedsMerger.h:21
L1MuonSeedsMerger::Action
Action
Definition: L1MuonSeedsMerger.h:21
L1MuonSeedsMerger::Less::operator()
bool operator()(const TrackAndHits &, const TrackAndHits &) const
Definition: L1MuonSeedsMerger.cc:46
L1MuonSeedsMerger::Less
Definition: L1MuonSeedsMerger.h:22
Track.h
L1MuonSeedsMerger::compare
Action compare(const TrackAndHits *, const TrackAndHits *) const
Definition: L1MuonSeedsMerger.cc:61
L1MuonSeedsMerger::L1MuonSeedsMerger
L1MuonSeedsMerger(const edm::ParameterSet &cfg)
Definition: L1MuonSeedsMerger.cc:6
L1MuonSeedsMerger::resolve
virtual void resolve(TracksAndHits &) const
Definition: L1MuonSeedsMerger.cc:11
b
double b
Definition: hdecay.h:118
edm::ParameterSet
Definition: ParameterSet.h:47
a
double a
Definition: hdecay.h:119
L1MuonSeedsMerger::TrackAndHits
std::pair< const reco::Track *, SeedingHitSet > TrackAndHits
Definition: L1MuonSeedsMerger.h:14
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:176
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
L1MuonSeedsMerger::merge
const TrackAndHits * merge(const TrackAndHits *, const TrackAndHits *) const
Definition: L1MuonSeedsMerger.cc:50
L1MuonSeedsMerger::killSecond
Definition: L1MuonSeedsMerger.h:21
looper.cfg
cfg
Definition: looper.py:296
TrackingRecHit
Definition: TrackingRecHit.h:21
L1MuonSeedsMerger::goAhead
Definition: L1MuonSeedsMerger.h:21
L1MuonSeedsMerger::theDeltaEtaCut
float theDeltaEtaCut
Definition: L1MuonSeedsMerger.h:29
SeedingHitSet::size
unsigned int size() const
Definition: SeedingHitSet.h:41
ParameterSet.h