CMS 3D CMS Logo

SiStripRecHitConverterAlgorithm.cc
Go to the documentation of this file.
4 
7 
10 
15 
16 
18  useQuality(conf.getParameter<bool>("useSiStripQuality")),
19  maskBad128StripBlocks( conf.existsAs<bool>("MaskBadAPVFibers") && conf.getParameter<bool>("MaskBadAPVFibers")),
20  tracker_cache_id(0),
21  cpe_cache_id(0),
22  quality_cache_id(0),
23  cpeTag(conf.getParameter<edm::ESInputTag>("StripCPE")),
24  matcherTag(conf.getParameter<edm::ESInputTag>("Matcher")),
25  qualityTag(conf.getParameter<edm::ESInputTag>("siStripQualityLabel"))
26 {}
27 
30 {
31  uint32_t tk_cache_id = es.get<TrackerDigiGeometryRecord>().cacheIdentifier();
32  uint32_t c_cache_id = es.get<TkStripCPERecord>().cacheIdentifier();
33  uint32_t q_cache_id = es.get<SiStripQualityRcd>().cacheIdentifier();
34 
35  if(tk_cache_id != tracker_cache_id) {
37  tracker_cache_id = tk_cache_id;
38  }
39  if(c_cache_id != cpe_cache_id) {
42  cpe_cache_id = c_cache_id;
43  }
44  if( useQuality && q_cache_id!=quality_cache_id) {
46  quality_cache_id = q_cache_id;
47  }
48 }
49 
52 { run(input, output, LocalVector(0.,0.,0.)); }
53 
54 
55 
58 {
59 
60  for (auto const & DS : *inputhandle) {
61  auto id = DS.id();
62  if(!useModule(id)) continue;
63 
64  Collector collector = StripSubdetector(id).stereo()
65  ? Collector(*output.stereo, id)
66  : Collector(*output.rphi, id);
67 
68  bool bad128StripBlocks[6]; fillBad128StripBlocks( id, bad128StripBlocks);
69 
70  GeomDetUnit const & du = *(tracker->idToDetUnit(id));
71  for(auto const & cluster : DS) {
72 
73  if(isMasked(cluster,bad128StripBlocks)) continue;
74 
76  collector.push_back(SiStripRecHit2D( parameters.first, parameters.second, du,
77  DS.makeRefTo(inputhandle, &cluster)
78  ));
79  }
80 
81  if (collector.empty()) collector.abort();
82  }
83  match(output,trackdirection);
84 }
85 
86 
87 namespace {
88 
89  struct CollectorHelper {
90  size_t nmatch;
91 
92 
93  typedef edm::OwnVector<SiStripMatchedRecHit2D> CollectorMatched;
95 
96  Collector & m_collector;
97  CollectorMatched & m_collectorMatched;
98  SiStripRecHit2DCollection::FastFiller & m_fillerRphiUnm;
99  std::vector<SiStripRecHit2D::ClusterRef::key_type> & m_matchedSteroClusters;
100 
101  static inline SiStripRecHit2D const & stereoHit(edmNew::DetSet<SiStripRecHit2D>::const_iterator iter) {
102  return *iter;
103  }
104 
105  static inline SiStripRecHit2D const & monoHit(edmNew::DetSet<SiStripRecHit2D>::const_iterator iter) {
106  return *iter;
107  }
108 
109  struct Add {
110  Add(CollectorHelper& ih) : h(ih){}
111  CollectorHelper& h;
112  void operator()(SiStripMatchedRecHit2D const & rh) { h.m_collectorMatched.push_back(rh);}
113  };
114 
115  CollectorHelper & collector() {
116  return *this;
117  }
118 
119  void operator()(SiStripMatchedRecHit2D const & rh) {m_collectorMatched.push_back(rh);}
120 
121 
122  CollectorHelper(
123  Collector & i_collector,
124  CollectorMatched & i_collectorMatched,
125  SiStripRecHit2DCollection::FastFiller & i_fillerRphiUnm,
126  std::vector<SiStripRecHit2D::ClusterRef::key_type> & i_matchedSteroClusters
127  ) : nmatch(0),
128  m_collector(i_collector),
129  m_collectorMatched(i_collectorMatched),
130  m_fillerRphiUnm(i_fillerRphiUnm),
131  m_matchedSteroClusters(i_matchedSteroClusters) {}
132 
134  if (!m_collectorMatched.empty()){
135  nmatch+=m_collectorMatched.size();
136  for (edm::OwnVector<SiStripMatchedRecHit2D>::const_iterator itm = m_collectorMatched.begin(),
137  edm = m_collectorMatched.end();
138  itm != edm;
139  ++itm) {
140  m_collector.push_back(*itm);
141  // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
142  m_matchedSteroClusters.push_back(itm->stereoClusterRef().key());
143  }
144  m_collectorMatched.clear();
145  } else {
146  // store a copy of this rphi hit as an unmatched rphi hit
147  m_fillerRphiUnm.push_back(*it);
148  }
149  }
150  };
151 }
152 
153 
155 match(products& output, LocalVector trackdirection) const
156 {
157  int nmatch=0;
158  edm::OwnVector<SiStripMatchedRecHit2D> collectorMatched; // gp/FIXME: avoid this
159 
160  // Remember the ends of the collections, as we will use them a lot
161  SiStripRecHit2DCollection::const_iterator edStereoDet = output.stereo->end();
162  SiStripRecHit2DCollection::const_iterator edRPhiDet = output.rphi->end();
163 
164  // two work vectors for bookeeping clusters used by the stereo part of the matched hits
165  std::vector<SiStripRecHit2D::ClusterRef::key_type> matchedSteroClusters;
166 
167  for (SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->begin(); itRPhiDet != edRPhiDet; ++itRPhiDet) {
168  edmNew::DetSet<SiStripRecHit2D> rphiHits = *itRPhiDet;
169  StripSubdetector specDetId(rphiHits.detId());
170  uint32_t partnerId = specDetId.partnerDetId();
171 
172  // if not part of a glued pair
173  if (partnerId == 0) {
174  // I must copy these as unmatched
175  if (!rphiHits.empty()) {
177  filler.resize(rphiHits.size());
178  std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
179  }
180  continue;
181  }
182 
183  SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->find(partnerId);
184 
185  // if the partner is not found (which probably can happen if it's empty)
186  if (itStereoDet == edStereoDet) {
187  // I must copy these as unmatched
188  if (!rphiHits.empty()) {
190  filler.resize(rphiHits.size());
191  std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
192  }
193  continue;
194  }
195 
196  edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
197 
198 
199  // Get ready for making glued hits
200  const GluedGeomDet* gluedDet = (const GluedGeomDet*)tracker->idToDet(DetId(specDetId.glued()));
202  Collector collector(*output.matched, specDetId.glued());
203 
204  // Prepare also the list for unmatched rphi hits
205  SiStripRecHit2DCollection::FastFiller fillerRphiUnm(*output.rphiUnmatched, rphiHits.detId());
206 
207  // a list of clusters used by the matched part of the stereo hits in this detector
208  matchedSteroClusters.clear(); // at the beginning, empty
209 
210 #ifdef DOUBLE_MATCH
211  CollectorHelper chelper(collector, collectorMatched,
212  fillerRphiUnm,
213  matchedSteroClusters
214  );
215  matcher->doubleMatch(rphiHits.begin(), rphiHits.end(),
216  stereoHits.begin(),stereoHits.end(),gluedDet,trackdirection,chelper);
217  nmatch+=chelper.nmatch;
218 #else
219  // Make simple collection of this (gp:FIXME: why do we need it?)
221  // gp:FIXME: use std::transform
222  stereoSimpleHits.reserve(stereoHits.size());
223  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
224  stereoSimpleHits.push_back(&*it);
225  }
226 
227  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = rphiHits.begin(), ed = rphiHits.end(); it != ed; ++it) {
228  matcher->match(&(*it),stereoSimpleHits.begin(),stereoSimpleHits.end(),collectorMatched,gluedDet,trackdirection);
229  if (collectorMatched.size()>0){
230  nmatch+=collectorMatched.size();
232  edm = collectorMatched.end();
233  itm != edm;
234  ++itm) {
235  collector.push_back(*itm);
236  // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
237  matchedSteroClusters.push_back(itm->stereoClusterRef().key());
238  }
239  collectorMatched.clear();
240  } else {
241  // store a copy of this rphi hit as an unmatched rphi hit
242  fillerRphiUnm.push_back(*it);
243  }
244  }
245 
246 #endif
247 
248 
249  // discard matched hits if the collection is empty
250  if (collector.empty()) collector.abort();
251 
252  // discard unmatched rphi hits if there are none
253  if (fillerRphiUnm.empty()) fillerRphiUnm.abort();
254 
255  // now look for unmatched stereo hits
256  SiStripRecHit2DCollection::FastFiller fillerStereoUnm(*output.stereoUnmatched, stereoHits.detId());
257  std::sort(matchedSteroClusters.begin(), matchedSteroClusters.end());
258  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
259  if (!std::binary_search(matchedSteroClusters.begin(), matchedSteroClusters.end(), it->cluster().key())) {
260  fillerStereoUnm.push_back(*it);
261  }
262  }
263  if (fillerStereoUnm.empty()) fillerStereoUnm.abort();
264 
265 
266  }
267 
268  for (SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->begin(); itStereoDet != edStereoDet; ++itStereoDet) {
269  edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
270  StripSubdetector specDetId(stereoHits.detId());
271  uint32_t partnerId = specDetId.partnerDetId();
272  if (partnerId == 0) continue;
273  SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->find(partnerId);
274  if (itRPhiDet == edRPhiDet) {
275  if (!stereoHits.empty()) {
277  filler.resize(stereoHits.size());
278  std::copy(stereoHits.begin(), stereoHits.end(), filler.begin());
279  }
280  }
281  }
282 
283  edm::LogInfo("SiStripRecHitConverter")
284  << "found\n"
285  << nmatch
286  << " matched RecHits\n";
287 }
288 
290 fillBad128StripBlocks(const uint32_t detid, bool bad128StripBlocks[6] ) const
291 {
293  short badApvs = quality->getBadApvs(detid);
294  short badFibers = quality->getBadFibers(detid);
295  for (int j = 0; j < 6; j++) {
296  bad128StripBlocks[j] = (badApvs & (1 << j));
297  }
298  for (int j = 0; j < 3; j++) {
299  if (badFibers & (1 << j)) {
300  bad128StripBlocks[2*j+0] = true;
301  bad128StripBlocks[2*j+1] = true;
302  }
303  }
304  }
305 }
306 
307 inline
309 isMasked(const SiStripCluster &cluster, bool bad128StripBlocks[6]) const
310 {
312  if ( bad128StripBlocks[cluster.firstStrip() >> 7] ) {
313  if ( bad128StripBlocks[(cluster.firstStrip()+cluster.amplitudes().size()) >> 7] ||
314  bad128StripBlocks[static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7] ) {
315  return true;
316  }
317  } else {
318  if ( bad128StripBlocks[(cluster.firstStrip()+cluster.amplitudes().size()) >> 7] &&
319  bad128StripBlocks[static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7] ) {
320  return true;
321  }
322  }
323  }
324  return false;
325 }
326 
327 inline
329 useModule(const uint32_t id) const
330 {
331  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tracker->idToDetUnit(id);
332  if(stripdet==nullptr) edm::LogWarning("SiStripRecHitConverter") << "Detid=" << id << " not found";
333  return stripdet!=nullptr && (!useQuality || quality->IsModuleUsable(id));
334 }
void push_back(data_type const &d)
std::unique_ptr< SiStripMatchedRecHit2DCollection > matched
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
short getBadFibers(const uint32_t &detid) const
Local3DVector LocalVector
Definition: LocalVector.h:12
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
void run(edm::Handle< edmNew::DetSetVector< SiStripCluster > > input, products &output)
edm::ESHandle< SiStripRecHitMatcher > matcher
edm::ESHandle< TrackerGeometry > tracker
virtual void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const
size_type size() const
Definition: OwnVector.h:264
std::unique_ptr< SiStripRecHit2DCollection > rphi
std::pair< LocalPoint, LocalError > LocalValues
uint16_t firstStrip() const
data_type const * const_iterator
Definition: DetSetNew.h:30
SiStripRecHitConverterAlgorithm(const edm::ParameterSet &)
void fillBad128StripBlocks(const uint32_t detid, bool bad128StripBlocks[6]) const
unsigned int partnerDetId() const
static std::string const input
Definition: EdmProvDump.cc:44
iterator begin()
Definition: OwnVector.h:244
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
std::unique_ptr< SiStripRecHit2DCollection > stereo
short getBadApvs(const uint32_t &detid) const
bool IsModuleUsable(const uint32_t &detid) const
std::unique_ptr< SiStripMatchedRecHit2D > match(const SiStripRecHit2D *monoRH, const SiStripRecHit2D *stereoRH, const GluedGeomDet *gluedDet, LocalVector trackdirection, bool force) const
bool empty() const
Definition: DetSetNew.h:90
void clear()
Definition: OwnVector.h:445
void match(products &output, LocalVector trackdirection) const
float barycenter() const
void doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend, StereoIterator seconditer, StereoIterator seconditerend, const GluedGeomDet *gluedDet, LocalVector trdir, CollectorHelper &collectorHelper) const
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:96
Definition: DetId.h:18
unsigned int stereo() const
stereo
iterator end()
Definition: OwnVector.h:249
std::vector< const SiStripRecHit2D * > SimpleHitCollection
const T & get() const
Definition: EventSetup.h:58
std::unique_ptr< SiStripRecHit2DCollection > stereoUnmatched
SiStripRecHit2DCollection::FastFiller Collector
id_type detId() const
Definition: DetSetNew.h:84
HLT enums.
iterator end()
Definition: DetSetNew.h:70
const TrackerGeomDet * idToDet(DetId) const override
edm::ESHandle< StripClusterParameterEstimator > parameterestimator
std::unique_ptr< SiStripRecHit2DCollection > rphiUnmatched
size_type size() const
Definition: DetSetNew.h:87
const std::vector< uint8_t > & amplitudes() const
bool isMasked(const SiStripCluster &cluster, bool bad128StripBlocks[6]) const
edm::ESHandle< SiStripQuality > quality
iterator begin()
Definition: DetSetNew.h:67