CMS 3D CMS Logo

MuonSeedOrcaPatternRecognition.cc
Go to the documentation of this file.
1 
13 
20 
23 
24 
25 // Geometry
28 
31 
32 // Framework
38 
39 // C++
40 #include <vector>
41 
42 using namespace std;
43 
44  const std::string metname = "Muon|RecoMuon|MuonSeedOrcaPatternRecognition";
45 
46 // Constructor
49  theCrackEtas(pset.getParameter<std::vector<double> >("crackEtas")),
50  theCrackWindow(pset.getParameter<double>("crackWindow")),
51  theDeltaPhiWindow(pset.existsAs<double>("deltaPhiSearchWindow") ? pset.getParameter<double>("deltaPhiSearchWindow") : 0.25),
52  theDeltaEtaWindow(pset.existsAs<double>("deltaEtaSearchWindow") ? pset.getParameter<double>("deltaEtaSearchWindow") : 0.2),
53 theDeltaCrackWindow(pset.existsAs<double>("deltaEtaCrackSearchWindow") ? pset.getParameter<double>("deltaEtaCrackSearchWindow") : 0.25)
54 {
55 
57  iC,
58  enableDTMeasurement,enableCSCMeasurement,false,false,false);
59 }
60 
61 
62 // reconstruct muon's seeds
64  std::vector<MuonRecHitContainer> & result)
65 {
66 
67  // divide the RecHits by DetLayer, in order to fill the
68  // RecHitContainer like it was in ORCA
69 
70  // Muon Geometry - DT, CSC and RPC
72  eSetup.get<MuonRecoGeometryRecord>().get(muonLayers);
73 
74  // get the DT layers
75  vector<const DetLayer*> dtLayers = muonLayers->allDTLayers();
76 
77  // get the CSC layers
78  vector<const DetLayer*> cscForwardLayers = muonLayers->forwardCSCLayers();
79  vector<const DetLayer*> cscBackwardLayers = muonLayers->backwardCSCLayers();
80 
81  // Backward (z<0) EndCap disk
82  const DetLayer* ME4Bwd = cscBackwardLayers[4];
83  const DetLayer* ME3Bwd = cscBackwardLayers[3];
84  const DetLayer* ME2Bwd = cscBackwardLayers[2];
85  const DetLayer* ME12Bwd = cscBackwardLayers[1];
86  const DetLayer* ME11Bwd = cscBackwardLayers[0];
87 
88  // Forward (z>0) EndCap disk
89  const DetLayer* ME11Fwd = cscForwardLayers[0];
90  const DetLayer* ME12Fwd = cscForwardLayers[1];
91  const DetLayer* ME2Fwd = cscForwardLayers[2];
92  const DetLayer* ME3Fwd = cscForwardLayers[3];
93  const DetLayer* ME4Fwd = cscForwardLayers[4];
94 
95  // barrel
96  const DetLayer* MB4DL = dtLayers[3];
97  const DetLayer* MB3DL = dtLayers[2];
98  const DetLayer* MB2DL = dtLayers[1];
99  const DetLayer* MB1DL = dtLayers[0];
100 
101  // instantiate the accessor
102  // Don not use RPC for seeding
103 
104  double barreldThetaCut = 0.2;
105  // still lose good muons to a tighter cut
106  double endcapdThetaCut = 1.0;
107  MuonRecHitContainer list9 = filterSegments(muonMeasurements->recHits(MB4DL,event), barreldThetaCut);
108  MuonRecHitContainer list6 = filterSegments(muonMeasurements->recHits(MB3DL,event), barreldThetaCut);
109  MuonRecHitContainer list7 = filterSegments(muonMeasurements->recHits(MB2DL,event), barreldThetaCut);
110  MuonRecHitContainer list8 = filterSegments(muonMeasurements->recHits(MB1DL,event), barreldThetaCut);
111 
112  dumpLayer("MB4 ", list9);
113  dumpLayer("MB3 ", list6);
114  dumpLayer("MB2 ", list7);
115  dumpLayer("MB1 ", list8);
116 
117  bool* MB1 = zero(list8.size());
118  bool* MB2 = zero(list7.size());
119  bool* MB3 = zero(list6.size());
120 
121  endcapPatterns(filterSegments(muonMeasurements->recHits(ME11Bwd,event), endcapdThetaCut),
122  filterSegments(muonMeasurements->recHits(ME12Bwd,event), endcapdThetaCut),
123  filterSegments(muonMeasurements->recHits(ME2Bwd,event), endcapdThetaCut),
124  filterSegments(muonMeasurements->recHits(ME3Bwd,event), endcapdThetaCut),
125  filterSegments(muonMeasurements->recHits(ME4Bwd,event), endcapdThetaCut),
126  list8, list7, list6,
127  MB1, MB2, MB3, result);
128 
129  endcapPatterns(filterSegments(muonMeasurements->recHits(ME11Fwd,event), endcapdThetaCut),
130  filterSegments(muonMeasurements->recHits(ME12Fwd,event), endcapdThetaCut),
131  filterSegments(muonMeasurements->recHits(ME2Fwd,event), endcapdThetaCut),
132  filterSegments(muonMeasurements->recHits(ME3Fwd,event), endcapdThetaCut),
133  filterSegments(muonMeasurements->recHits(ME4Fwd,event), endcapdThetaCut),
134  list8, list7, list6,
135  MB1, MB2, MB3, result);
136 
137 
138  // ---------- Barrel only
139 
140  unsigned int counter = 0;
141  if ( list9.size() < 100 ) { // +v
142  for (MuonRecHitContainer::iterator iter=list9.begin(); iter!=list9.end(); iter++ ){
143  MuonRecHitContainer seedSegments;
144  seedSegments.push_back(*iter);
145  complete(seedSegments, list6, MB3);
146  complete(seedSegments, list7, MB2);
147  complete(seedSegments, list8, MB1);
148  if(check(seedSegments)) result.push_back(seedSegments);
149  }
150  }
151 
152 
153  if ( list6.size() < 100 ) { // +v
154  for ( counter = 0; counter<list6.size(); counter++ ){
155  if ( !MB3[counter] ) {
156  MuonRecHitContainer seedSegments;
157  seedSegments.push_back(list6[counter]);
158  complete(seedSegments, list7, MB2);
159  complete(seedSegments, list8, MB1);
160  complete(seedSegments, list9);
161  if(check(seedSegments)) result.push_back(seedSegments);
162  }
163  }
164  }
165 
166 
167  if ( list7.size() < 100 ) { // +v
168  for ( counter = 0; counter<list7.size(); counter++ ){
169  if ( !MB2[counter] ) {
170  MuonRecHitContainer seedSegments;
171  seedSegments.push_back(list7[counter]);
172  complete(seedSegments, list8, MB1);
173  complete(seedSegments, list9);
174  complete(seedSegments, list6, MB3);
175  if (seedSegments.size()>1 ||
176  (seedSegments.size()==1 && seedSegments[0]->dimension()==4) )
177  {
178  result.push_back(seedSegments);
179  }
180  }
181  }
182  }
183 
184 
185  if ( list8.size() < 100 ) { // +v
186  for ( counter = 0; counter<list8.size(); counter++ ){
187  if ( !MB1[counter] ) {
188  MuonRecHitContainer seedSegments;
189  seedSegments.push_back(list8[counter]);
190  complete(seedSegments, list9);
191  complete(seedSegments, list6, MB3);
192  complete(seedSegments, list7, MB2);
193  if (seedSegments.size()>1 ||
194  (seedSegments.size()==1 && seedSegments[0]->dimension()==4) )
195  {
196  result.push_back(seedSegments);
197  }
198  }
199  }
200  }
201 
202  if ( MB3 ) delete [] MB3;
203  if ( MB2 ) delete [] MB2;
204  if ( MB1 ) delete [] MB1;
205 
206  if(result.empty())
207  {
208  // be a little stricter with single segment seeds
209  barreldThetaCut = 0.2;
210  endcapdThetaCut = 0.2;
211 
213  MuonRecHitContainer tmp = filterSegments(muonMeasurements->recHits(ME3Bwd,event), endcapdThetaCut);
214  copy(tmp.begin(),tmp.end(),back_inserter(all));
215 
216  tmp = filterSegments(muonMeasurements->recHits(ME2Bwd,event), endcapdThetaCut);
217  copy(tmp.begin(),tmp.end(),back_inserter(all));
218 
219  tmp = filterSegments(muonMeasurements->recHits(ME12Bwd,event), endcapdThetaCut);
220  copy(tmp.begin(),tmp.end(),back_inserter(all));
221 
222  tmp = filterSegments(muonMeasurements->recHits(ME11Bwd,event), endcapdThetaCut);
223  copy(tmp.begin(),tmp.end(),back_inserter(all));
224 
225  tmp = filterSegments(muonMeasurements->recHits(ME11Fwd,event), endcapdThetaCut);
226  copy(tmp.begin(),tmp.end(),back_inserter(all));
227 
228  tmp = filterSegments(muonMeasurements->recHits(ME12Fwd,event), endcapdThetaCut);
229  copy(tmp.begin(),tmp.end(),back_inserter(all));
230 
231  tmp = filterSegments(muonMeasurements->recHits(ME2Fwd,event), endcapdThetaCut);
232  copy(tmp.begin(),tmp.end(),back_inserter(all));
233 
234  tmp = filterSegments(muonMeasurements->recHits(ME3Fwd,event), endcapdThetaCut);
235  copy(tmp.begin(),tmp.end(),back_inserter(all));
236 
237  tmp = filterSegments(muonMeasurements->recHits(ME4Fwd,event), endcapdThetaCut);
238  copy(tmp.begin(),tmp.end(),back_inserter(all));
239 
240  tmp = filterSegments(muonMeasurements->recHits(MB4DL,event), barreldThetaCut);
241  copy(tmp.begin(),tmp.end(),back_inserter(all));
242 
243  tmp = filterSegments(muonMeasurements->recHits(MB3DL,event), barreldThetaCut);
244  copy(tmp.begin(),tmp.end(),back_inserter(all));
245 
246  tmp = filterSegments(muonMeasurements->recHits(MB2DL,event), barreldThetaCut);
247  copy(tmp.begin(),tmp.end(),back_inserter(all));
248 
249  tmp = filterSegments(muonMeasurements->recHits(MB1DL,event), barreldThetaCut);
250  copy(tmp.begin(),tmp.end(),back_inserter(all));
251 
252  LogTrace(metname)<<"Number of segments: "<<all.size();
253 
254  for(MuonRecHitContainer::const_iterator segmentItr = all.begin();
255  segmentItr != all.end(); ++segmentItr)
256  {
257  MuonRecHitContainer singleSegmentContainer;
258  singleSegmentContainer.push_back(*segmentItr);
259  result.push_back(singleSegmentContainer);
260  }
261  }
262 
263 }
264 
265 
266 bool * MuonSeedOrcaPatternRecognition::zero(unsigned listSize)
267 {
268  bool * result = 0;
269  if (listSize) {
270  result = new bool[listSize];
271  for ( size_t i=0; i<listSize; i++ ) result[i]=false;
272  }
273  return result;
274 }
275 
276 
278  const MuonRecHitContainer & me11, const MuonRecHitContainer & me12,
279  const MuonRecHitContainer & me2, const MuonRecHitContainer & me3,
280  const MuonRecHitContainer & me4, const MuonRecHitContainer & mb1,
282  bool * MB1, bool * MB2, bool * MB3,
283  std::vector<MuonRecHitContainer> & result)
284 {
285  dumpLayer("ME4 ", me4);
286  dumpLayer("ME3 ", me3);
287  dumpLayer("ME2 ", me2);
288  dumpLayer("ME12 ", me12);
289  dumpLayer("ME11 ", me11);
290 
291  std::vector<MuonRecHitContainer> patterns;
292  MuonRecHitContainer crackSegments;
293  rememberCrackSegments(me11, crackSegments);
294  rememberCrackSegments(me12, crackSegments);
295  rememberCrackSegments(me2, crackSegments);
296  rememberCrackSegments(me3, crackSegments);
297  rememberCrackSegments(me4, crackSegments);
298 
299 
300  MuonRecHitContainer list24 = me4;
301  MuonRecHitContainer list23 = me3;
302 
303  MuonRecHitContainer list12 = me2;
304 
305  MuonRecHitContainer list22 = me12;
306  MuonRecHitContainer list21 = me11;
307 
308  MuonRecHitContainer list11 = list21;
309  MuonRecHitContainer list5 = list22;
310  MuonRecHitContainer list13 = list23;
311  MuonRecHitContainer list4 = list24;
312 
313  if ( list21.size() == 0 ) {
314  list11 = list22; list5 = list21;
315  }
316 
317  if ( list24.size() < list23.size() && list24.size() > 0 ) {
318  list13 = list24; list4 = list23;
319  }
320 
321  if ( list23.size() == 0 ) {
322  list13 = list24; list4 = list23;
323  }
324 
325  MuonRecHitContainer list1 = list11;
326  MuonRecHitContainer list2 = list12;
327  MuonRecHitContainer list3 = list13;
328 
329 
330  if ( list12.size() == 0 ) {
331  list3 = list12;
332  if ( list11.size() <= list13.size() && list11.size() > 0 ) {
333  list1 = list11; list2 = list13;}
334  else { list1 = list13; list2 = list11;}
335  }
336 
337  if ( list13.size() == 0 ) {
338  if ( list11.size() <= list12.size() && list11.size() > 0 ) {
339  list1 = list11; list2 = list12;}
340  else { list1 = list12; list2 = list11;}
341  }
342 
343  if ( list12.size() != 0 && list13.size() != 0 ) {
344  if ( list11.size()<=list12.size() && list11.size()<=list13.size() && list11.size()>0 ) { // ME 1
345  if ( list12.size() > list13.size() ) {
346  list2 = list13; list3 = list12;}
347  }
348  else if ( list12.size() <= list13.size() ) { // start with ME 2
349  list1 = list12;
350  if ( list11.size() <= list13.size() && list11.size() > 0 ) {
351  list2 = list11; list3 = list13;}
352  else { list2 = list13; list3 = list11;}
353  }
354  else { // start with ME 3
355  list1 = list13;
356  if ( list11.size() <= list12.size() && list11.size() > 0 ) {
357  list2 = list11; list3 = list12;}
358  else { list2 = list12; list3 = list11;}
359  }
360  }
361 
362 
363  bool* ME2 = zero(list2.size());
364  bool* ME3 = zero(list3.size());
365  bool* ME4 = zero(list4.size());
366  bool* ME5 = zero(list5.size());
367 
368 
369  // creates list of compatible track segments
370 
371  for (MuonRecHitContainer::iterator iter = list1.begin(); iter!=list1.end(); iter++ ){
372  if ( (*iter)->recHits().size() < 4 && list3.size() > 0 ) continue; // 3p.tr-seg. are not so good for starting
373  MuonRecHitContainer seedSegments;
374  seedSegments.push_back(*iter);
375  complete(seedSegments, list2, ME2);
376  complete(seedSegments, list3, ME3);
377  complete(seedSegments, list4, ME4);
378  complete(seedSegments, list5, ME5);
379  complete(seedSegments, mb3, MB3);
380  complete(seedSegments, mb2, MB2);
381  complete(seedSegments, mb1, MB1);
382  if(check(seedSegments)) patterns.push_back(seedSegments);
383  }
384 
385 
386  unsigned int counter;
387 
388  for ( counter = 0; counter<list2.size(); counter++ ){
389 
390  if ( !ME2[counter] ) {
391  MuonRecHitContainer seedSegments;
392  seedSegments.push_back(list2[counter]);
393  complete(seedSegments, list3, ME3);
394  complete(seedSegments, list4, ME4);
395  complete(seedSegments, list5, ME5);
396  complete(seedSegments, mb3, MB3);
397  complete(seedSegments, mb2, MB2);
398  complete(seedSegments, mb1, MB1);
399  if(check(seedSegments)) patterns.push_back(seedSegments);
400  }
401  }
402 
403 
404  if ( list3.size() < 20 ) { // +v
405  for ( counter = 0; counter<list3.size(); counter++ ){
406  if ( !ME3[counter] ) {
407  MuonRecHitContainer seedSegments;
408  seedSegments.push_back(list3[counter]);
409  complete(seedSegments, list4, ME4);
410  complete(seedSegments, list5, ME5);
411  complete(seedSegments, mb3, MB3);
412  complete(seedSegments, mb2, MB2);
413  complete(seedSegments, mb1, MB1);
414  if(check(seedSegments)) patterns.push_back(seedSegments);
415  }
416  }
417  }
418 
419  if ( list4.size() < 20 ) { // +v
420  for ( counter = 0; counter<list4.size(); counter++ ){
421  if ( !ME4[counter] ) {
422  MuonRecHitContainer seedSegments;
423  seedSegments.push_back(list4[counter]);
424  complete(seedSegments, list5, ME5);
425  complete(seedSegments, mb3, MB3);
426  complete(seedSegments, mb2, MB2);
427  complete(seedSegments, mb1, MB1);
428  if(check(seedSegments)) patterns.push_back(seedSegments);
429  }
430  }
431  }
432 
433  if ( ME5 ) delete [] ME5;
434  if ( ME4 ) delete [] ME4;
435  if ( ME3 ) delete [] ME3;
436  if ( ME2 ) delete [] ME2;
437 
438  if(!patterns.empty())
439  {
440  result.insert(result.end(), patterns.begin(), patterns.end());
441  }
442  else
443  {
444  if(!crackSegments.empty())
445  {
446  // make some single-segment seeds
447  for(MuonRecHitContainer::const_iterator crackSegmentItr = crackSegments.begin();
448  crackSegmentItr != crackSegments.end(); ++crackSegmentItr)
449  {
450  MuonRecHitContainer singleSegmentPattern;
451  singleSegmentPattern.push_back(*crackSegmentItr);
452  result.push_back(singleSegmentPattern);
453  }
454  }
455  }
456 }
457 
458 
459 
461  const MuonRecHitContainer &recHits, bool* used) const {
462 
463  MuonRecHitContainer good_rhit;
464  MuonPatternRecoDumper theDumper;
465  //+v get all rhits compatible with the seed on dEta/dPhi Glob.
466 
467  ConstMuonRecHitPointer first = seedSegments[0]; // first rechit of seed
468 
469  GlobalPoint ptg2 = first->globalPosition(); // its global pos +v
470 
471  for (unsigned nr = 0; nr < recHits.size(); ++nr ){
472  MuonRecHitPointer recHit(recHits[nr]);
473  GlobalPoint ptg1(recHit->globalPosition());
474  float deta = fabs (ptg1.eta()-ptg2.eta());
475  // Geom::Phi should keep it in the range [-pi, pi]
476  float dphi = fabs( deltaPhi(ptg1.barePhi(), ptg2.barePhi()) );
477  // be a little more lenient in cracks
478  bool crack = isCrack(recHit) || isCrack(first);
479  //float detaWindow = 0.3;
480  float detaWindow = crack ? theDeltaCrackWindow : theDeltaEtaWindow;
481  if ( deta > detaWindow || dphi > theDeltaPhiWindow ) {
482  continue;
483  } // +vvp!!!
484 
485  good_rhit.push_back(recHit);
486  if (used) markAsUsed(nr, recHits, used);
487  } // recHits iter
488 
489  // select the best rhit among the compatible ones (based on Dphi Glob & Dir)
490  MuonRecHitPointer best=bestMatch(first, good_rhit);
491  if(best && best->isValid() ) seedSegments.push_back(best);
492 }
493 
494 
495 
498  MuonRecHitContainer & good_rhit) const
499 {
500  MuonRecHitPointer best = 0;
501  if(good_rhit.size() == 1) return good_rhit[0];
502  double bestDiscrim = 10000.;
503  for (MuonRecHitContainer::iterator iter=good_rhit.begin();
504  iter!=good_rhit.end(); iter++)
505  {
506  double discrim = discriminator(first, *iter);
507  if(discrim < bestDiscrim)
508  {
509  bestDiscrim = discrim;
510  best = *iter;
511  }
512  }
513  return best;
514 }
515 
516 
518 {
519  GlobalPoint gp1= first->globalPosition();
520  GlobalPoint gp2= other->globalPosition();
521  GlobalVector gd1 = first->globalDirection();
522  GlobalVector gd2 = other->globalDirection();
523  if(first->isDT() || other->isDT()) {
524  return fabs(deltaPhi(gd1.barePhi(), gd2.barePhi()));
525  }
526 
527  // penalize those 3-hit segments
528  int nhits = other->recHits().size();
529  int penalty = std::max(nhits-2, 1);
530  float dphig = deltaPhi(gp1.barePhi(), gp2.barePhi());
531  // ME1A has slanted wires, so matching theta position doesn't work well.
532  if(isME1A(first) || isME1A(other)) {
533  return fabs(dphig/penalty);
534  }
535 
536  float dthetag = gp1.theta()-gp2.theta();
537  float dphid2 = fabs(deltaPhi(gd2.barePhi(), gp2.barePhi()));
538  if (dphid2 > M_PI*.5) dphid2 = M_PI - dphid2; //+v
539  float dthetad2 = gp2.theta()-gd2.theta();
540  // for CSC, make a big chi-squared of relevant variables
541  // FIXME for 100 GeV mnd above muons, this doesn't perform as well as
542  // previous methods. needs investigation.
543  float chisq = ((dphig/0.02)*(dphig/0.02)
544  + (dthetag/0.003)*(dthetag/0.003)
545  + (dphid2/0.06)*(dphid2/0.06)
546  + (dthetad2/0.08)*(dthetad2/0.08)
547  );
548  return chisq / penalty;
549 }
550 
551 
553 {
554  return (segments.size() > 1);
555 }
556 
557 
558 void MuonSeedOrcaPatternRecognition::markAsUsed(int nr, const MuonRecHitContainer &recHits, bool* used) const
559 {
560  used[nr] = true;
561  // if it's ME1A with two other segments in the container, mark the ghosts as used, too.
562  if(recHits[nr]->isCSC())
563  {
564  CSCDetId detId(recHits[nr]->geographicalId().rawId());
565  if(detId.ring() == 4)
566  {
567  std::vector<unsigned> chamberHitNs;
568  for(unsigned i = 0; i < recHits.size(); ++i)
569  {
570  if(recHits[i]->geographicalId().rawId() == detId.rawId())
571  {
572  chamberHitNs.push_back(i);
573  }
574  }
575  if(chamberHitNs.size() == 3)
576  {
577  for(unsigned i = 0; i < 3; ++i)
578  {
579  used[chamberHitNs[i]] = true;
580  }
581  }
582  }
583  }
584 }
585 
586 
588 {
589  bool result = false;
590  double absEta = fabs(segment->globalPosition().eta());
591  for(std::vector<double>::const_iterator crackItr = theCrackEtas.begin();
592  crackItr != theCrackEtas.end(); ++crackItr)
593  {
594  if(fabs(absEta-*crackItr) < theCrackWindow) {
595  result = true;
596  }
597  }
598  return result;
599 }
600 
601 
603  MuonRecHitContainer & crackSegments) const
604 {
605  for(MuonRecHitContainer::const_iterator segmentItr = segments.begin();
606  segmentItr != segments.end(); ++segmentItr)
607  {
608  if((**segmentItr).hit()->dimension() == 4 && isCrack(*segmentItr))
609  {
610  crackSegments.push_back(*segmentItr);
611  }
612  }
613 }
614 
615 
616 
617 void MuonSeedOrcaPatternRecognition::dumpLayer(const char * name, const MuonRecHitContainer & segments) const
618 {
619  MuonPatternRecoDumper theDumper;
620 
621  LogTrace(metname) << name << std::endl;
622  for(MuonRecHitContainer::const_iterator segmentItr = segments.begin();
623  segmentItr != segments.end(); ++segmentItr)
624  {
625  LogTrace(metname) << theDumper.dumpMuonId((**segmentItr).geographicalId());
626  }
627 }
628 
629 
632 {
633 MuonPatternRecoDumper theDumper;
635  for(MuonRecHitContainer::const_iterator segmentItr = segments.begin();
636  segmentItr != segments.end(); ++segmentItr)
637  {
638  double dtheta = (*segmentItr)->globalDirection().theta() - (*segmentItr)->globalPosition().theta();
639  if((*segmentItr)->isDT())
640  {
641  // only apply the cut to 4D segments
642  if((*segmentItr)->dimension() == 2 || fabs(dtheta) < dThetaCut)
643  {
644  result.push_back(*segmentItr);
645  }
646  else
647  {
648  LogTrace(metname) << "Cutting segment " << theDumper.dumpMuonId((**segmentItr).geographicalId()) << " because dtheta = " << dtheta;
649  }
650 
651  }
652  else if((*segmentItr)->isCSC())
653  {
654  if(fabs(dtheta) < dThetaCut)
655  {
656  result.push_back(*segmentItr);
657  }
658  else
659  {
660  LogTrace(metname) << "Cutting segment " << theDumper.dumpMuonId((**segmentItr).geographicalId()) << " because dtheta = " << dtheta;
661  }
662  }
663  }
665  return result;
666 }
667 
668 
670 {
671  if(segments.empty()) return;
672  MuonPatternRecoDumper theDumper;
673  // need to optimize cuts
674  double dphiCut = 0.05;
675  double detaCut = 0.05;
676  std::vector<unsigned> toKill;
677  std::vector<unsigned> me1aOverlaps;
678  // loop over all segment pairs to see if there are two that match up in eta and phi
679  // but from different chambers
680  unsigned nseg = segments.size();
681  for(unsigned i = 0; i < nseg-1; ++i)
682  {
683  GlobalPoint pg1 = segments[i]->globalPosition();
684  for(unsigned j = i+1; j < nseg; ++j)
685  {
686  GlobalPoint pg2 = segments[j]->globalPosition();
687  if(segments[i]->geographicalId().rawId() != segments[j]->geographicalId().rawId()
688  && fabs(deltaPhi(pg1.barePhi(), pg2.barePhi())) < dphiCut
689  && fabs(pg1.eta()-pg2.eta()) < detaCut)
690  {
691  LogTrace(metname) << "OVERLAP " << theDumper.dumpMuonId(segments[i]->geographicalId()) << " " <<
692  theDumper.dumpMuonId(segments[j]->geographicalId());
693  // see which one is best
694  toKill.push_back( (countHits(segments[i]) < countHits(segments[j])) ? i : j);
695  if(isME1A(segments[i]))
696  {
697  me1aOverlaps.push_back(i);
698  me1aOverlaps.push_back(j);
699  }
700  }
701  }
702  }
703 
704  if(toKill.empty()) return;
705 
706  // try to kill ghosts assigned to overlaps
707  for(unsigned i = 0; i < me1aOverlaps.size(); ++i)
708  {
709  DetId detId(segments[me1aOverlaps[i]]->geographicalId());
710  vector<unsigned> inSameChamber;
711  for(unsigned j = 0; j < nseg; ++j)
712  {
713  if(i != j && segments[j]->geographicalId() == detId)
714  {
715  inSameChamber.push_back(j);
716  }
717  }
718  if(inSameChamber.size() == 2)
719  {
720  toKill.push_back(inSameChamber[0]);
721  toKill.push_back(inSameChamber[1]);
722  }
723  }
724  // now kill the killable
726  for(unsigned i = 0; i < nseg; ++i)
727  {
728  if(std::find(toKill.begin(), toKill.end(), i) == toKill.end())
729  {
730  result.push_back(segments[i]);
731  }
732 
733  }
734  segments.swap(result);
735 }
736 
737 
739 {
740  return segment->isCSC() && CSCDetId(segment->geographicalId()).ring() == 4;
741 }
742 
743 
745  int count = 0;
746  vector<TrackingRecHit*> components = (*segment).recHits();
747  for(vector<TrackingRecHit*>::const_iterator component = components.begin();
748  component != components.end(); ++component)
749  {
750  int componentSize = (**component).recHits().size();
751  count += (componentSize == 0) ? 1 : componentSize;
752  }
753  return count;
754 }
755 
756 
MuonTransientTrackingRecHit::ConstMuonRecHitPointer ConstMuonRecHitPointer
MuonRecHitPointer bestMatch(const ConstMuonRecHitPointer &first, MuonRecHitContainer &good_rhit) const
MuonTransientTrackingRecHit::MuonRecHitPointer MuonRecHitPointer
const std::string metname
void filterOverlappingChambers(MuonRecHitContainer &segments) const
bool isCrack(const ConstMuonRecHitPointer &segment) const
bool check(const MuonRecHitContainer &segments)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
const PhiMemoryImage patterns[9]
std::string dumpMuonId(const DetId &id) const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
T barePhi() const
Definition: PV3DBase.h:68
bool enableCSCMeasurement
Enable the CSC measurement.
int countHits(const MuonRecHitPointer &segment) const
MuonRecHitContainer recHits(const DetLayer *layer, const edm::Event &iEvent)
returns the rechits which are on the layer
void produce(const edm::Event &event, const edm::EventSetup &eSetup, std::vector< MuonRecHitContainer > &result)
void endcapPatterns(const MuonRecHitContainer &me11, const MuonRecHitContainer &me12, const MuonRecHitContainer &me2, const MuonRecHitContainer &me3, const MuonRecHitContainer &me4, const MuonRecHitContainer &mb1, const MuonRecHitContainer &mb2, const MuonRecHitContainer &mb3, bool *MB1, bool *MB2, bool *MB3, std::vector< MuonRecHitContainer > &result)
void rememberCrackSegments(const MuonRecHitContainer &segments, MuonRecHitContainer &crackSegments) const
#define LogTrace(id)
#define M_PI
int ring() const
Definition: CSCDetId.h:75
Definition: DetId.h:18
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
const std::vector< const DetLayer * > & allDTLayers() const
return the DT DetLayers (barrel), inside-out
bool isME1A(const ConstMuonRecHitPointer &segment) const
void complete(MuonRecHitContainer &seedSegments, const MuonRecHitContainer &recHits, bool *used=0) const
const std::vector< const DetLayer * > & forwardCSCLayers() const
return the forward (+Z) CSC DetLayers, inside-out
const T & get() const
Definition: EventSetup.h:56
std::string const & label() const
Definition: InputTag.h:36
T eta() const
Definition: PV3DBase.h:76
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
void markAsUsed(int nr, const MuonRecHitContainer &recHits, bool *used) const
MuonSeedOrcaPatternRecognition(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
static std::atomic< unsigned int > counter
bool enableDTMeasurement
Enable the DT measurement.
const std::vector< const DetLayer * > & backwardCSCLayers() const
return the backward (-Z) CSC DetLayers, inside-out
void dumpLayer(const char *name, const MuonRecHitContainer &segments) const
double discriminator(const ConstMuonRecHitPointer &first, MuonRecHitPointer &other) const
bool isCSC(const GeomDetEnumerators::SubDetector m)
edm::InputTag theCSCRecSegmentLabel
the name of the CSC rec hits collection
Definition: event.py:1
edm::InputTag theDTRecSegmentLabel
the name of the DT rec hits collection
MuonRecHitContainer filterSegments(const MuonRecHitContainer &segments, double dThetaCut) const
apply some cuts to segments before using them