CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
62  DS = inputhandle->begin(); DS != dse; ++DS ) {
63  edmNew::det_id_type id = (*DS).id();
64  if(!useModule(id)) continue;
65 
66  Collector collector = StripSubdetector(id).stereo()
67  ? Collector(*output.stereo, id)
68  : Collector(*output.rphi, id);
69 
70  bool bad128StripBlocks[6]; fillBad128StripBlocks( id, bad128StripBlocks);
71 
72  GeomDetUnit const & du = *(tracker->idToDetUnit(id));
75  cluster = (*DS).begin(); cluster != cle; ++cluster ) {
76 
77  if(isMasked(*cluster,bad128StripBlocks)) continue;
78 
80  collector.push_back(SiStripRecHit2D( parameters.first, parameters.second, du, edmNew::makeRefTo(inputhandle,cluster) ));
81  }
82 
83  if (collector.empty()) collector.abort();
84  }
85  match(output,trackdirection);
86 }
87 
88 
89 namespace {
90 
91  struct CollectorHelper {
92  size_t nmatch;
93 
94 
95  typedef edm::OwnVector<SiStripMatchedRecHit2D> CollectorMatched;
97 
98  Collector & m_collector;
99  CollectorMatched & m_collectorMatched;
100  SiStripRecHit2DCollection::FastFiller & m_fillerRphiUnm;
101  std::vector<SiStripRecHit2D::ClusterRef::key_type> & m_matchedSteroClusters;
102 
103  static inline SiStripRecHit2D const & stereoHit(edmNew::DetSet<SiStripRecHit2D>::const_iterator iter) {
104  return *iter;
105  }
106 
108  return *iter;
109  }
110 
111  struct Add {
112  Add(CollectorHelper& ih) : h(ih){}
113  CollectorHelper& h;
114  void operator()(SiStripMatchedRecHit2D const & rh) { h.m_collectorMatched.push_back(rh);}
115  };
116 
117  CollectorHelper & collector() {
118  return *this;
119  }
120 
121  void operator()(SiStripMatchedRecHit2D const & rh) {m_collectorMatched.push_back(rh);}
122 
123 
124  CollectorHelper(
125  Collector & i_collector,
126  CollectorMatched & i_collectorMatched,
127  SiStripRecHit2DCollection::FastFiller & i_fillerRphiUnm,
128  std::vector<SiStripRecHit2D::ClusterRef::key_type> & i_matchedSteroClusters
129  ) : nmatch(0),
130  m_collector(i_collector),
131  m_collectorMatched(i_collectorMatched),
132  m_fillerRphiUnm(i_fillerRphiUnm),
133  m_matchedSteroClusters(i_matchedSteroClusters) {}
134 
136  if (!m_collectorMatched.empty()){
137  nmatch+=m_collectorMatched.size();
138  for (edm::OwnVector<SiStripMatchedRecHit2D>::const_iterator itm = m_collectorMatched.begin(),
139  edm = m_collectorMatched.end();
140  itm != edm;
141  ++itm) {
142  m_collector.push_back(*itm);
143  // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
144  m_matchedSteroClusters.push_back(itm->stereoClusterRef().key());
145  }
146  m_collectorMatched.clear();
147  } else {
148  // store a copy of this rphi hit as an unmatched rphi hit
149  m_fillerRphiUnm.push_back(*it);
150  }
151  }
152  };
153 }
154 
155 
157 match(products& output, LocalVector trackdirection) const
158 {
159  int nmatch=0;
160  edm::OwnVector<SiStripMatchedRecHit2D> collectorMatched; // gp/FIXME: avoid this
161 
162  // Remember the ends of the collections, as we will use them a lot
163  SiStripRecHit2DCollection::const_iterator edStereoDet = output.stereo->end();
164  SiStripRecHit2DCollection::const_iterator edRPhiDet = output.rphi->end();
165 
166  // two work vectors for bookeeping clusters used by the stereo part of the matched hits
167  std::vector<SiStripRecHit2D::ClusterRef::key_type> matchedSteroClusters;
168 
169  for (SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->begin(); itRPhiDet != edRPhiDet; ++itRPhiDet) {
170  edmNew::DetSet<SiStripRecHit2D> rphiHits = *itRPhiDet;
171  StripSubdetector specDetId(rphiHits.detId());
172  uint32_t partnerId = specDetId.partnerDetId();
173 
174  // if not part of a glued pair
175  if (partnerId == 0) {
176  // I must copy these as unmatched
177  if (!rphiHits.empty()) {
178  SiStripRecHit2DCollection::FastFiller filler(*output.rphiUnmatched, rphiHits.detId());
179  filler.resize(rphiHits.size());
180  std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
181  }
182  continue;
183  }
184 
185  SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->find(partnerId);
186 
187  // if the partner is not found (which probably can happen if it's empty)
188  if (itStereoDet == edStereoDet) {
189  // I must copy these as unmatched
190  if (!rphiHits.empty()) {
191  SiStripRecHit2DCollection::FastFiller filler(*output.rphiUnmatched, rphiHits.detId());
192  filler.resize(rphiHits.size());
193  std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
194  }
195  continue;
196  }
197 
198  edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
199 
200 
201  // Get ready for making glued hits
202  const GluedGeomDet* gluedDet = (const GluedGeomDet*)tracker->idToDet(DetId(specDetId.glued()));
204  Collector collector(*output.matched, specDetId.glued());
205 
206  // Prepare also the list for unmatched rphi hits
207  SiStripRecHit2DCollection::FastFiller fillerRphiUnm(*output.rphiUnmatched, rphiHits.detId());
208 
209  // a list of clusters used by the matched part of the stereo hits in this detector
210  matchedSteroClusters.clear(); // at the beginning, empty
211 
212 #ifdef DOUBLE_MATCH
213  CollectorHelper chelper(collector, collectorMatched,
214  fillerRphiUnm,
215  matchedSteroClusters
216  );
217  matcher->doubleMatch(rphiHits.begin(), rphiHits.end(),
218  stereoHits.begin(),stereoHits.end(),gluedDet,trackdirection,chelper);
219  nmatch+=chelper.nmatch;
220 #else
221  // Make simple collection of this (gp:FIXME: why do we need it?)
223  // gp:FIXME: use std::transform
224  stereoSimpleHits.reserve(stereoHits.size());
225  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
226  stereoSimpleHits.push_back(&*it);
227  }
228 
229  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = rphiHits.begin(), ed = rphiHits.end(); it != ed; ++it) {
230  matcher->match(&(*it),stereoSimpleHits.begin(),stereoSimpleHits.end(),collectorMatched,gluedDet,trackdirection);
231  if (collectorMatched.size()>0){
232  nmatch+=collectorMatched.size();
234  edm = collectorMatched.end();
235  itm != edm;
236  ++itm) {
237  collector.push_back(*itm);
238  // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
239  matchedSteroClusters.push_back(itm->stereoClusterRef().key());
240  }
241  collectorMatched.clear();
242  } else {
243  // store a copy of this rphi hit as an unmatched rphi hit
244  fillerRphiUnm.push_back(*it);
245  }
246  }
247 
248 #endif
249 
250 
251  // discard matched hits if the collection is empty
252  if (collector.empty()) collector.abort();
253 
254  // discard unmatched rphi hits if there are none
255  if (fillerRphiUnm.empty()) fillerRphiUnm.abort();
256 
257  // now look for unmatched stereo hits
258  SiStripRecHit2DCollection::FastFiller fillerStereoUnm(*output.stereoUnmatched, stereoHits.detId());
259  std::sort(matchedSteroClusters.begin(), matchedSteroClusters.end());
260  for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
261  if (!std::binary_search(matchedSteroClusters.begin(), matchedSteroClusters.end(), it->cluster().key())) {
262  fillerStereoUnm.push_back(*it);
263  }
264  }
265  if (fillerStereoUnm.empty()) fillerStereoUnm.abort();
266 
267 
268  }
269 
270  for (SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->begin(); itStereoDet != edStereoDet; ++itStereoDet) {
271  edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
272  StripSubdetector specDetId(stereoHits.detId());
273  uint32_t partnerId = specDetId.partnerDetId();
274  if (partnerId == 0) continue;
275  SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->find(partnerId);
276  if (itRPhiDet == edRPhiDet) {
277  if (!stereoHits.empty()) {
278  SiStripRecHit2DCollection::FastFiller filler(*output.stereoUnmatched, stereoHits.detId());
279  filler.resize(stereoHits.size());
280  std::copy(stereoHits.begin(), stereoHits.end(), filler.begin());
281  }
282  }
283  }
284 
285  edm::LogInfo("SiStripRecHitConverter")
286  << "found\n"
287  << nmatch
288  << " matched RecHits\n";
289 }
290 
292 fillBad128StripBlocks(const uint32_t detid, bool bad128StripBlocks[6] ) const
293 {
295  short badApvs = quality->getBadApvs(detid);
296  short badFibers = quality->getBadFibers(detid);
297  for (int j = 0; j < 6; j++) {
298  bad128StripBlocks[j] = (badApvs & (1 << j));
299  }
300  for (int j = 0; j < 3; j++) {
301  if (badFibers & (1 << j)) {
302  bad128StripBlocks[2*j+0] = true;
303  bad128StripBlocks[2*j+1] = true;
304  }
305  }
306  }
307 }
308 
309 inline
311 isMasked(const SiStripCluster &cluster, bool bad128StripBlocks[6]) const
312 {
314  if ( bad128StripBlocks[cluster.firstStrip() >> 7] ) {
315  if ( bad128StripBlocks[(cluster.firstStrip()+cluster.amplitudes().size()) >> 7] ||
316  bad128StripBlocks[static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7] ) {
317  return true;
318  }
319  } else {
320  if ( bad128StripBlocks[(cluster.firstStrip()+cluster.amplitudes().size()) >> 7] &&
321  bad128StripBlocks[static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7] ) {
322  return true;
323  }
324  }
325  }
326  return false;
327 }
328 
329 inline
331 useModule(const uint32_t id) const
332 {
333  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tracker->idToDetUnit(id);
334  if(stripdet==0) edm::LogWarning("SiStripRecHitConverter") << "Detid=" << id << " not found";
335  return stripdet!=0 && (!useQuality || quality->IsModuleUsable(id));
336 }
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
void push_back(data_type const &d)
dictionary parameters
Definition: Parameters.py:2
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Local3DVector LocalVector
Definition: LocalVector.h:12
void run(edm::Handle< edmNew::DetSetVector< SiStripCluster > > input, products &output)
unsigned int det_id_type
Definition: DetSetNew.h:10
std::auto_ptr< SiStripRecHit2DCollection > rphiUnmatched
edm::ESHandle< SiStripRecHitMatcher > matcher
edm::ESHandle< TrackerGeometry > tracker
size_type size() const
Definition: OwnVector.h:254
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:43
iterator begin()
Definition: OwnVector.h:234
bool empty() const
Definition: DetSetNew.h:89
void clear()
Definition: OwnVector.h:377
void match(products &output, LocalVector trackdirection) const
int j
Definition: DBlmapReader.cc:9
float barycenter() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::auto_ptr< SiStripRecHit2DCollection > stereoUnmatched
std::auto_ptr< SiStripRecHit2DCollection > rphi
tuple conf
Definition: dbtoconf.py:185
#define TYPELOOKUP_DATA_REG(_dataclass_)
Definition: typelookup.h:96
Definition: DetId.h:18
unsigned int stereo() const
stereo
iterator end()
Definition: OwnVector.h:239
std::vector< const SiStripRecHit2D * > SimpleHitCollection
const T & get() const
Definition: EventSetup.h:55
SiStripRecHit2DCollection::FastFiller Collector
id_type detId() const
Definition: DetSetNew.h:83
std::pair< LocalPoint, LocalError > LocalValues
iterator end()
Definition: DetSetNew.h:70
edm::ESHandle< StripClusterParameterEstimator > parameterestimator
size_type size() const
Definition: DetSetNew.h:86
std::auto_ptr< SiStripMatchedRecHit2DCollection > matched
const std::vector< uint8_t > & amplitudes() const
std::auto_ptr< SiStripRecHit2DCollection > stereo
bool isMasked(const SiStripCluster &cluster, bool bad128StripBlocks[6]) const
edm::ESHandle< SiStripQuality > quality
iterator begin()
Definition: DetSetNew.h:67