CMS 3D CMS Logo

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