CMS 3D CMS Logo

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