CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HitExtractor.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkSeedingLayers_HitExtractor_H
2 #define RecoTracker_TkSeedingLayers_HitExtractor_H
3 
4 #include <vector>
5 #include <iterator>
8 
9 namespace edm { class Event; class EventSetup; }
10 namespace ctfseeding { class SeedingLayer; }
11 
12 namespace ctfseeding {
13 
14 class HitExtractor {
15 public:
16  typedef std::vector<TransientTrackingRecHit::ConstRecHitPointer> Hits;
17  virtual ~HitExtractor(){}
19  skipClusters=false;}
20  virtual Hits hits(const SeedingLayer & sl, const edm::Event& , const edm::EventSetup& ) const =0;
21 
22  //skip clusters
23  void useSkipClusters( const edm::InputTag & m) {
24  skipClusters=true;
26  }
29 };
30 
31 class HitConv {
32 public:
33  HitConv(const SeedingLayer &sl, const edm::EventSetup &es) : sl_(sl), es_(es) {}
34  template<typename H>
36  const TrackingRecHit* trh = &hit;
37  return sl_.hitBuilder()->build(trh); }
38 private:
39  const SeedingLayer &sl_;
41 
42 };
43 
44  template <typename DSTV, typename A, typename B>
45  inline void range2SeedingHits(DSTV const & dstv,
47  std::pair<A,B> const & sel,
48  const SeedingLayer &sl, const edm::EventSetup &es) {
49  typename DSTV::Range range = dstv.equal_range(sel.first,sel.second);
50  size_t ts = v.size();
51  for(typename DSTV::const_iterator id=range.first; id!=range.second; id++)
52  ts += std::distance((*id).begin(), (*id).end());
53  v.reserve(ts);
54  for(typename DSTV::const_iterator id=range.first; id!=range.second; id++){
55  std::transform((*id).begin(), (*id).end(), std::back_inserter(v), HitConv(sl,es));
56  }
57  }
58 
59 }
60 
61 #endif
virtual Hits hits(const SeedingLayer &sl, const edm::Event &, const edm::EventSetup &) const =0
void range2SeedingHits(DSTV const &dstv, HitExtractor::Hits &v, std::pair< A, B > const &sel, const SeedingLayer &sl, const edm::EventSetup &es)
Definition: HitExtractor.h:45
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
void useSkipClusters(const edm::InputTag &m)
Definition: HitExtractor.h:23
const SeedingLayer & sl_
Definition: HitExtractor.h:39
const TransientTrackingRecHitBuilder * hitBuilder() const
Definition: SeedingLayer.cc:85
PixelRecoRange< float > Range
HitConv(const SeedingLayer &sl, const edm::EventSetup &es)
Definition: HitExtractor.h:33
edm::InputTag theSkipClusters
Definition: HitExtractor.h:28
TransientTrackingRecHit::ConstRecHitPointer operator()(const H &hit)
Definition: HitExtractor.h:35
const edm::EventSetup & es_
Definition: HitExtractor.h:40
mathSSE::Vec4< T > v
std::vector< TransientTrackingRecHit::ConstRecHitPointer > Hits
Definition: HitExtractor.h:16