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 */
15 #include <boost/unordered_set.hpp>
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 
28 /* C++ Headers */
29 #include <vector>
30 #include <deque>
31 #include <utility>
32 
36 
37 /* ====================================================================== */
38 
39 /* Class DTCombinatorialPatternReco Interface */
40 
42 public:
45 
47  ~DTCombinatorialPatternReco() override;
48 
49  /* Operations */
50 
53  const std::vector<DTRecHit1DPair>& hits) override;
54 
56  std::string algoName() const override { return theAlgoName; }
57 
60  void setES(const edm::EventSetup& setup) override;
61 
62 protected:
63 private:
65 
66  // typedef std::pair<DTHitPairForFit*, DTEnums::DTCellSide> AssPoint;
67 
68  // create the DTHitPairForFit from the pairs for easy use
69  std::vector<std::shared_ptr<DTHitPairForFit>> initHits(const DTSuperLayer* sl,
70  const std::vector<DTRecHit1DPair>& hits);
71 
72  // search for candidate, starting from pairs of hits in different layers
73  std::vector<DTSegmentCand*> buildSegments(const DTSuperLayer* sl,
74  const std::vector<std::shared_ptr<DTHitPairForFit>>& hits);
75 
76  // find all the hits compatible with the candidate
77  std::vector<DTSegmentCand::AssPoint> findCompatibleHits(const LocalPoint& pos,
78  const LocalVector& dir,
79  const std::vector<std::shared_ptr<DTHitPairForFit>>& hits);
80 
81  // build segments from hits collection
82  DTSegmentCand* buildBestSegment(std::vector<DTSegmentCand::AssPoint>& assHits, const DTSuperLayer* sl);
83 
84  bool checkDoubleCandidates(std::vector<DTSegmentCand*>& segs, DTSegmentCand* seg);
85 
88  void buildPointsCollection(std::vector<DTSegmentCand::AssPoint>& points,
89  std::deque<std::shared_ptr<DTHitPairForFit>>& pointsNoLR,
90  std::vector<DTSegmentCand*>& candidates,
91  const DTSuperLayer* sl);
92 
93 private:
95  unsigned int theMaxAllowedHits;
98  bool debug;
99  bool usePairs;
100  DTSegmentUpdator* theUpdator; // the updator and fitter
101  DTSegmentCleaner* theCleaner; // the cleaner
102 
104 
105 public:
106  // The type must be public, as otherwise the global 'hash_value' function can't locate it
107  class TriedPattern {
108  public:
109  typedef std::vector<short unsigned int> values;
110 
111  // empty costructor
112  TriedPattern() : hash_(1) { values_.reserve(8); }
113 
114  // equality operator
115  bool operator==(const TriedPattern& other) const {
116  return (hash_ == other.hash_) && // cheap
117  (values_ == other.values_); // expensive last resort
118  }
119 
121  void push_back(short unsigned int i) {
122  boost::hash_combine(hash_, i);
123  values_.push_back(i);
124  }
127  size_t hash() const { return hash_; }
128 
129  // some extra methods to look like a std::vector
130  typedef values::const_iterator const_iterator;
131  const_iterator begin() const { return values_.begin(); }
132  const_iterator end() const { return values_.end(); }
133  values::size_type size() const { return values_.size(); }
134 
135  private:
137  size_t hash_;
138  };
139  typedef boost::unordered_set<TriedPattern> TriedPatterns;
140 
141 private:
143 };
144 
145 inline std::size_t hash_value(const DTCombinatorialPatternReco::TriedPattern& t) { return t.hash(); }
146 #endif // DTSegment_DTCombinatorialPatternReco_h
Vector3DBase< float, LocalTag >
HLT_2018_cff.points
points
Definition: HLT_2018_cff.py:20125
DTCombinatorialPatternReco::theCleaner
DTSegmentCleaner * theCleaner
Definition: DTCombinatorialPatternReco.h:101
DTRecSegment2DBaseAlgo.h
mps_fire.i
i
Definition: mps_fire.py:355
DTCombinatorialPatternReco::theDTGeometry
edm::ESHandle< DTGeometry > theDTGeometry
Definition: DTCombinatorialPatternReco.h:103
hash_value
std::size_t hash_value(const DTCombinatorialPatternReco::TriedPattern &t)
Definition: DTCombinatorialPatternReco.h:145
DTCombinatorialPatternReco::algoName
std::string algoName() const override
return the algo name
Definition: DTCombinatorialPatternReco.h:56
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
ESHandle.h
DTCombinatorialPatternReco::TriedPatterns
boost::unordered_set< TriedPattern > TriedPatterns
Definition: DTCombinatorialPatternReco.h:139
DTCombinatorialPatternReco::TriedPattern::hash_
size_t hash_
Definition: DTCombinatorialPatternReco.h:137
edm
HLT enums.
Definition: AlignableModifier.h:19
pos
Definition: PixelAliasList.h:18
DTCombinatorialPatternReco::usePairs
bool usePairs
Definition: DTCombinatorialPatternReco.h:99
DTCombinatorialPatternReco
Definition: DTCombinatorialPatternReco.h:41
DTCombinatorialPatternReco::TriedPattern::end
const_iterator end() const
Definition: DTCombinatorialPatternReco.h:132
DTSuperLayer
Definition: DTSuperLayer.h:24
DTCombinatorialPatternReco::checkDoubleCandidates
bool checkDoubleCandidates(std::vector< DTSegmentCand * > &segs, DTSegmentCand *seg)
Definition: DTCombinatorialPatternReco.cc:408
DTCombinatorialPatternReco::buildBestSegment
DTSegmentCand * buildBestSegment(std::vector< DTSegmentCand::AssPoint > &assHits, const DTSuperLayer *sl)
Definition: DTCombinatorialPatternReco.cc:289
singleTopDQM_cfi.setup
setup
Definition: singleTopDQM_cfi.py:37
DTCombinatorialPatternReco::initHits
std::vector< std::shared_ptr< DTHitPairForFit > > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
Definition: DTCombinatorialPatternReco.cc:83
DTSegmentUpdator
Definition: DTSegmentUpdator.h:43
DTSegmentCand.h
DTSegmentCleaner
Definition: DTSegmentCleaner.h:31
trigger::size_type
uint16_t size_type
Definition: TriggerTypeDefs.h:18
DTCombinatorialPatternReco::TriedPattern
Definition: DTCombinatorialPatternReco.h:107
edm::ESHandle< DTGeometry >
trackingPlots.other
other
Definition: trackingPlots.py:1465
DTSegmentCand
Definition: DTSegmentCand.h:34
DTCombinatorialPatternReco::buildPointsCollection
void buildPointsCollection(std::vector< DTSegmentCand::AssPoint > &points, std::deque< std::shared_ptr< DTHitPairForFit >> &pointsNoLR, std::vector< DTSegmentCand * > &candidates, const DTSuperLayer *sl)
Definition: DTCombinatorialPatternReco.cc:353
DTCombinatorialPatternReco::theMaxAllowedHits
unsigned int theMaxAllowedHits
Definition: DTCombinatorialPatternReco.h:95
Point3DBase< float, LocalTag >
OrderedSet.t
t
Definition: OrderedSet.py:90
DTCombinatorialPatternReco::TriedPattern::operator==
bool operator==(const TriedPattern &other) const
Definition: DTCombinatorialPatternReco.h:115
DTCombinatorialPatternReco::theAlgoName
std::string theAlgoName
Definition: DTCombinatorialPatternReco.h:94
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
DTGeometry.h
DTCombinatorialPatternReco::DTCombinatorialPatternReco
DTCombinatorialPatternReco(const edm::ParameterSet &pset)
Constructor.
Definition: DTCombinatorialPatternReco.cc:33
DTCombinatorialPatternReco4D
Definition: DTCombinatorialPatternReco4D.h:40
edm::ParameterSet
Definition: ParameterSet.h:36
ParameterSet
Definition: Functions.h:16
DTCombinatorialPatternReco::TriedPattern::TriedPattern
TriedPattern()
Definition: DTCombinatorialPatternReco.h:112
DTCombinatorialPatternReco::TriedPattern::size
values::size_type size() const
Definition: DTCombinatorialPatternReco.h:133
DTCombinatorialPatternReco::TriedPattern::begin
const_iterator begin() const
Definition: DTCombinatorialPatternReco.h:131
DTCombinatorialPatternReco::theTriedPattern
TriedPatterns theTriedPattern
Definition: DTCombinatorialPatternReco.h:142
edm::EventSetup
Definition: EventSetup.h:57
DTCombinatorialPatternReco::theUpdator
DTSegmentUpdator * theUpdator
Definition: DTCombinatorialPatternReco.h:100
DTCombinatorialPatternReco::~DTCombinatorialPatternReco
~DTCombinatorialPatternReco() override
Destructor.
Definition: DTCombinatorialPatternReco.cc:46
DTCombinatorialPatternReco::TriedPattern::values_
values values_
Definition: DTCombinatorialPatternReco.h:136
DTCombinatorialPatternReco::TriedPattern::values
std::vector< short unsigned int > values
Definition: DTCombinatorialPatternReco.h:109
DTCombinatorialPatternReco::setES
void setES(const edm::EventSetup &setup) override
Definition: DTCombinatorialPatternReco.cc:77
HLT_2018_cff.candidates
candidates
Definition: HLT_2018_cff.py:53513
DTCombinatorialPatternReco::findCompatibleHits
std::vector< DTSegmentCand::AssPoint > findCompatibleHits(const LocalPoint &pos, const LocalVector &dir, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
Definition: DTCombinatorialPatternReco.cc:216
DTCombinatorialPatternReco::TriedPattern::push_back
void push_back(short unsigned int i)
push back value, and update the hash
Definition: DTCombinatorialPatternReco.h:121
DTHitPairForFit
Definition: DTHitPairForFit.h:31
DTCombinatorialPatternReco::reconstruct
edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits) override
this function is called in the producer
Definition: DTCombinatorialPatternReco.cc:52
DTCombinatorialPatternReco::theAlphaMaxTheta
double theAlphaMaxTheta
Definition: DTCombinatorialPatternReco.h:96
DTCombinatorialPatternReco::buildSegments
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
Definition: DTCombinatorialPatternReco.cc:92
EventSetup
DTCombinatorialPatternReco::theAlphaMaxPhi
double theAlphaMaxPhi
Definition: DTCombinatorialPatternReco.h:97
DTCombinatorialPatternReco::debug
bool debug
Definition: DTCombinatorialPatternReco.h:98
DTCombinatorialPatternReco::TriedPattern::hash
size_t hash() const
Definition: DTCombinatorialPatternReco.h:127
DTCombinatorialPatternReco::TriedPattern::const_iterator
values::const_iterator const_iterator
Definition: DTCombinatorialPatternReco.h:130
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
DTRecSegment2DBaseAlgo
Definition: DTRecSegment2DBaseAlgo.h:32
edm::OwnVector
Definition: OwnVector.h:24
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23