CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonSelectors.cc
Go to the documentation of this file.
7 namespace muon {
9 {
10  static SelectionTypeStringToEnum selectionTypeStringToEnumMap[] = {
11  { "All", All },
12  { "AllGlobalMuons", AllGlobalMuons },
13  { "AllStandAloneMuons", AllStandAloneMuons },
14  { "AllTrackerMuons", AllTrackerMuons },
15  { "TrackerMuonArbitrated", TrackerMuonArbitrated },
16  { "AllArbitrated", AllArbitrated },
17  { "GlobalMuonPromptTight", GlobalMuonPromptTight },
18  { "TMLastStationLoose", TMLastStationLoose },
19  { "TMLastStationTight", TMLastStationTight },
20  { "TM2DCompatibilityLoose", TM2DCompatibilityLoose },
21  { "TM2DCompatibilityTight", TM2DCompatibilityTight },
22  { "TMOneStationLoose", TMOneStationLoose },
23  { "TMOneStationTight", TMOneStationTight },
24  { "TMLastStationOptimizedLowPtLoose", TMLastStationOptimizedLowPtLoose },
25  { "TMLastStationOptimizedLowPtTight", TMLastStationOptimizedLowPtTight },
26  { "GMTkChiCompatibility", GMTkChiCompatibility },
27  { "GMStaChiCompatibility", GMStaChiCompatibility},
28  { "GMTkKinkTight", GMTkKinkTight},
29  { "TMLastStationAngLoose", TMLastStationAngLoose },
30  { "TMLastStationAngTight", TMLastStationAngTight },
31  { "TMOneStationAngLoose", TMOneStationAngLoose },
32  { "TMOneStationAngTight", TMOneStationAngTight },
33  { "TMLastStationOptimizedBarrelLowPtLoose", TMLastStationOptimizedBarrelLowPtLoose },
34  { "TMLastStationOptimizedBarrelLowPtTight", TMLastStationOptimizedBarrelLowPtTight },
35  { 0, (SelectionType)-1 }
36  };
37 
39  bool found = false;
40  for(int i = 0; selectionTypeStringToEnumMap[i].label && (! found); ++i)
41  if (! strcmp(label.c_str(), selectionTypeStringToEnumMap[i].label)) {
42  found = true;
43  value = selectionTypeStringToEnumMap[i].value;
44  }
45 
46  // in case of unrecognized selection type
47  if (! found) throw cms::Exception("MuonSelectorError") << label << " is not a recognized SelectionType";
48  return value;
49 }
50 }
51 
53  double maxChamberDist,
54  double maxChamberDistPull,
55  reco::Muon::ArbitrationType arbitrationType )
56 {
57  unsigned int theMask = 0;
58 
59  for(int stationIdx = 1; stationIdx < 5; ++stationIdx)
60  for(int detectorIdx = 1; detectorIdx < 3; ++detectorIdx)
61  if(muon.trackDist(stationIdx,detectorIdx,arbitrationType) < maxChamberDist &&
62  muon.trackDist(stationIdx,detectorIdx,arbitrationType)/muon.trackDistErr(stationIdx,detectorIdx,arbitrationType) < maxChamberDistPull)
63  theMask += 1<<((stationIdx-1)+4*(detectorIdx-1));
64 
65  return theMask;
66 }
67 
68 // ------------ method to calculate the calo compatibility for a track with matched muon info ------------
70  return muon.caloCompatibility();
71 }
72 
73 // ------------ method to calculate the segment compatibility for a track with matched muon info ------------
75  bool use_weight_regain_at_chamber_boundary = true;
76  bool use_match_dist_penalty = true;
77 
78  int nr_of_stations_crossed = 0;
79  int nr_of_stations_with_segment = 0;
80  std::vector<int> stations_w_track(8);
81  std::vector<int> station_has_segmentmatch(8);
82  std::vector<int> station_was_crossed(8);
83  std::vector<float> stations_w_track_at_boundary(8);
84  std::vector<float> station_weight(8);
85  int position_in_stations = 0;
86  float full_weight = 0.;
87 
88  for(int i = 1; i<=8; ++i) {
89  // ********************************************************;
90  // *** fill local info for this muon (do some counting) ***;
91  // ************** begin ***********************************;
92  if(i<=4) { // this is the section for the DTs
93  if( muon.trackDist(i,1,arbitrationType) < 999999 ) { //current "raw" info that a track is close to a chamber
94  ++nr_of_stations_crossed;
95  station_was_crossed[i-1] = 1;
96  if(muon.trackDist(i,1,arbitrationType) > -10. ) stations_w_track_at_boundary[i-1] = muon.trackDist(i,1,arbitrationType);
97  else stations_w_track_at_boundary[i-1] = 0.;
98  }
99  if( muon.segmentX(i,1,arbitrationType) < 999999 ) { //current "raw" info that a segment is matched to the current track
100  ++nr_of_stations_with_segment;
101  station_has_segmentmatch[i-1] = 1;
102  }
103  }
104  else { // this is the section for the CSCs
105  if( muon.trackDist(i-4,2,arbitrationType) < 999999 ) { //current "raw" info that a track is close to a chamber
106  ++nr_of_stations_crossed;
107  station_was_crossed[i-1] = 1;
108  if(muon.trackDist(i-4,2,arbitrationType) > -10. ) stations_w_track_at_boundary[i-1] = muon.trackDist(i-4,2,arbitrationType);
109  else stations_w_track_at_boundary[i-1] = 0.;
110  }
111  if( muon.segmentX(i-4,2,arbitrationType) < 999999 ) { //current "raw" info that a segment is matched to the current track
112  ++nr_of_stations_with_segment;
113  station_has_segmentmatch[i-1] = 1;
114  }
115  }
116  // rough estimation of chamber border efficiency (should be parametrized better, this is just a quick guess):
117  // TF1 * merf = new TF1("merf","-0.5*(TMath::Erf(x/6.)-1)",-100,100);
118  // use above value to "unpunish" missing segment if close to border, i.e. rather than not adding any weight, add
119  // the one from the function. Only for dist ~> -10 cm, else full punish!.
120 
121  // ********************************************************;
122  // *** fill local info for this muon (do some counting) ***;
123  // ************** end *************************************;
124  }
125 
126  // ********************************************************;
127  // *** calculate weights for each station *****************;
128  // ************** begin ***********************************;
129  // const float slope = 0.5;
130  // const float attenuate_weight_regain = 1.;
131  // if attenuate_weight_regain < 1., additional punishment if track is close to boundary and no segment
132  const float attenuate_weight_regain = 0.5;
133 
134  for(int i = 1; i<=8; ++i) { // loop over all possible stations
135 
136  // first set all weights if a station has been crossed
137  // later penalize if a station did not have a matching segment
138 
139  //old logic if(station_has_segmentmatch[i-1] > 0 ) { // the track has an associated segment at the current station
140  if( station_was_crossed[i-1] > 0 ) { // the track crossed this chamber (or was nearby)
141  // - Apply a weight depending on the "depth" of the muon passage.
142  // - The station_weight is later reduced for stations with badly matched segments.
143  // - Even if there is no segment but the track passes close to a chamber boundary, the
144  // weight is set non zero and can go up to 0.5 of the full weight if the track is quite
145  // far from any station.
146  ++position_in_stations;
147 
148  switch ( nr_of_stations_crossed ) { // define different weights depending on how many stations were crossed
149  case 1 :
150  station_weight[i-1] = 1.;
151  break;
152  case 2 :
153  if ( position_in_stations == 1 ) station_weight[i-1] = 0.33;
154  else station_weight[i-1] = 0.67;
155  break;
156  case 3 :
157  if ( position_in_stations == 1 ) station_weight[i-1] = 0.23;
158  else if( position_in_stations == 2 ) station_weight[i-1] = 0.33;
159  else station_weight[i-1] = 0.44;
160  break;
161  case 4 :
162  if ( position_in_stations == 1 ) station_weight[i-1] = 0.10;
163  else if( position_in_stations == 2 ) station_weight[i-1] = 0.20;
164  else if( position_in_stations == 3 ) station_weight[i-1] = 0.30;
165  else station_weight[i-1] = 0.40;
166  break;
167 
168  default :
169 // LogTrace("MuonIdentification")<<" // Message: A muon candidate track has more than 4 stations with matching segments.";
170 // LogTrace("MuonIdentification")<<" // Did not expect this - please let me know: ibloch@fnal.gov";
171  // for all other cases
172  station_weight[i-1] = 1./nr_of_stations_crossed;
173  }
174 
175  if( use_weight_regain_at_chamber_boundary ) { // reconstitute some weight if there is no match but the segment is close to a boundary:
176  if(station_has_segmentmatch[i-1] <= 0 && stations_w_track_at_boundary[i-1] != 0. ) {
177  // if segment is not present but track in inefficient region, do not count as "missing match" but add some reduced weight.
178  // original "match weight" is currently reduced by at least attenuate_weight_regain, variing with an error function down to 0 if the track is
179  // inside the chamber.
180  station_weight[i-1] = station_weight[i-1]*attenuate_weight_regain*0.5*(TMath::Erf(stations_w_track_at_boundary[i-1]/6.)+1.); // remark: the additional scale of 0.5 normalizes Err to run from 0 to 1 in y
181  }
182  else if(station_has_segmentmatch[i-1] <= 0 && stations_w_track_at_boundary[i-1] == 0.) { // no segment match and track well inside chamber
183  // full penalization
184  station_weight[i-1] = 0.;
185  }
186  }
187  else { // always fully penalize tracks with no matching segment, whether the segment is close to the boundary or not.
188  if(station_has_segmentmatch[i-1] <= 0) station_weight[i-1] = 0.;
189  }
190 
191  if( station_has_segmentmatch[i-1] > 0 && 42 == 42 ) { // if track has matching segment, but the matching is not high quality, penalize
192  if(i<=4) { // we are in the DTs
193  if( muon.dY(i,1,arbitrationType) < 999999 && muon.dX(i,1,arbitrationType) < 999999) { // have both X and Y match
194  if(
195  TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.))> 1. ) {
196  // reduce weight
197  if(use_match_dist_penalty) {
198  // only use pull if 3 sigma is not smaller than 3 cm
199  if(TMath::Sqrt(TMath::Power(muon.dX(i,1,arbitrationType),2.)+TMath::Power(muon.dY(i,1,arbitrationType),2.)) < 3. && TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.)) > 3. ) {
200  station_weight[i-1] *= 1./TMath::Power(
201  TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i,1,arbitrationType),2.)+TMath::Power(muon.dY(i,1,arbitrationType),2.)),(double)1.),.25);
202  }
203  else {
204  station_weight[i-1] *= 1./TMath::Power(
205  TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.)),.25);
206  }
207  }
208  }
209  }
210  else if (muon.dY(i,1,arbitrationType) >= 999999) { // has no match in Y
211  if( muon.pullX(i,1,arbitrationType) > 1. ) { // has a match in X. Pull larger that 1 to avoid increasing the weight (just penalize, don't anti-penalize)
212  // reduce weight
213  if(use_match_dist_penalty) {
214  // only use pull if 3 sigma is not smaller than 3 cm
215  if( muon.dX(i,1,arbitrationType) < 3. && muon.pullX(i,1,arbitrationType) > 3. ) {
216  station_weight[i-1] *= 1./TMath::Power(TMath::Max((double)muon.dX(i,1,arbitrationType),(double)1.),.25);
217  }
218  else {
219  station_weight[i-1] *= 1./TMath::Power(muon.pullX(i,1,arbitrationType),.25);
220  }
221  }
222  }
223  }
224  else { // has no match in X
225  if( muon.pullY(i,1,arbitrationType) > 1. ) { // has a match in Y. Pull larger that 1 to avoid increasing the weight (just penalize, don't anti-penalize)
226  // reduce weight
227  if(use_match_dist_penalty) {
228  // only use pull if 3 sigma is not smaller than 3 cm
229  if( muon.dY(i,1,arbitrationType) < 3. && muon.pullY(i,1,arbitrationType) > 3. ) {
230  station_weight[i-1] *= 1./TMath::Power(TMath::Max((double)muon.dY(i,1,arbitrationType),(double)1.),.25);
231  }
232  else {
233  station_weight[i-1] *= 1./TMath::Power(muon.pullY(i,1,arbitrationType),.25);
234  }
235  }
236  }
237  }
238  }
239  else { // We are in the CSCs
240  if(
241  TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)) > 1. ) {
242  // reduce weight
243  if(use_match_dist_penalty) {
244  // only use pull if 3 sigma is not smaller than 3 cm
245  if(TMath::Sqrt(TMath::Power(muon.dX(i-4,2,arbitrationType),2.)+TMath::Power(muon.dY(i-4,2,arbitrationType),2.)) < 3. && TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)) > 3. ) {
246  station_weight[i-1] *= 1./TMath::Power(
247  TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i-4,2,arbitrationType),2.)+TMath::Power(muon.dY(i-4,2,arbitrationType),2.)),(double)1.),.25);
248  }
249  else {
250  station_weight[i-1] *= 1./TMath::Power(
251  TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)),.25);
252  }
253  }
254  }
255  }
256  }
257 
258  // Thoughts:
259  // - should penalize if the segment has only x OR y info
260  // - should also use the segment direction, as it now works!
261 
262  }
263  else { // track did not pass a chamber in this station - just reset weight
264  station_weight[i-1] = 0.;
265  }
266 
267  //increment final weight for muon:
268  full_weight += station_weight[i-1];
269  }
270 
271  // if we don't expect any matches, we set the compatibility to
272  // 0.5 as the track is as compatible with a muon as it is with
273  // background - we should maybe rather set it to -0.5!
274  if( nr_of_stations_crossed == 0 ) {
275  // full_weight = attenuate_weight_regain*0.5;
276  full_weight = 0.5;
277  }
278 
279  // ********************************************************;
280  // *** calculate weights for each station *****************;
281  // ************** end *************************************;
282 
283  return full_weight;
284 
285 }
286 
289  double minCompatibility,
290  reco::Muon::ArbitrationType arbitrationType ) {
291  if (!muon.isMatchesValid()) return false;
292  bool goodMuon = false;
293 
294  switch( type ) {
295  case TM2DCompatibility:
296  // Simplistic first cut in the 2D segment- vs calo-compatibility plane. Will have to be refined!
297  if( ( (0.8*caloCompatibility( muon ))+(1.2*segmentCompatibility( muon, arbitrationType )) ) > minCompatibility ) goodMuon = true;
298  else goodMuon = false;
299  return goodMuon;
300  break;
301  default :
302  // LogTrace("MuonIdentification")<<" // Invalid Algorithm Type called!";
303  goodMuon = false;
304  return goodMuon;
305  break;
306  }
307 }
308 
311  int minNumberOfMatches,
312  double maxAbsDx,
313  double maxAbsPullX,
314  double maxAbsDy,
315  double maxAbsPullY,
316  double maxChamberDist,
317  double maxChamberDistPull,
318  reco::Muon::ArbitrationType arbitrationType,
319  bool syncMinNMatchesNRequiredStationsInBarrelOnly,
320  bool applyAlsoAngularCuts)
321 {
322  if (!muon.isMatchesValid()) return false;
323  bool goodMuon = false;
324 
325  if (type == TMLastStation) {
326  // To satisfy my own paranoia, if the user specifies that the
327  // minimum number of matches is zero, then return true.
328  if(minNumberOfMatches == 0) return true;
329 
330  unsigned int theStationMask = muon.stationMask(arbitrationType);
331  unsigned int theRequiredStationMask = RequiredStationMask(muon, maxChamberDist, maxChamberDistPull, arbitrationType);
332 
333  // Require that there be at least a minimum number of segments
334  int numSegs = 0;
335  int numRequiredStations = 0;
336  for(int it = 0; it < 8; ++it) {
337  if(theStationMask & 1<<it) ++numSegs;
338  if(theRequiredStationMask & 1<<it) ++numRequiredStations;
339  }
340 
341  // Make sure the minimum number of matches is not greater than
342  // the number of required stations but still greater than zero
343  if (syncMinNMatchesNRequiredStationsInBarrelOnly) {
344  // Note that we only do this in the barrel region!
345  if (fabs(muon.eta()) < 1.2) {
346  if(minNumberOfMatches > numRequiredStations)
347  minNumberOfMatches = numRequiredStations;
348  if(minNumberOfMatches < 1) //SK: this only happens for negative values
349  minNumberOfMatches = 1;
350  }
351  } else {
352  if(minNumberOfMatches > numRequiredStations)
353  minNumberOfMatches = numRequiredStations;
354  if(minNumberOfMatches < 1) //SK: this only happens for negative values
355  minNumberOfMatches = 1;
356  }
357 
358  if(numSegs >= minNumberOfMatches) goodMuon = 1;
359 
360  // Require that last required station have segment
361  // If there are zero required stations keep track
362  // of the last station with a segment so that we may
363  // apply the quality cuts below to it instead
364  int lastSegBit = 0;
365  if(theRequiredStationMask) {
366  for(int stationIdx = 7; stationIdx >= 0; --stationIdx)
367  if(theRequiredStationMask & 1<<stationIdx){
368  if(theStationMask & 1<<stationIdx) {
369  lastSegBit = stationIdx;
370  goodMuon &= 1;
371  break;
372  } else {
373  goodMuon = false;
374  break;
375  }
376  }
377  } else {
378  for(int stationIdx = 7; stationIdx >= 0; --stationIdx)
379  if(theStationMask & 1<<stationIdx) {
380  lastSegBit = stationIdx;
381  break;
382  }
383  }
384 
385  if(!goodMuon) return false;
386 
387  // Impose pull cuts on last segment
388  int station = 0, detector = 0;
389  station = lastSegBit < 4 ? lastSegBit+1 : lastSegBit-3;
390  detector = lastSegBit < 4 ? 1 : 2;
391 
392  // Check x information
393  if(fabs(muon.pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
394  fabs(muon.dX(station,detector,arbitrationType)) > maxAbsDx)
395  return false;
396 
397  if(applyAlsoAngularCuts && fabs(muon.pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX)
398  return false;
399 
400  // Is this a tight algorithm, i.e. do we bother to check y information?
401  if (maxAbsDy < 999999) { // really if maxAbsDy < 1E9 as currently defined
402 
403  // Check y information
404  if (detector == 2) { // CSC
405  if(fabs(muon.pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
406  fabs(muon.dY(station,2,arbitrationType)) > maxAbsDy)
407  return false;
408 
409  if(applyAlsoAngularCuts && fabs(muon.pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY)
410  return false;
411  } else {
412  //
413  // In DT, if this is a "Tight" algorithm and the last segment is
414  // missing y information (always the case in station 4!!!), impose
415  // respective cuts on the next station in the stationMask that has
416  // a segment with y information. If there are no segments with y
417  // information then there is nothing to penalize. Should we
418  // penalize in Tight for having zero segments with y information?
419  // That is the fundamental question. Of course I am being uber
420  // paranoid; if this is a good muon then there will probably be at
421  // least one segment with y information but not always. Suppose
422  // somehow a muon only creates segments in station 4, then we
423  // definitely do not want to require that there be at least one
424  // segment with y information because we will lose it completely.
425  //
426 
427  for (int stationIdx = station; stationIdx > 0; --stationIdx) {
428  if(! (theStationMask & 1<<(stationIdx-1))) // don't bother if the station is not in the stationMask
429  continue;
430 
431  if(muon.dY(stationIdx,1,arbitrationType) > 999998) // no y-information
432  continue;
433 
434  if(fabs(muon.pullY(stationIdx,1,arbitrationType,1)) > maxAbsPullY &&
435  fabs(muon.dY(stationIdx,1,arbitrationType)) > maxAbsDy) {
436  return false;
437  }
438 
439  if(applyAlsoAngularCuts && fabs(muon.pullDyDz(stationIdx,1,arbitrationType,1)) > maxAbsPullY)
440  return false;
441 
442  // If we get this far then great this is a good muon
443  return true;
444  }
445  }
446  }
447 
448  return goodMuon;
449  } // TMLastStation
450 
451  // TMOneStation requires only that there be one "good" segment, regardless
452  // of the required stations. We do not penalize if there are absolutely zero
453  // segments with y information in the Tight algorithm. Maybe I'm being
454  // paranoid but so be it. If it's really a good muon then we will probably
455  // find at least one segment with both x and y information but you never
456  // know, and I don't want to deal with a potential inefficiency in the DT
457  // like we did with the original TMLastStation. Incidentally, not penalizing
458  // for total lack of y information in the Tight algorithm is what is done in
459  // the new TMLastStation
460  //
461  if (type == TMOneStation) {
462  unsigned int theStationMask = muon.stationMask(arbitrationType);
463 
464  // Of course there must be at least one segment
465  if (! theStationMask) return false;
466 
467  int station = 0, detector = 0;
468  // Keep track of whether or not there is a DT segment with y information.
469  // In the end, if it turns out there are absolutely zero DT segments with
470  // y information, require only that there was a segment with good x info.
471  // This of course only applies to the Tight algorithms.
472  bool existsGoodDTSegX = false;
473  bool existsDTSegY = false;
474 
475  // Impose cuts on the segments in the station mask until we find a good one
476  // Might as well start with the lowest bit to speed things up.
477  for(int stationIdx = 0; stationIdx <= 7; ++stationIdx)
478  if(theStationMask & 1<<stationIdx) {
479  station = stationIdx < 4 ? stationIdx+1 : stationIdx-3;
480  detector = stationIdx < 4 ? 1 : 2;
481 
482  if((fabs(muon.pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
483  fabs(muon.dX(station,detector,arbitrationType)) > maxAbsDx) ||
484  (applyAlsoAngularCuts && fabs(muon.pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX))
485  continue;
486  else if (detector == 1)
487  existsGoodDTSegX = true;
488 
489  // Is this a tight algorithm? If yes, use y information
490  if (maxAbsDy < 999999) {
491  if (detector == 2) { // CSC
492  if((fabs(muon.pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
493  fabs(muon.dY(station,2,arbitrationType)) > maxAbsDy) ||
494  (applyAlsoAngularCuts && fabs(muon.pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY))
495  continue;
496  } else {
497 
498  if(muon.dY(station,1,arbitrationType) > 999998) // no y-information
499  continue;
500  else
501  existsDTSegY = true;
502 
503  if((fabs(muon.pullY(station,1,arbitrationType,1)) > maxAbsPullY &&
504  fabs(muon.dY(station,1,arbitrationType)) > maxAbsDy) ||
505  (applyAlsoAngularCuts && fabs(muon.pullDyDz(station,1,arbitrationType,1)) > maxAbsPullY)) {
506  continue;
507  }
508  }
509  }
510 
511  // If we get this far then great this is a good muon
512  return true;
513  }
514 
515  // If we get this far then for sure there are no "good" CSC segments. For
516  // DT, check if there were any segments with y information. If there
517  // were none, but there was a segment with good x, then we're happy. If
518  // there WERE segments with y information, then they must have been shit
519  // since we are here so fail it. Of course, if this is a Loose algorithm
520  // then fail immediately since if we had good x we would already have
521  // returned true
522  if (maxAbsDy < 999999) {
523  if (existsDTSegY)
524  return false;
525  else if (existsGoodDTSegX)
526  return true;
527  } else
528  return false;
529  } // TMOneStation
530 
531  return goodMuon;
532 }
533 
535  reco::Muon::ArbitrationType arbitrationType)
536 {
537  switch (type)
538  {
539  case muon::All:
540  return true;
541  break;
543  return muon.isGlobalMuon();
544  break;
546  return muon.isTrackerMuon();
547  break;
549  return muon.isStandAloneMuon();
550  break;
552  return muon.isTrackerMuon() && muon.numberOfMatches(arbitrationType)>0;
553  break;
554  case muon::AllArbitrated:
555  return ! muon.isTrackerMuon() || muon.numberOfMatches(arbitrationType)>0;
556  break;
558  return muon.isGlobalMuon() && muon.globalTrack()->normalizedChi2()<10. && muon.globalTrack()->hitPattern().numberOfValidMuonHits() >0;
559  break;
560  // For "Loose" algorithms we choose maximum y quantity cuts of 1E9 instead of
561  // 9999 as before. We do this because the muon methods return 999999 (note
562  // there are six 9's) when the requested information is not available. For
563  // example, if a muon fails to traverse the z measuring superlayer in a station
564  // in the DT, then all methods involving segmentY in this station return
565  // 999999 to demonstrate that the information is missing. In order to not
566  // penalize muons for missing y information in Loose algorithms where we do
567  // not care at all about y information, we raise these limits. In the
568  // TMLastStation and TMOneStation algorithms we actually use this huge number
569  // to determine whether to consider y information at all.
571  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,true,false);
572  break;
574  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,true,false);
575  break;
577  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
578  break;
580  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
581  break;
583  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
584  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
585  else
586  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,false,false);
587  break;
589  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
590  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
591  else
592  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,false,false);
593  break;
594  //compatibility loose
596  return muon.isTrackerMuon() && isGoodMuon(muon,TM2DCompatibility,0.7,arbitrationType);
597  break;
598  //compatibility tight
600  return muon.isTrackerMuon() && isGoodMuon(muon,TM2DCompatibility,1.0,arbitrationType);
601  break;
603  return muon.isGlobalMuon() && muon.isQualityValid() && fabs(muon.combinedQuality().trkRelChi2 - muon.innerTrack()->normalizedChi2()) < 2.0;
604  break;
606  return muon.isGlobalMuon() && muon.isQualityValid() && fabs(muon.combinedQuality().staRelChi2 - muon.outerTrack()->normalizedChi2()) < 2.0;
607  break;
608  case muon::GMTkKinkTight:
609  return muon.isGlobalMuon() && muon.isQualityValid() && muon.combinedQuality().trkKink < 100.0;
610  break;
612  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,false,true);
613  break;
615  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,false,true);
616  break;
618  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,true);
619  break;
621  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,true);
622  break;
624  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
625  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
626  else
627  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,true,false);
628  break;
630  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
631  return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
632  else
633  return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,true,false);
634  break;
635  default:
636  return false;
637  }
638 }
639 
640 bool muon::overlap( const reco::Muon& muon1, const reco::Muon& muon2,
641  double pullX, double pullY, bool checkAdjacentChambers)
642 {
643  unsigned int nMatches1 = muon1.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
644  unsigned int nMatches2 = muon2.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
645  unsigned int betterMuon = ( muon1.pt() > muon2.pt() ? 1 : 2 );
646  for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber1 = muon1.matches().begin();
647  chamber1 != muon1.matches().end(); ++chamber1 )
648  for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber2 = muon2.matches().begin();
649  chamber2 != muon2.matches().end(); ++chamber2 )
650  {
651 
652  // if ( (chamber1->segmentMatches.empty() || chamber2->segmentMatches.empty()) ) continue;
653 
654  // handle case where both muons have information about the same chamber
655  // here we know how close they are
656  if ( chamber1->id == chamber2->id ){
657  // found the same chamber
658  if ( fabs(chamber1->x-chamber2->x) <
659  pullX * sqrt(chamber1->xErr*chamber1->xErr+chamber2->xErr*chamber2->xErr) )
660  {
661  if ( betterMuon == 1 )
662  nMatches2--;
663  else
664  nMatches1--;
665  if ( nMatches1==0 || nMatches2==0 ) return true;
666  continue;
667  }
668  if ( fabs(chamber1->y-chamber2->y) <
669  pullY * sqrt(chamber1->yErr*chamber1->yErr+chamber2->yErr*chamber2->yErr) )
670  {
671  if ( betterMuon == 1 )
672  nMatches2--;
673  else
674  nMatches1--;
675  if ( nMatches1==0 || nMatches2==0 ) return true;
676  }
677  } else {
678  if ( ! checkAdjacentChambers ) continue;
679  // check if tracks are pointing into overlaping region of the CSC detector
680  if ( chamber1->id.subdetId() != MuonSubdetId::CSC ||
681  chamber2->id.subdetId() != MuonSubdetId::CSC ) continue;
682  CSCDetId id1(chamber1->id);
683  CSCDetId id2(chamber2->id);
684  if ( id1.endcap() != id2.endcap() ) continue;
685  if ( id1.station() != id2.station() ) continue;
686  if ( id1.ring() != id2.ring() ) continue;
687  if ( abs(id1.chamber() - id2.chamber())>1 ) continue;
688  // FIXME: we don't handle 18->1; 36->1 transitions since
689  // I don't know how to check for sure how many chambers
690  // are there. Probably need to hard code some checks.
691 
692  // Now we have to make sure that both tracks are close to an edge
693  // FIXME: ignored Y coordinate for now
694  if ( fabs(chamber1->edgeX) > chamber1->xErr*pullX ) continue;
695  if ( fabs(chamber2->edgeX) > chamber2->xErr*pullX ) continue;
696  if ( chamber1->x * chamber2->x < 0 ) { // check if the same edge
697  if ( betterMuon == 1 )
698  nMatches2--;
699  else
700  nMatches1--;
701  if ( nMatches1==0 || nMatches2==0 ) return true;
702  }
703  }
704  }
705  return false;
706 }
707 
708 
710 
711  if(!muon.isPFMuon() || !muon.isGlobalMuon() ) return false;
712 
713  bool muID = isGoodMuon(muon,GlobalMuonPromptTight) && (muon.numberOfMatchedStations() > 1);
714 
715 
716  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
717  muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
718 
719 
720  bool ip = fabs(muon.muonBestTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.muonBestTrack()->dz(vtx.position())) < 0.5;
721 
722  return muID && hits && ip;
723 }
724 
726  return muon.isPFMuon() && ( muon.isGlobalMuon() || muon.isTrackerMuon());
727 }
728 
729 bool muon::isSoftMuon(const reco::Muon& muon, const reco::Vertex& vtx){
730 
731  bool muID = muon::isGoodMuon(muon, TMOneStationTight);
732  if(!muID) return false;
733 
734  bool layers = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 && muon.innerTrack()->hitPattern().pixelLayersWithMeasurement() > 1;
735 
736  bool chi2 = muon.innerTrack()->normalizedChi2() < 1.8;
737 
738  bool ip = fabs(muon.innerTrack()->dxy(vtx.position())) < 3. && fabs(muon.innerTrack()->dz(vtx.position())) < 30.;
739 
740  return muID && layers && ip && chi2;
741 }
742 
743 bool muon::isHighPtMuon(const reco::Muon& muon, const reco::Vertex& vtx, reco::TunePType tunePType){
744 
745  bool muID = muon.isGlobalMuon() && muon.globalTrack()->hitPattern().numberOfValidMuonHits() >0 && (muon.numberOfMatchedStations() > 1);
746  if(!muID) return false;
747 
748  if(tunePType == reco::improvedTuneP){
749  // Get the optimized track
750  reco::TrackRef cktTrack = (muon::tevOptimized(muon, 200, 17., 40., 0.25)).first;
751 
752 
753  bool momQuality = cktTrack->ptError()/cktTrack->pt() < 0.3;
754 
755 
756  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 && muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
757 
758 
759 
760  bool ip = fabs(cktTrack->dxy(vtx.position())) < 0.2 && fabs(cktTrack->dz(vtx.position())) < 0.5;
761 
762 
763 
764  return muID && hits && momQuality && ip;}
765 
766  else if(tunePType == reco::defaultTuneP){
767  // Get the optimized track
768  reco::TrackRef cktTrack = (muon::tevOptimized(muon, 200, 4., 6., -1)).first;
769 
770  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 8 && muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
771 
772  bool ip = fabs(cktTrack->dxy(vtx.position())) < 0.2 && fabs(cktTrack->dz(vtx.position())) < 0.5;
773 
774  return muID && hits && ip;}
775 
776  else return false;
777 }
778 
779 
780 int muon::sharedSegments( const reco::Muon& mu, const reco::Muon& mu2, unsigned int segmentArbitrationMask ) {
781  int ret = 0;
782 
783  // Will do with a stupid double loop, since creating and filling a map is probably _more_ inefficient for a single lookup.
784  for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.matches().begin();
785  chamberMatch != mu.matches().end(); ++chamberMatch) {
786  if (chamberMatch->segmentMatches.empty()) continue;
787  for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch2 = mu2.matches().begin();
788  chamberMatch2 != mu2.matches().end(); ++chamberMatch2) {
789  if (chamberMatch2->segmentMatches.empty()) continue;
790  if (chamberMatch2->id() != chamberMatch->id()) continue;
791  for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
792  segmentMatch != chamberMatch->segmentMatches.end(); ++segmentMatch) {
793  if (!segmentMatch->isMask(segmentArbitrationMask)) continue;
794  for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch2 = chamberMatch2->segmentMatches.begin();
795  segmentMatch2 != chamberMatch2->segmentMatches.end(); ++segmentMatch2) {
796  if (!segmentMatch2->isMask(segmentArbitrationMask)) continue;
797  if ((segmentMatch->cscSegmentRef.isNonnull() && segmentMatch->cscSegmentRef == segmentMatch2->cscSegmentRef) ||
798  (segmentMatch-> dtSegmentRef.isNonnull() && segmentMatch-> dtSegmentRef == segmentMatch2-> dtSegmentRef) ) {
799  ++ret;
800  } // is the same
801  } // segment of mu2 in chamber
802  } // segment of mu1 in chamber
803  } // chamber of mu2
804  } // chamber of mu1
805 
806  return ret;
807 }
int chamber() const
Definition: CSCDetId.h:70
type
Definition: HCALResponse.h:22
int i
Definition: DBlmapReader.cc:9
AlgorithmType
Definition: MuonSelectors.h:64
float trackDistErr(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:707
virtual TrackRef innerTrack() const
Definition: Muon.h:49
bool isTrackerMuon() const
Definition: Muon.h:212
#define abs(x)
Definition: mlp_lapack.h:159
bool isGlobalMuon() const
Definition: Muon.h:211
float trkKink
value of the kink algorithm applied to the inner track stub
Definition: MuonQuality.h:11
bool isMatchesValid() const
Definition: Muon.h:139
float dY(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:344
bool isStandAloneMuon() const
Definition: Muon.h:213
float trkRelChi2
chi2 value for the inner track stub with respect to the global track
Definition: MuonQuality.h:15
const Point & position() const
position
Definition: Vertex.h:93
virtual double eta() const
momentum pseudorapidity
float caloCompatibility(const reco::Muon &muon)
bool isLooseMuon(const reco::Muon &)
SelectionType
Selector type.
Definition: MuonSelectors.h:19
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
bool isSoftMuon(const reco::Muon &, const reco::Vertex &)
int endcap() const
Definition: CSCDetId.h:95
static const int CSC
Definition: MuonSubdetId.h:15
virtual TrackRef muonBestTrack() const
Definition: Muon.h:64
ArbitrationType
define arbitration schemes
Definition: Muon.h:173
T sqrt(T t)
Definition: SSEVec.h:46
MuonQuality combinedQuality() const
get energy deposition information
Definition: Muon.h:122
float pullDxDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
Definition: Muon.cc:391
float segmentCompatibility(const reco::Muon &muon, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
const int mu
Definition: Constants.h:23
bool isQualityValid() const
Definition: Muon.h:120
bool first
Definition: L1TdeRCT.cc:94
unsigned int stationMask(ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:113
virtual TrackRef outerTrack() const
reference to Track reconstructed in the muon detector only
Definition: Muon.h:52
float pullY(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
Definition: Muon.cc:380
float segmentX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:412
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
Definition: Muon.cc:56
int ring() const
Definition: CSCDetId.h:77
float staRelChi2
chi2 value for the outer track stub with respect to the global track
Definition: MuonQuality.h:17
int sharedSegments(const reco::Muon &muon1, const reco::Muon &muon2, unsigned int segmentArbitrationMask=reco::MuonSegmentMatch::BestInChamberByDR)
SelectionType selectionTypeFromString(const std::string &label)
Definition: MuonSelectors.cc:8
virtual double pt() const
transverse momentum
a lightweight &quot;map&quot; for selection type string label and enum value
Definition: MuonSelectors.h:54
unsigned int RequiredStationMask(const reco::Muon &muon, double maxChamberDist, double maxChamberDistPull, reco::Muon::ArbitrationType arbitrationType)
std::vector< MuonChamberMatch > & matches()
get muon matching information
Definition: Muon.h:141
float trackDist(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:690
float dX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:336
int numberOfMatchedStations(ArbitrationType type=SegmentAndTrackArbitration) const
Definition: Muon.cc:100
bool isPFMuon() const
Definition: Muon.h:215
reco::Muon::MuonTrackTypePair tevOptimized(const reco::TrackRef &combinedTrack, const reco::TrackRef &trackerTrack, const reco::TrackRef &tpfmsTrack, const reco::TrackRef &pickyTrack, const double ptThreshold=200., const double tune1=4., const double tune2=6., const double dptcut=-1.)
Definition: MuonCocktails.cc:9
int station() const
Definition: CSCDetId.h:88
bool isHighPtMuon(const reco::Muon &, const reco::Vertex &, reco::TunePType)
bool isTightMuon(const reco::Muon &, const reco::Vertex &)
TunePType
Definition: MuonFwd.h:32
float caloCompatibility() const
Definition: Muon.h:151
float pullDyDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
Definition: Muon.cc:401
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:55
float pullX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
Definition: Muon.cc:370