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 
10 
12 
17 
20 
21 #include<tuple>
22 
23 #include<iostream>
24 
25 using namespace ctfseeding;
26 using namespace std;
27 using namespace edm;
28 
30  theLayerSubDet(subdet), theSide(side), theIdLayer(idLayer),
31  minAbsZ(0), theMinRing(1), theMaxRing(0),
32  hasMatchedHits(false), hasRPhiHits(false), hasStereoHits(false),
33  hasRingSelector(false), hasSimpleRphiHitsCleaner(true)
34 {}
35 
38 }
39 
41 {
42  hasRingSelector=true;
45 }
46 
48 {
49  if (!hasRingSelector) return true;
50  return (ring >= theMinRing) & (ring <= theMaxRing);
51 }
52 
55  return stripClusterMask->mask(clus.key());
56 }
57 
58 
59 
60 std::pair<bool,ProjectedSiStripRecHit2D *>
62  TkHitRef matched,
64  const SiStripMatchedRecHit2D & hit = (SiStripMatchedRecHit2D const&)(matched);
65 
66  assert(dynamic_cast<SiStripMatchedRecHit2D const*>(&matched));
67 
68  ProjectedSiStripRecHit2D * replaceMe = nullptr;
69  bool rejectSt = skipThis(hit.stereoClusterRef(), stripClusterMask);
70  bool rejectMono = skipThis(hit.monoClusterRef(), stripClusterMask);
71 
72  if ((!rejectSt)&(!rejectMono)){
73  // keepit
74  return std::make_pair(false,replaceMe);
75  }
76 
77  if (failProjection || (rejectSt&rejectMono) ){
78  //only skip if both hits are done
79  return std::make_pair(true,replaceMe);
80  }
81 
82  // replace with one
83 
84  auto cloner = ttrhBuilder.cloner();
85  replaceMe = cloner.project(hit, rejectSt, TrajectoryStateOnSurface());
86  if (rejectSt)
87  LogDebug("HitExtractorSTRP")<<"a matched hit is partially masked, and the mono hit got projected onto: "<<replaceMe->geographicalId().rawId()<<" key: "<<hit.monoClusterRef().key();
88  else
89  LogDebug("HitExtractorSTRP")<<"a matched hit is partially masked, and the stereo hit got projected onto: "<<replaceMe->geographicalId().rawId()<<" key: "<<hit.stereoClusterRef().key();
90 
91  return std::make_pair(true,replaceMe);
92 }
93 
94 
96  const edm::Event& ev, HitExtractor::Hits & hits,
97  bool matched,
98  unsigned int cleanFrom) const{
99  LogDebug("HitExtractorPIX")<<"getting: "<<hits.size()<<" in input.";
100  edm::Handle<SkipClustersCollection> stripClusterMask;
101  ev.getByToken(theSkipClusters,stripClusterMask);
102  unsigned int skipped=0;
103  unsigned int projected=0;
104  for (unsigned int iH=cleanFrom;iH<hits.size();++iH){
105  assert(hits[iH]->isValid());
106  if (matched) {
107  bool replace; ProjectedSiStripRecHit2D * replaceMe; std::tie(replace,replaceMe) = skipThis(ttrhBuilder, *hits[iH],stripClusterMask);
108  if (replace) {
109  if (!replaceMe) {
110  LogDebug("HitExtractorSTRP")<<"skipping a matched hit on :"<<hits[iH]->geographicalId().rawId();
111  skipped++;
112  } else projected++;
113  hits[iH].reset(replaceMe);
114  if (replaceMe==nullptr) assert(hits[iH].empty());
115  else assert(hits[iH].isOwn());
116  }
117  }
118  else if (skipThis(hits[iH]->firstClusterRef(),stripClusterMask)){
119  LogDebug("HitExtractorSTRP")<<"skipping a hit on :"<<hits[iH]->geographicalId().rawId()<<" key: ";
120  skipped++;
121  hits[iH].reset();
122  }
123  }
124  // remove empty elements...
125  auto last = std::remove_if(hits.begin()+cleanFrom,hits.end(),[]( HitPointer const & p) {return p.empty();});
126  hits.resize(last-hits.begin());
127 
128  // std::cout << "HitExtractorSTRP " <<"skipped :"<<skipped<<" strip rechits because of clusters and projected: "<<projected << std::endl;
129  LogDebug("HitExtractorSTRP")<<"skipped :"<<skipped<<" strip rechits because of clusters and projected: "<<projected;
130 }
131 
133 {
135  TrackerLayerIdAccessor accessor;
136  unsigned int cleanFrom=0;
137  //
138  // TIB
139  //
141  if (hasMatchedHits) {
143  ev.getByToken( theMatchedHits, matchedHits);
144  if (skipClusters) cleanFrom=result.size();
145  range2SeedingHits( *matchedHits, result, accessor.stripTIBLayer(theIdLayer));
146  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,true,cleanFrom);
147  }
148  if (hasRPhiHits) {
150  ev.getByToken( theRPhiHits, rphiHits);
151  if (hasMatchedHits){
152  if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
153  if (skipClusters) cleanFrom=result.size();
154  range2SeedingHits( *rphiHits, result, accessor.stripTIBLayer(theIdLayer));
155  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
156  }
157  } else {
158  if (skipClusters) cleanFrom=result.size();
159  range2SeedingHits( *rphiHits, result, accessor.stripTIBLayer(theIdLayer));
160  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
161  }
162  }
163  if (hasStereoHits) {
165  ev.getByToken( theStereoHits, stereoHits);
166  if (skipClusters) cleanFrom=result.size();
167  range2SeedingHits( *stereoHits, result, accessor.stripTIBLayer(theIdLayer));
168  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
169  }
170  }
171 
172  //
173  // TID
174  //
176  if (hasMatchedHits) {
178  ev.getByToken( theMatchedHits, matchedHits);
179  if (skipClusters) cleanFrom=result.size();
180  std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
181  SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
182  for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
183  int ring = TIDDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
184  for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
185  result.emplace_back(*hit);
186  }
187  }
188  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,true,cleanFrom);
189  }
190  if (hasRPhiHits) {
192  ev.getByToken( theRPhiHits, rphiHits);
193  if (skipClusters) cleanFrom=result.size();
194  std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
195  SiStripRecHit2DCollection::Range range = rphiHits->equal_range(getter.first, getter.second);
196  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
197  int ring = TIDDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
198  if ((SiStripDetId(it->detId()).partnerDetId() != 0) && hasSimpleRphiHitsCleaner) continue; // this is a brutal "cleaning". Add something smarter in the future
199  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
200  result.emplace_back(*hit);
201  }
202  }
203  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
204  }
205  if (hasStereoHits) {
207  ev.getByToken( theStereoHits, stereoHits);
208  if (skipClusters) cleanFrom=result.size();
209  std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
210  SiStripRecHit2DCollection::Range range = stereoHits->equal_range(getter.first, getter.second);
211  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
212  int ring = TIDDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
213  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
214  result.emplace_back(*hit);
215  }
216  }
217  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
218  }
219  }
220  //
221  // TOB
222  //
224  if (hasMatchedHits) {
226  ev.getByToken( theMatchedHits, matchedHits);
227  if (skipClusters) cleanFrom=result.size();
228  if (minAbsZ>0.) {
229  std::pair<DetId,DetIdTOBSameLayerComparator> getter = accessor.stripTOBLayer(theIdLayer);
230  SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
231  for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
232  for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
233  if (fabs(hit->globalPosition().z())>=minAbsZ) result.emplace_back(*hit);
234  }
235  }
236  } else {
237  range2SeedingHits( *matchedHits, result, accessor.stripTOBLayer(theIdLayer));
238  }
239  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,true,cleanFrom);
240  }
241  if (hasRPhiHits) {
243  ev.getByToken( theRPhiHits, rphiHits);
244  if (hasMatchedHits){
245  if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
246  if (skipClusters) cleanFrom=result.size();
247  range2SeedingHits( *rphiHits, result, accessor.stripTOBLayer(theIdLayer));
248  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
249  }
250  } else {
251  if (skipClusters) cleanFrom=result.size();
252  range2SeedingHits( *rphiHits, result, accessor.stripTOBLayer(theIdLayer));
253  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
254  }
255  }
256  if (hasStereoHits) {
258  ev.getByToken( theStereoHits, stereoHits);
259  if (skipClusters) cleanFrom=result.size();
260  range2SeedingHits( *stereoHits, result, accessor.stripTOBLayer(theIdLayer));
261  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
262  }
263  }
264 
265  //
266  // TEC
267  //
269  if (hasMatchedHits) {
271  ev.getByToken( theMatchedHits, matchedHits);
272  if (skipClusters) cleanFrom=result.size();
273  std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
274  SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
275  for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
276  int ring = TECDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
277  for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
278  result.emplace_back(*hit);
279  }
280  }
281  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,true,cleanFrom);
282  }
283  if (hasRPhiHits) {
285  ev.getByToken( theRPhiHits, rphiHits);
286  if (skipClusters) cleanFrom=result.size();
287  std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
288  SiStripRecHit2DCollection::Range range = rphiHits->equal_range(getter.first, getter.second);
289  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
290  int ring = TECDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
291  if ((SiStripDetId(it->detId()).partnerDetId() != 0) && hasSimpleRphiHitsCleaner) continue; // this is a brutal "cleaning". Add something smarter in the future
292  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
293  result.emplace_back(*hit);
294  }
295  }
296  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
297 
298  }
299  if (hasStereoHits) {
301  ev.getByToken( theStereoHits, stereoHits);
302  if (skipClusters) cleanFrom=result.size();
303  std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
304  SiStripRecHit2DCollection::Range range = stereoHits->equal_range(getter.first, getter.second);
305  for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
306  int ring = TECDetId( it->detId() ).ring(); if (!ringRange(ring)) continue;
307  for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
308  result.emplace_back(*hit);
309  }
310  }
311  if (skipClusters) cleanedOfClusters(ttrhBuilder, ev,result,false,cleanFrom);
312  }
313  }
314  /* done in each skipCluster...
315  // std::cout << "HitExtractorSTRP before cleanup "<<" giving: "<<result.size()<< std::endl;
316  // remove empty elements...
317  auto last = std::remove_if(result.begin(),result.end(),[]( HitPointer const & p) {return p.empty();});
318  result.resize(last-result.begin());
319  */
320  LogDebug("HitExtractorSTRP")<<" giving: "<<result.size()<<" out";
321  // std::cout << "HitExtractorSTRP "<<" giving: "<<result.size()<< std::endl;
322  return result;
323 }
324 
325 
#define LogDebug(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
HitExtractorSTRP(GeomDetEnumerators::SubDetector subdet, SeedingLayer::Side &side, int idLayer)
void range2SeedingHits(DSTV const &dstv, HitExtractor::Hits &v, std::pair< A, B > const &sel)
Definition: HitExtractor.h:45
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
std::pair< DetId, DetIdTIDSameDiskComparator > stripTIDDisk(int side, int disk)
std::vector< HitPointer > Hits
Definition: HitExtractor.h:25
bool ringRange(int ring) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
BaseTrackerRecHit const & TkHitRef
Definition: HitExtractor.h:23
OmniClusterRef const & stereoClusterRef() const
std::pair< DetId, DetIdTIBSameLayerComparator > stripTIBLayer(int layer)
bool ev
std::pair< const_iterator, const_iterator > Range
edm::EDGetTokenT< SiStripMatchedRecHit2DCollection > theMatchedHits
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
std::pair< bool, ProjectedSiStripRecHit2D * > skipThis(const TkTransientTrackingRecHitBuilder &ttrhBuilder, TkHitRef matched, edm::Handle< edm::ContainerMask< edmNew::DetSetVector< SiStripCluster > > > &stripClusterMask) const
std::pair< DetId, DetIdTOBSameLayerComparator > stripTOBLayer(int layer)
const GeomDetEnumerators::SubDetector theLayerSubDet
edm::EDGetTokenT< SiStripRecHit2DCollection > theStereoHits
tuple result
Definition: query.py:137
std::pair< DetId, DetIdTECSameDiskComparator > stripTECDisk(int side, int disk)
void cleanedOfClusters(const TkTransientTrackingRecHitBuilder &ttrhBuilder, const edm::Event &ev, HitExtractor::Hits &hits, bool matched, unsigned int cleanFrom=0) const
#define end
Definition: vmac.h:37
OmniClusterRef const & monoClusterRef() const
virtual HitExtractor::Hits hits(const TkTransientTrackingRecHitBuilder &ttrhBuilder, const edm::Event &, const edm::EventSetup &) const override
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
edm::EDGetTokenT< SiStripRecHit2DCollection > theRPhiHits
DetId geographicalId() const
volatile std::atomic< bool > shutdown_flag false
unsigned int key() const
void useSkipClusters_(const edm::InputTag &m, edm::ConsumesCollector &iC) override
edm::EDGetTokenT< SkipClustersCollection > theSkipClusters