CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
SeedFinder Class Reference

#include <SeedFinder.h>

Public Member Functions

void addHitSelector (SeedFinderSelector *seedFinderSelector, unsigned int nHits)
 
std::vector< unsigned int > getSeed (const std::vector< const FastTrackerRecHit *> &trackerRecHits) const
 
const SeedingNode< TrackingLayer > * insertHit (const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
 
std::vector< unsigned int > iterateHits (unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
 
 SeedFinder (const SeedingTree< TrackingLayer > &seedingTree, const TrackerTopology &trackerTopology)
 

Private Attributes

const SeedingTree< TrackingLayer > & _seedingTree
 
std::vector< std::vector< SeedFinderSelector * > > _selectorFunctionsByHits
 
const TrackerTopology_trackerTopology
 

Detailed Description

Definition at line 17 of file SeedFinder.h.

Constructor & Destructor Documentation

◆ SeedFinder()

SeedFinder::SeedFinder ( const SeedingTree< TrackingLayer > &  seedingTree,
const TrackerTopology trackerTopology 
)
inline

Definition at line 25 of file SeedFinder.h.

26  : _seedingTree(seedingTree), _trackerTopology(&trackerTopology) {}
const SeedingTree< TrackingLayer > & _seedingTree
Definition: SeedFinder.h:21
const TrackerTopology * _trackerTopology
Definition: SeedFinder.h:22

Member Function Documentation

◆ addHitSelector()

void SeedFinder::addHitSelector ( SeedFinderSelector seedFinderSelector,
unsigned int  nHits 
)
inline

Definition at line 28 of file SeedFinder.h.

References _selectorFunctionsByHits, nHits, and TrajectorySeedProducer_cfi::seedFinderSelector.

28  {
29  if (_selectorFunctionsByHits.size() < nHits) {
32  }
33  //shift indices by -1 so that _selectorFunctionsByHits[0] tests 1 hit
35  }
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
std::vector< std::vector< SeedFinderSelector * > > _selectorFunctionsByHits
Definition: SeedFinder.h:20

◆ getSeed()

std::vector<unsigned int> SeedFinder::getSeed ( const std::vector< const FastTrackerRecHit *> &  trackerRecHits) const
inline

Definition at line 37 of file SeedFinder.h.

References _seedingTree, _trackerTopology, iterateHits(), and SeedingTree< DATA >::numberOfNodes().

37  {
38  std::vector<int> hitIndicesInTree(_seedingTree.numberOfNodes(), -1);
39  //A SeedingNode is associated by its index to this list. The list stores the indices of the hits in 'trackerRecHits'
40  /* example
41  SeedingNode | hit index | hit
42  -------------------------------------------------------------------------------
43  index= 0: [BPix1] | hitIndicesInTree[0] (=1) | trackerRecHits[1]
44  index= 1: -- [BPix2] | hitIndicesInTree[1] (=3) | trackerRecHits[3]
45  index= 2: -- -- [BPix3] | hitIndicesInTree[2] (=4) | trackerRecHits[4]
46  index= 3: -- -- [FPix1_pos] | hitIndicesInTree[3] (=6) | trackerRecHits[6]
47  index= 4: -- -- [FPix1_neg] | hitIndicesInTree[4] (=7) | trackerRecHits[7]
48 
49  The implementation has been chosen such that the tree only needs to be build once upon construction.
50  */
51  std::vector<TrajectorySeedHitCandidate> seedHitCandidates;
52  for (const FastTrackerRecHit* trackerRecHit : trackerRecHits) {
53  TrajectorySeedHitCandidate seedHitCandidate(trackerRecHit, _trackerTopology);
54  seedHitCandidates.push_back(seedHitCandidate);
55  }
56  return iterateHits(0, seedHitCandidates, hitIndicesInTree, true);
57 
58  //TODO: create pairs of TrackingLayer -> remove TrajectorySeedHitCandidate class
59  }
std::vector< unsigned int > iterateHits(unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
Definition: SeedFinder.h:117
unsigned int numberOfNodes() const
Definition: SeedingTree.h:157
const SeedingTree< TrackingLayer > & _seedingTree
Definition: SeedFinder.h:21
const TrackerTopology * _trackerTopology
Definition: SeedFinder.h:22

◆ insertHit()

const SeedingNode<TrackingLayer>* SeedFinder::insertHit ( const std::vector< TrajectorySeedHitCandidate > &  trackerRecHits,
std::vector< int > &  hitIndicesInTree,
const SeedingNode< TrackingLayer > *  node,
unsigned int  trackerHit 
) const
inline

Definition at line 63 of file SeedFinder.h.

References _selectorFunctionsByHits, relativeConstraints::empty, SeedingNode< DATA >::getChild(), SeedingNode< DATA >::getChildrenSize(), SeedingNode< DATA >::getData(), SeedingNode< DATA >::getDepth(), SeedingNode< DATA >::getIndex(), SeedingNode< DATA >::getParent(), TrajectorySeedHitCandidate::getTrackingLayer(), TrajectorySeedHitCandidate::hit(), wplusjetsAnalysis_cfi::NHits, and fileCollector::seed.

Referenced by iterateHits().

66  {
67  if (!node->getParent() || hitIndicesInTree[node->getParent()->getIndex()] >= 0) {
68  if (hitIndicesInTree[node->getIndex()] < 0) {
69  const TrajectorySeedHitCandidate& currentTrackerHit = trackerRecHits[trackerHit];
70  if (currentTrackerHit.getTrackingLayer() != node->getData()) {
71  return nullptr;
72  }
73 
74  const unsigned int NHits = node->getDepth() + 1;
75  if (_selectorFunctionsByHits.size() >= NHits) {
76  //are there any selector functions stored for NHits?
77  if (!_selectorFunctionsByHits[NHits - 1].empty()) {
78  //fill vector of Hits from node to root to be passed to the selector function
79  std::vector<const FastTrackerRecHit*> seedCandidateHitList(node->getDepth() + 1);
80  seedCandidateHitList[node->getDepth()] = currentTrackerHit.hit();
81  const SeedingNode<TrackingLayer>* parentNode = node->getParent();
82  while (parentNode != nullptr) {
83  seedCandidateHitList[parentNode->getDepth()] =
84  trackerRecHits[hitIndicesInTree[parentNode->getIndex()]].hit();
85  parentNode = parentNode->getParent();
86  }
87 
88  //loop over selector functions
89  for (SeedFinderSelector* selectorFunction : _selectorFunctionsByHits[NHits - 1]) {
90  if (!selectorFunction->pass(seedCandidateHitList)) {
91  return nullptr;
92  }
93  }
94  }
95  }
96 
97  //the hit was not rejected by all selector functions -> insert it into the tree
98  hitIndicesInTree[node->getIndex()] = trackerHit;
99  if (node->getChildrenSize() == 0) {
100  return node;
101  }
102 
103  return nullptr;
104  } else {
105  for (unsigned int ichild = 0; ichild < node->getChildrenSize(); ++ichild) {
107  insertHit(trackerRecHits, hitIndicesInTree, node->getChild(ichild), trackerHit);
108  if (seed) {
109  return seed;
110  }
111  }
112  }
113  }
114  return nullptr;
115  }
const SeedingNode< DATA > * getChild(unsigned int ichild) const
Definition: SeedingTree.h:92
const FastTrackerRecHit * hit() const
The Hit itself.
unsigned int getIndex() const
Definition: SeedingTree.h:81
const SeedingNode< TrackingLayer > * insertHit(const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
Definition: SeedFinder.h:63
unsigned int getDepth() const
Definition: SeedingTree.h:65
const TrackingLayer & getTrackingLayer() const
const SeedingNode * getParent() const
Definition: SeedingTree.h:83
unsigned int getChildrenSize() const
Definition: SeedingTree.h:90
const DATA & getData() const
Definition: SeedingTree.h:96
std::vector< std::vector< SeedFinderSelector * > > _selectorFunctionsByHits
Definition: SeedFinder.h:20

◆ iterateHits()

std::vector<unsigned int> SeedFinder::iterateHits ( unsigned int  start,
const std::vector< TrajectorySeedHitCandidate > &  trackerRecHits,
std::vector< int >  hitIndicesInTree,
bool  processSkippedHits 
) const
inline

Definition at line 117 of file SeedFinder.h.

References _seedingTree, SeedingNode< DATA >::getDepth(), SeedingNode< DATA >::getIndex(), SeedingNode< DATA >::getParent(), SeedingTree< DATA >::getRoot(), SeedingTree< DATA >::getSingleSet(), insertHit(), and SeedingTree< DATA >::numberOfRoots().

Referenced by getSeed().

120  {
121  for (unsigned int irecHit = start; irecHit < trackerRecHits.size(); ++irecHit) {
122  // only accept hits that are on one of the requested layers
123  if (_seedingTree.getSingleSet().find(trackerRecHits[irecHit].getTrackingLayer()) ==
124  _seedingTree.getSingleSet().end()) {
125  continue;
126  }
127 
128  unsigned int currentHitIndex = irecHit;
129 
130  for (unsigned int inext = currentHitIndex + 1; inext < trackerRecHits.size(); ++inext) {
131  //if multiple hits are on the same layer -> follow all possibilities by recusion
132  if (trackerRecHits[currentHitIndex].getTrackingLayer() == trackerRecHits[inext].getTrackingLayer()) {
133  if (processSkippedHits) {
134  //recusively call the method again with hit 'inext' but skip all following on the same layer though 'processSkippedHits=false'
135  std::vector<unsigned int> seedHits = iterateHits(inext, trackerRecHits, hitIndicesInTree, false);
136  if (!seedHits.empty()) {
137  return seedHits;
138  }
139  }
140  irecHit += 1;
141  } else {
142  break;
143  }
144  }
145 
146  //processSkippedHits=true
147 
148  const SeedingNode<TrackingLayer>* seedNode = nullptr;
149  for (unsigned int iroot = 0; seedNode == nullptr && iroot < _seedingTree.numberOfRoots(); ++iroot) {
150  seedNode = insertHit(trackerRecHits, hitIndicesInTree, _seedingTree.getRoot(iroot), currentHitIndex);
151  }
152  if (seedNode) {
153  std::vector<unsigned int> seedIndices(seedNode->getDepth() + 1);
154  while (seedNode) {
155  seedIndices[seedNode->getDepth()] = hitIndicesInTree[seedNode->getIndex()];
156  seedNode = seedNode->getParent();
157  }
158  return seedIndices;
159  }
160  }
161 
162  return std::vector<unsigned int>();
163  }
Definition: start.py:1
std::vector< unsigned int > iterateHits(unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
Definition: SeedFinder.h:117
const SeedingNode< DATA > * getRoot(unsigned int i) const
Definition: SeedingTree.h:159
unsigned int getIndex() const
Definition: SeedingTree.h:81
const SeedingTree< TrackingLayer > & _seedingTree
Definition: SeedFinder.h:21
const SeedingNode< TrackingLayer > * insertHit(const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
Definition: SeedFinder.h:63
unsigned int getDepth() const
Definition: SeedingTree.h:65
const SeedingNode * getParent() const
Definition: SeedingTree.h:83
const SingleSet & getSingleSet() const
Definition: SeedingTree.h:144
unsigned int numberOfRoots() const
Definition: SeedingTree.h:155

Member Data Documentation

◆ _seedingTree

const SeedingTree<TrackingLayer>& SeedFinder::_seedingTree
private

Definition at line 21 of file SeedFinder.h.

Referenced by getSeed(), and iterateHits().

◆ _selectorFunctionsByHits

std::vector<std::vector<SeedFinderSelector*> > SeedFinder::_selectorFunctionsByHits
private

Definition at line 20 of file SeedFinder.h.

Referenced by addHitSelector(), and insertHit().

◆ _trackerTopology

const TrackerTopology* SeedFinder::_trackerTopology
private

Definition at line 22 of file SeedFinder.h.

Referenced by getSeed().