CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HitExtractorSTRP.cc
Go to the documentation of this file.
1 #include "HitExtractorSTRP.h"
4 
7 
12 
17 
18 using namespace ctfseeding;
19 using namespace std;
20 using namespace edm;
21 
23  SeedingLayer::Side & side, int idLayer)
24  : theLayer(detLayer), theSide(side), theIdLayer(idLayer),
25  hasMatchedHits(false), hasRPhiHits(false), hasStereoHits(false),
26  hasRingSelector(false), theMinRing(1), theMaxRing(0), hasSimpleRphiHitsCleaner(true)
27 { }
28 
29 void HitExtractorSTRP::useRingSelector(int minRing, int maxRing)
30 {
31  hasRingSelector=true;
32  theMinRing=minRing;
33  theMaxRing=maxRing;
34 }
35 
37 {
38  if (!hasRingSelector) return true;
39  else if ( ring >= theMinRing && ring <= theMaxRing) return true;
40  else return false;
41 }
42 
44 {
45  TrackerLayerIdAccessor accessor;
47 
48  //
49  // TIB
50  //
52  if (hasMatchedHits) {
54  ev.getByLabel( theMatchedHits, matchedHits);
55  range2SeedingHits( *matchedHits, result, accessor.stripTIBLayer(theIdLayer), sl, es);
56  }
57  if (hasRPhiHits) {
59  ev.getByLabel( theRPhiHits, rphiHits);
60  if (hasMatchedHits){
61  if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
62  range2SeedingHits( *rphiHits, result, accessor.stripTIBLayer(theIdLayer), sl, es);
63  }
64  } else {
65  range2SeedingHits( *rphiHits, result, accessor.stripTIBLayer(theIdLayer), sl, es);
66  }
67  }
68  if (hasStereoHits) {
70  ev.getByLabel( theStereoHits, stereoHits);
71  range2SeedingHits( *stereoHits, result, accessor.stripTIBLayer(theIdLayer), sl, es);
72  }
73  }
74 
75  //
76  // TID
77  //
79  if (hasMatchedHits) {
81  ev.getByLabel( theMatchedHits, matchedHits);
82  std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
83  SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
84  for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
85  int ring = TIDDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
86  for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
87  result.push_back( sl.hitBuilder()->build(hit) );
88  }
89  }
90  }
91  if (hasRPhiHits) {
93  ev.getByLabel( theRPhiHits, rphiHits);
94  std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
95  SiStripRecHit2DCollection::Range range = rphiHits->equal_range(getter.first, getter.second);
96  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
97  int ring = TIDDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
98  if ((SiStripDetId(it->detId()).partnerDetId() != 0) && hasSimpleRphiHitsCleaner) continue; // this is a brutal "cleaning". Add something smarter in the future
99  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
100  result.push_back( sl.hitBuilder()->build(hit) );
101  }
102  }
103  }
104  if (hasStereoHits) {
106  ev.getByLabel( theStereoHits, stereoHits);
107  std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
108  SiStripRecHit2DCollection::Range range = stereoHits->equal_range(getter.first, getter.second);
109  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
110  int ring = TIDDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
111  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
112  result.push_back( sl.hitBuilder()->build(hit) );
113  }
114  }
115  }
116  }
117  //
118  // TOB
119  //
121  if (hasMatchedHits) {
123  ev.getByLabel( theMatchedHits, matchedHits);
124  range2SeedingHits( *matchedHits, result, accessor.stripTOBLayer(theIdLayer), sl, es);
125  }
126  if (hasRPhiHits) {
128  ev.getByLabel( theRPhiHits, rphiHits);
129  if (hasMatchedHits){
130  if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
131  range2SeedingHits( *rphiHits, result, accessor.stripTOBLayer(theIdLayer), sl, es);
132  }
133  } else {
134  range2SeedingHits( *rphiHits, result, accessor.stripTOBLayer(theIdLayer), sl, es);
135  }
136  }
137  if (hasStereoHits) {
139  ev.getByLabel( theStereoHits, stereoHits);
140  range2SeedingHits( *stereoHits, result, accessor.stripTOBLayer(theIdLayer), sl, es);
141  }
142  }
143 
144  //
145  // TEC
146  //
148  if (hasMatchedHits) {
150  ev.getByLabel( theMatchedHits, matchedHits);
151  std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
152  SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
153  for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
154  int ring = TECDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
155  for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
156  result.push_back( sl.hitBuilder()->build(hit) );
157  }
158  }
159  }
160  if (hasRPhiHits) {
162  ev.getByLabel( theRPhiHits, rphiHits);
163  std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
164  SiStripRecHit2DCollection::Range range = rphiHits->equal_range(getter.first, getter.second);
165  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
166  int ring = TECDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
167  if ((SiStripDetId(it->detId()).partnerDetId() != 0) && hasSimpleRphiHitsCleaner) continue; // this is a brutal "cleaning". Add something smarter in the future
168  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
169  result.push_back( sl.hitBuilder()->build(hit) );
170  }
171  }
172 
173  }
174  if (hasStereoHits) {
176  ev.getByLabel( theStereoHits, stereoHits);
177  std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
178  SiStripRecHit2DCollection::Range range = stereoHits->equal_range(getter.first, getter.second);
179  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
180  int ring = TECDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
181  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
182  result.push_back( sl.hitBuilder()->build(hit) );
183  }
184  }
185  }
186  }
187 
188 
189  return result;
190 }
191 
192 
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
HitExtractorSTRP(const DetLayer *detLayer, SeedingLayer::Side &side, int idLayer)
std::pair< DetId, DetIdTIDSameDiskComparator > stripTIDDisk(int side, int disk)
bool ringRange(int ring) const
std::pair< DetId, DetIdTIBSameLayerComparator > stripTIBLayer(int layer)
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
std::pair< const_iterator, const_iterator > Range
void range2SeedingHits(DSTV const &dstv, HitExtractor::Hits &v, std::pair< A, B > const &sel, const SeedingLayer &sl, const edm::EventSetup &es)
Definition: HitExtractor.h:34
std::pair< DetId, DetIdTOBSameLayerComparator > stripTOBLayer(int layer)
virtual HitExtractor::Hits hits(const SeedingLayer &sl, const edm::Event &, const edm::EventSetup &) const
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
tuple result
Definition: query.py:137
std::pair< DetId, DetIdTECSameDiskComparator > stripTECDisk(int side, int disk)
#define end
Definition: vmac.h:38
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
void useRingSelector(int minRing, int maxRing)
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
const TransientTrackingRecHitBuilder * hitBuilder() const
Definition: SeedingLayer.cc:85
std::vector< TransientTrackingRecHit::ConstRecHitPointer > Hits
Definition: HitExtractor.h:16