CMS 3D CMS Logo

DTCombinatorialPatternReco.h
Go to the documentation of this file.
1 #ifndef DTSegment_DTCombinatorialPatternReco_h
2 #define DTSegment_DTCombinatorialPatternReco_h
3 
13 /* Base Class Headers */
16 
17 /* Collaborating Class Declarations */
18 namespace edm {
19  class ParameterSet;
20  class EventSetup;
21  // class ESHandle;
22 } // namespace edm
23 class DTSegmentUpdator;
24 class DTSegmentCleaner;
25 class DTHitPairForFit;
26 class DTSegmentCand;
27 class MuonGeometryRecord;
28 
29 /* C++ Headers */
30 #include <deque>
31 #include <functional>
32 #include <unordered_set>
33 #include <utility>
34 #include <vector>
35 
40 
41 /* ====================================================================== */
42 
43 /* Class DTCombinatorialPatternReco Interface */
44 
46 public:
49 
51  ~DTCombinatorialPatternReco() override;
52 
53  /* Operations */
54 
57  const std::vector<DTRecHit1DPair>& hits) override;
58 
60  std::string algoName() const override { return theAlgoName; }
61 
64  void setES(const edm::EventSetup& setup) override;
65 
66 protected:
67 private:
69 
70  // typedef std::pair<DTHitPairForFit*, DTEnums::DTCellSide> AssPoint;
71 
72  // create the DTHitPairForFit from the pairs for easy use
73  std::vector<std::shared_ptr<DTHitPairForFit>> initHits(const DTSuperLayer* sl,
74  const std::vector<DTRecHit1DPair>& hits);
75 
76  // search for candidate, starting from pairs of hits in different layers
77  std::vector<DTSegmentCand*> buildSegments(const DTSuperLayer* sl,
78  const std::vector<std::shared_ptr<DTHitPairForFit>>& hits);
79 
80  // find all the hits compatible with the candidate
81  std::vector<DTSegmentCand::AssPoint> findCompatibleHits(const LocalPoint& pos,
82  const LocalVector& dir,
83  const std::vector<std::shared_ptr<DTHitPairForFit>>& hits);
84 
85  // build segments from hits collection
86  DTSegmentCand* buildBestSegment(std::vector<DTSegmentCand::AssPoint>& assHits, const DTSuperLayer* sl);
87 
88  bool checkDoubleCandidates(std::vector<DTSegmentCand*>& segs, DTSegmentCand* seg);
89 
92  void buildPointsCollection(std::vector<DTSegmentCand::AssPoint>& points,
93  std::deque<std::shared_ptr<DTHitPairForFit>>& pointsNoLR,
94  std::vector<DTSegmentCand*>& candidates,
95  const DTSuperLayer* sl);
96 
97 private:
99  unsigned int theMaxAllowedHits;
102  bool debug;
103  bool usePairs;
104  DTSegmentUpdator* theUpdator; // the updator and fitter
105  DTSegmentCleaner* theCleaner; // the cleaner
106 
109 
110 public:
111  // The type must be public, as otherwise the global 'hash_value' function can't locate it
112  class TriedPattern {
113  public:
114  typedef std::vector<short unsigned int> values;
115 
116  // empty costructor
117  TriedPattern() : hash_(1) { values_.reserve(8); }
118 
119  // equality operator
120  bool operator==(const TriedPattern& other) const {
121  return (hash_ == other.hash_) && // cheap
122  (values_ == other.values_); // expensive last resort
123  }
124 
125  //Same implementation as boost::hash_combine
126  template <class T>
127  inline void hash_combine(std::size_t& seed, const T& v) {
128  std::hash<T> hasher;
129  seed ^= hasher(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
130  }
131 
133  void push_back(short unsigned int i) {
134  hash_combine(hash_, i);
135  values_.push_back(i);
136  }
139  size_t hash() const { return hash_; }
140 
141  // some extra methods to look like a std::vector
142  typedef values::const_iterator const_iterator;
143  const_iterator begin() const { return values_.begin(); }
144  const_iterator end() const { return values_.end(); }
145  values::size_type size() const { return values_.size(); }
146 
147  //Custom hash functor for std::unordered_set
148  class HashFunction {
149  public:
150  size_t operator()(const TriedPattern& p) const { return p.hash(); }
151  };
152 
153  private:
155  size_t hash_;
156  };
157  typedef std::unordered_set<TriedPattern, TriedPattern::HashFunction> TriedPatterns;
158 
159 private:
161 };
162 
163 inline std::size_t hash_value(const DTCombinatorialPatternReco::TriedPattern& t) { return t.hash(); }
164 
165 #endif // DTSegment_DTCombinatorialPatternReco_h
~DTCombinatorialPatternReco() override
Destructor.
std::vector< std::shared_ptr< DTHitPairForFit > > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits) override
this function is called in the producer
void hash_combine(std::size_t &seed, const T &v)
uint16_t size_type
edm::ESHandle< DTGeometry > theDTGeometry
void buildPointsCollection(std::vector< DTSegmentCand::AssPoint > &points, std::deque< std::shared_ptr< DTHitPairForFit >> &pointsNoLR, std::vector< DTSegmentCand *> &candidates, const DTSuperLayer *sl)
DTCombinatorialPatternReco(const edm::ParameterSet &pset, edm::ConsumesCollector)
Constructor.
std::size_t hash_value(const DTCombinatorialPatternReco::TriedPattern &t)
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > theDTGeometryToken
void push_back(short unsigned int i)
push back value, and update the hash
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
bool checkDoubleCandidates(std::vector< DTSegmentCand *> &segs, DTSegmentCand *seg)
DTSegmentCand * buildBestSegment(std::vector< DTSegmentCand::AssPoint > &assHits, const DTSuperLayer *sl)
std::string algoName() const override
return the algo name
HLT enums.
std::vector< DTSegmentCand::AssPoint > findCompatibleHits(const LocalPoint &pos, const LocalVector &dir, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
std::unordered_set< TriedPattern, TriedPattern::HashFunction > TriedPatterns
long double T
void setES(const edm::EventSetup &setup) override
bool operator==(const TriedPattern &other) const