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 
8 namespace muon {
10  static const 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  {"RPCMuLoose", RPCMuLoose},
36  {"AllME0Muons", AllME0Muons},
37  {"ME0MuonArbitrated", ME0MuonArbitrated},
38  {"AllGEMMuons", AllGEMMuons},
39  {"GEMMuonArbitrated", GEMMuonArbitrated},
40  {"TriggerIdLoose", TriggerIdLoose},
41  {nullptr, (SelectionType)-1}};
42 
44  bool found = false;
45  for (int i = 0; selectionTypeStringToEnumMap[i].label && (!found); ++i)
46  if (!strcmp(label.c_str(), selectionTypeStringToEnumMap[i].label)) {
47  found = true;
48  value = selectionTypeStringToEnumMap[i].value;
49  }
50 
51  // in case of unrecognized selection type
52  if (!found)
53  throw cms::Exception("MuonSelectorError") << label << " is not a recognized SelectionType";
54  return value;
55  }
56 } // namespace muon
57 
59  double maxChamberDist,
60  double maxChamberDistPull,
61  reco::Muon::ArbitrationType arbitrationType) {
62  unsigned int theMask = 0;
63 
64  for (int stationIdx = 1; stationIdx < 5; ++stationIdx)
65  for (int detectorIdx = 1; detectorIdx < 3; ++detectorIdx)
66  if (muon.trackDist(stationIdx, detectorIdx, arbitrationType) < maxChamberDist &&
67  muon.trackDist(stationIdx, detectorIdx, arbitrationType) /
68  muon.trackDistErr(stationIdx, detectorIdx, arbitrationType) <
69  maxChamberDistPull)
70  theMask += 1 << ((stationIdx - 1) + 4 * (detectorIdx - 1));
71 
72  return theMask;
73 }
74 
75 // ------------ method to calculate the calo compatibility for a track with matched muon info ------------
76 float muon::caloCompatibility(const reco::Muon& muon) { return muon.caloCompatibility(); }
77 
78 // ------------ method to calculate the segment compatibility for a track with matched muon info ------------
80  bool use_weight_regain_at_chamber_boundary = true;
81  bool use_match_dist_penalty = true;
82 
83  int nr_of_stations_crossed = 0;
84  int nr_of_stations_with_segment = 0;
85  std::vector<int> stations_w_track(8);
86  std::vector<int> station_has_segmentmatch(8);
87  std::vector<int> station_was_crossed(8);
88  std::vector<float> stations_w_track_at_boundary(8);
89  std::vector<float> station_weight(8);
90  int position_in_stations = 0;
91  float full_weight = 0.;
92 
93  for (int i = 1; i <= 8; ++i) {
94  // ********************************************************;
95  // *** fill local info for this muon (do some counting) ***;
96  // ************** begin ***********************************;
97  if (i <= 4) { // this is the section for the DTs
98  if (muon.trackDist(i, 1, arbitrationType) < 999999) { //current "raw" info that a track is close to a chamber
99  ++nr_of_stations_crossed;
100  station_was_crossed[i - 1] = 1;
101  if (muon.trackDist(i, 1, arbitrationType) > -10.)
102  stations_w_track_at_boundary[i - 1] = muon.trackDist(i, 1, arbitrationType);
103  else
104  stations_w_track_at_boundary[i - 1] = 0.;
105  }
106  //current "raw" info that a segment is matched to the current track
107  if (muon.segmentX(i, 1, arbitrationType) < 999999) {
108  ++nr_of_stations_with_segment;
109  station_has_segmentmatch[i - 1] = 1;
110  }
111  } else { // this is the section for the CSCs
112  if (muon.trackDist(i - 4, 2, arbitrationType) < 999999) { //current "raw" info that a track is close to a chamber
113  ++nr_of_stations_crossed;
114  station_was_crossed[i - 1] = 1;
115  if (muon.trackDist(i - 4, 2, arbitrationType) > -10.)
116  stations_w_track_at_boundary[i - 1] = muon.trackDist(i - 4, 2, arbitrationType);
117  else
118  stations_w_track_at_boundary[i - 1] = 0.;
119  }
120  //current "raw" info that a segment is matched to the current track
121  if (muon.segmentX(i - 4, 2, arbitrationType) < 999999) {
122  ++nr_of_stations_with_segment;
123  station_has_segmentmatch[i - 1] = 1;
124  }
125  }
126  // rough estimation of chamber border efficiency (should be parametrized better, this is just a quick guess):
127  // TF1 * merf = new TF1("merf","-0.5*(TMath::Erf(x/6.)-1)",-100,100);
128  // use above value to "unpunish" missing segment if close to border, i.e. rather than not adding any weight, add
129  // the one from the function. Only for dist ~> -10 cm, else full punish!.
130 
131  // ********************************************************;
132  // *** fill local info for this muon (do some counting) ***;
133  // ************** end *************************************;
134  }
135 
136  // ********************************************************;
137  // *** calculate weights for each station *****************;
138  // ************** begin ***********************************;
139  // const float slope = 0.5;
140  // const float attenuate_weight_regain = 1.;
141  // if attenuate_weight_regain < 1., additional punishment if track is close to boundary and no segment
142  const float attenuate_weight_regain = 0.5;
143 
144  for (int i = 1; i <= 8; ++i) { // loop over all possible stations
145 
146  // first set all weights if a station has been crossed
147  // later penalize if a station did not have a matching segment
148 
149  //old logic if(station_has_segmentmatch[i-1] > 0 ) { // the track has an associated segment at the current station
150  if (station_was_crossed[i - 1] > 0) { // the track crossed this chamber (or was nearby)
151  // - Apply a weight depending on the "depth" of the muon passage.
152  // - The station_weight is later reduced for stations with badly matched segments.
153  // - Even if there is no segment but the track passes close to a chamber boundary, the
154  // weight is set non zero and can go up to 0.5 of the full weight if the track is quite
155  // far from any station.
156  ++position_in_stations;
157 
158  switch (nr_of_stations_crossed) { // define different weights depending on how many stations were crossed
159  case 1:
160  station_weight[i - 1] = 1.;
161  break;
162  case 2:
163  if (position_in_stations == 1)
164  station_weight[i - 1] = 0.33;
165  else
166  station_weight[i - 1] = 0.67;
167  break;
168  case 3:
169  if (position_in_stations == 1)
170  station_weight[i - 1] = 0.23;
171  else if (position_in_stations == 2)
172  station_weight[i - 1] = 0.33;
173  else
174  station_weight[i - 1] = 0.44;
175  break;
176  case 4:
177  if (position_in_stations == 1)
178  station_weight[i - 1] = 0.10;
179  else if (position_in_stations == 2)
180  station_weight[i - 1] = 0.20;
181  else if (position_in_stations == 3)
182  station_weight[i - 1] = 0.30;
183  else
184  station_weight[i - 1] = 0.40;
185  break;
186 
187  default:
188  // LogTrace("MuonIdentification")<<" // Message: A muon candidate track has more than 4 stations with matching segments.";
189  // LogTrace("MuonIdentification")<<" // Did not expect this - please let me know: ibloch@fnal.gov";
190  // for all other cases
191  station_weight[i - 1] = 1. / nr_of_stations_crossed;
192  }
193 
194  if (use_weight_regain_at_chamber_boundary) { // reconstitute some weight if there is no match but the segment is close to a boundary:
195  if (station_has_segmentmatch[i - 1] <= 0 && stations_w_track_at_boundary[i - 1] != 0.) {
196  // if segment is not present but track in inefficient region, do not count as "missing match" but add some reduced weight.
197  // original "match weight" is currently reduced by at least attenuate_weight_regain, variing with an error function down to 0 if the track is
198  // inside the chamber.
199  // remark: the additional scale of 0.5 normalizes Err to run from 0 to 1 in y
200  station_weight[i - 1] = station_weight[i - 1] * attenuate_weight_regain * 0.5 *
201  (TMath::Erf(stations_w_track_at_boundary[i - 1] / 6.) + 1.);
202  } else if (station_has_segmentmatch[i - 1] <= 0 &&
203  stations_w_track_at_boundary[i - 1] == 0.) { // no segment match and track well inside chamber
204  // full penalization
205  station_weight[i - 1] = 0.;
206  }
207  } else { // always fully penalize tracks with no matching segment, whether the segment is close to the boundary or not.
208  if (station_has_segmentmatch[i - 1] <= 0)
209  station_weight[i - 1] = 0.;
210  }
211 
212  // if track has matching segment, but the matching is not high quality, penalize
213  if (station_has_segmentmatch[i - 1] > 0 && 42 == 42) {
214  if (i <= 4) { // we are in the DTs
215  if (muon.dY(i, 1, arbitrationType) < 999999 &&
216  muon.dX(i, 1, arbitrationType) < 999999) { // have both X and Y match
217  if (TMath::Sqrt(TMath::Power(muon.pullX(i, 1, arbitrationType), 2.) +
218  TMath::Power(muon.pullY(i, 1, arbitrationType), 2.)) > 1.) {
219  // reduce weight
220  if (use_match_dist_penalty) {
221  // only use pull if 3 sigma is not smaller than 3 cm
222  if (TMath::Sqrt(TMath::Power(muon.dX(i, 1, arbitrationType), 2.) +
223  TMath::Power(muon.dY(i, 1, arbitrationType), 2.)) < 3. &&
224  TMath::Sqrt(TMath::Power(muon.pullX(i, 1, arbitrationType), 2.) +
225  TMath::Power(muon.pullY(i, 1, arbitrationType), 2.)) > 3.) {
226  station_weight[i - 1] *=
227  1. /
228  TMath::Power(TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i, 1, arbitrationType), 2.) +
229  TMath::Power(muon.dY(i, 1, arbitrationType), 2.)),
230  (double)1.),
231  .25);
232  } else {
233  station_weight[i - 1] *=
234  1. / TMath::Power(TMath::Sqrt(TMath::Power(muon.pullX(i, 1, arbitrationType), 2.) +
235  TMath::Power(muon.pullY(i, 1, arbitrationType), 2.)),
236  .25);
237  }
238  }
239  }
240  } else if (muon.dY(i, 1, arbitrationType) >= 999999) { // has no match in Y
241  // has a match in X. Pull larger that 1 to avoid increasing the weight (just penalize, don't anti-penalize)
242  if (muon.pullX(i, 1, arbitrationType) > 1.) {
243  // reduce weight
244  if (use_match_dist_penalty) {
245  // only use pull if 3 sigma is not smaller than 3 cm
246  if (muon.dX(i, 1, arbitrationType) < 3. && muon.pullX(i, 1, arbitrationType) > 3.) {
247  station_weight[i - 1] *=
248  1. / TMath::Power(TMath::Max((double)muon.dX(i, 1, arbitrationType), (double)1.), .25);
249  } else {
250  station_weight[i - 1] *= 1. / TMath::Power(muon.pullX(i, 1, arbitrationType), .25);
251  }
252  }
253  }
254  } else { // has no match in X
255  // has a match in Y. Pull larger that 1 to avoid increasing the weight (just penalize, don't anti-penalize)
256  if (muon.pullY(i, 1, arbitrationType) > 1.) {
257  // reduce weight
258  if (use_match_dist_penalty) {
259  // only use pull if 3 sigma is not smaller than 3 cm
260  if (muon.dY(i, 1, arbitrationType) < 3. && muon.pullY(i, 1, arbitrationType) > 3.) {
261  station_weight[i - 1] *=
262  1. / TMath::Power(TMath::Max((double)muon.dY(i, 1, arbitrationType), (double)1.), .25);
263  } else {
264  station_weight[i - 1] *= 1. / TMath::Power(muon.pullY(i, 1, arbitrationType), .25);
265  }
266  }
267  }
268  }
269  } else { // We are in the CSCs
270  if (TMath::Sqrt(TMath::Power(muon.pullX(i - 4, 2, arbitrationType), 2.) +
271  TMath::Power(muon.pullY(i - 4, 2, arbitrationType), 2.)) > 1.) {
272  // reduce weight
273  if (use_match_dist_penalty) {
274  // only use pull if 3 sigma is not smaller than 3 cm
275  if (TMath::Sqrt(TMath::Power(muon.dX(i - 4, 2, arbitrationType), 2.) +
276  TMath::Power(muon.dY(i - 4, 2, arbitrationType), 2.)) < 3. &&
277  TMath::Sqrt(TMath::Power(muon.pullX(i - 4, 2, arbitrationType), 2.) +
278  TMath::Power(muon.pullY(i - 4, 2, arbitrationType), 2.)) > 3.) {
279  station_weight[i - 1] *=
280  1. /
281  TMath::Power(TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i - 4, 2, arbitrationType), 2.) +
282  TMath::Power(muon.dY(i - 4, 2, arbitrationType), 2.)),
283  (double)1.),
284  .25);
285  } else {
286  station_weight[i - 1] *=
287  1. / TMath::Power(TMath::Sqrt(TMath::Power(muon.pullX(i - 4, 2, arbitrationType), 2.) +
288  TMath::Power(muon.pullY(i - 4, 2, arbitrationType), 2.)),
289  .25);
290  }
291  }
292  }
293  }
294  }
295 
296  // Thoughts:
297  // - should penalize if the segment has only x OR y info
298  // - should also use the segment direction, as it now works!
299 
300  } else { // track did not pass a chamber in this station - just reset weight
301  station_weight[i - 1] = 0.;
302  }
303 
304  //increment final weight for muon:
305  full_weight += station_weight[i - 1];
306  }
307 
308  // if we don't expect any matches, we set the compatibility to
309  // 0.5 as the track is as compatible with a muon as it is with
310  // background - we should maybe rather set it to -0.5!
311  if (nr_of_stations_crossed == 0) {
312  // full_weight = attenuate_weight_regain*0.5;
313  full_weight = 0.5;
314  }
315 
316  // ********************************************************;
317  // *** calculate weights for each station *****************;
318  // ************** end *************************************;
319 
320  return full_weight;
321 }
322 
325  double minCompatibility,
326  reco::Muon::ArbitrationType arbitrationType) {
327  if (!muon.isMatchesValid())
328  return false;
329  bool goodMuon = false;
330 
331  switch (type) {
332  case TM2DCompatibility:
333  // Simplistic first cut in the 2D segment- vs calo-compatibility plane. Will have to be refined!
334  if (((0.8 * caloCompatibility(muon)) + (1.2 * segmentCompatibility(muon, arbitrationType))) > minCompatibility)
335  goodMuon = true;
336  else
337  goodMuon = false;
338  return goodMuon;
339  break;
340  default:
341  // LogTrace("MuonIdentification")<<" // Invalid Algorithm Type called!";
342  goodMuon = false;
343  return goodMuon;
344  break;
345  }
346 }
347 
350  int minNumberOfMatches,
351  double maxAbsDx,
352  double maxAbsPullX,
353  double maxAbsDy,
354  double maxAbsPullY,
355  double maxChamberDist,
356  double maxChamberDistPull,
357  reco::Muon::ArbitrationType arbitrationType,
358  bool syncMinNMatchesNRequiredStationsInBarrelOnly,
359  bool applyAlsoAngularCuts) {
360  if (!muon.isMatchesValid())
361  return false;
362  bool goodMuon = false;
363 
364  if (type == TMLastStation) {
365  // To satisfy my own paranoia, if the user specifies that the
366  // minimum number of matches is zero, then return true.
367  if (minNumberOfMatches == 0)
368  return true;
369 
370  unsigned int theStationMask = muon.stationMask(arbitrationType);
371  unsigned int theRequiredStationMask =
372  RequiredStationMask(muon, maxChamberDist, maxChamberDistPull, arbitrationType);
373 
374  // Require that there be at least a minimum number of segments
375  int numSegs = 0;
376  int numRequiredStations = 0;
377  for (int it = 0; it < 8; ++it) {
378  if (theStationMask & 1 << it)
379  ++numSegs;
380  if (theRequiredStationMask & 1 << it)
381  ++numRequiredStations;
382  }
383 
384  // Make sure the minimum number of matches is not greater than
385  // the number of required stations but still greater than zero
386  if (syncMinNMatchesNRequiredStationsInBarrelOnly) {
387  // Note that we only do this in the barrel region!
388  if (fabs(muon.eta()) < 1.2) {
389  if (minNumberOfMatches > numRequiredStations)
390  minNumberOfMatches = numRequiredStations;
391  if (minNumberOfMatches < 1) //SK: this only happens for negative values
392  minNumberOfMatches = 1;
393  }
394  } else {
395  if (minNumberOfMatches > numRequiredStations)
396  minNumberOfMatches = numRequiredStations;
397  if (minNumberOfMatches < 1) //SK: this only happens for negative values
398  minNumberOfMatches = 1;
399  }
400 
401  if (numSegs >= minNumberOfMatches)
402  goodMuon = true;
403 
404  // Require that last required station have segment
405  // If there are zero required stations keep track
406  // of the last station with a segment so that we may
407  // apply the quality cuts below to it instead
408  int lastSegBit = 0;
409  if (theRequiredStationMask) {
410  for (int stationIdx = 7; stationIdx >= 0; --stationIdx)
411  if (theRequiredStationMask & 1 << stationIdx) {
412  if (theStationMask & 1 << stationIdx) {
413  lastSegBit = stationIdx;
414  goodMuon &= 1;
415  break;
416  } else {
417  goodMuon = false;
418  break;
419  }
420  }
421  } else {
422  for (int stationIdx = 7; stationIdx >= 0; --stationIdx)
423  if (theStationMask & 1 << stationIdx) {
424  lastSegBit = stationIdx;
425  break;
426  }
427  }
428 
429  if (!goodMuon)
430  return false;
431 
432  // Impose pull cuts on last segment
433  int station = 0, detector = 0;
434  station = lastSegBit < 4 ? lastSegBit + 1 : lastSegBit - 3;
435  detector = lastSegBit < 4 ? 1 : 2;
436 
437  // Check x information
438  if (fabs(muon.pullX(station, detector, arbitrationType, true)) > maxAbsPullX &&
439  fabs(muon.dX(station, detector, arbitrationType)) > maxAbsDx)
440  return false;
441 
442  if (applyAlsoAngularCuts && fabs(muon.pullDxDz(station, detector, arbitrationType, true)) > maxAbsPullX)
443  return false;
444 
445  // Is this a tight algorithm, i.e. do we bother to check y information?
446  if (maxAbsDy < 999999) { // really if maxAbsDy < 1E9 as currently defined
447 
448  // Check y information
449  if (detector == 2) { // CSC
450  if (fabs(muon.pullY(station, 2, arbitrationType, true)) > maxAbsPullY &&
451  fabs(muon.dY(station, 2, arbitrationType)) > maxAbsDy)
452  return false;
453 
454  if (applyAlsoAngularCuts && fabs(muon.pullDyDz(station, 2, arbitrationType, true)) > maxAbsPullY)
455  return false;
456  } else {
457  //
458  // In DT, if this is a "Tight" algorithm and the last segment is
459  // missing y information (always the case in station 4!!!), impose
460  // respective cuts on the next station in the stationMask that has
461  // a segment with y information. If there are no segments with y
462  // information then there is nothing to penalize. Should we
463  // penalize in Tight for having zero segments with y information?
464  // That is the fundamental question. Of course I am being uber
465  // paranoid; if this is a good muon then there will probably be at
466  // least one segment with y information but not always. Suppose
467  // somehow a muon only creates segments in station 4, then we
468  // definitely do not want to require that there be at least one
469  // segment with y information because we will lose it completely.
470  //
471 
472  for (int stationIdx = station; stationIdx > 0; --stationIdx) {
473  if (!(theStationMask & 1 << (stationIdx - 1))) // don't bother if the station is not in the stationMask
474  continue;
475 
476  if (muon.dY(stationIdx, 1, arbitrationType) > 999998) // no y-information
477  continue;
478 
479  if (fabs(muon.pullY(stationIdx, 1, arbitrationType, true)) > maxAbsPullY &&
480  fabs(muon.dY(stationIdx, 1, arbitrationType)) > maxAbsDy) {
481  return false;
482  }
483 
484  if (applyAlsoAngularCuts && fabs(muon.pullDyDz(stationIdx, 1, arbitrationType, true)) > maxAbsPullY)
485  return false;
486 
487  // If we get this far then great this is a good muon
488  return true;
489  }
490  }
491  }
492 
493  return goodMuon;
494  } // TMLastStation
495 
496  // TMOneStation requires only that there be one "good" segment, regardless
497  // of the required stations. We do not penalize if there are absolutely zero
498  // segments with y information in the Tight algorithm. Maybe I'm being
499  // paranoid but so be it. If it's really a good muon then we will probably
500  // find at least one segment with both x and y information but you never
501  // know, and I don't want to deal with a potential inefficiency in the DT
502  // like we did with the original TMLastStation. Incidentally, not penalizing
503  // for total lack of y information in the Tight algorithm is what is done in
504  // the new TMLastStation
505  //
506  if (type == TMOneStation) {
507  unsigned int theStationMask = muon.stationMask(arbitrationType);
508 
509  // Of course there must be at least one segment
510  if (!theStationMask)
511  return false;
512 
513  int station = 0, detector = 0;
514  // Keep track of whether or not there is a DT segment with y information.
515  // In the end, if it turns out there are absolutely zero DT segments with
516  // y information, require only that there was a segment with good x info.
517  // This of course only applies to the Tight algorithms.
518  bool existsGoodDTSegX = false;
519  bool existsDTSegY = false;
520 
521  // Impose cuts on the segments in the station mask until we find a good one
522  // Might as well start with the lowest bit to speed things up.
523  for (int stationIdx = 0; stationIdx <= 7; ++stationIdx)
524  if (theStationMask & 1 << stationIdx) {
525  station = stationIdx < 4 ? stationIdx + 1 : stationIdx - 3;
526  detector = stationIdx < 4 ? 1 : 2;
527 
528  if ((fabs(muon.pullX(station, detector, arbitrationType, true)) > maxAbsPullX &&
529  fabs(muon.dX(station, detector, arbitrationType)) > maxAbsDx) ||
530  (applyAlsoAngularCuts && fabs(muon.pullDxDz(station, detector, arbitrationType, true)) > maxAbsPullX))
531  continue;
532  else if (detector == 1)
533  existsGoodDTSegX = true;
534 
535  // Is this a tight algorithm? If yes, use y information
536  if (maxAbsDy < 999999) {
537  if (detector == 2) { // CSC
538  if ((fabs(muon.pullY(station, 2, arbitrationType, true)) > maxAbsPullY &&
539  fabs(muon.dY(station, 2, arbitrationType)) > maxAbsDy) ||
540  (applyAlsoAngularCuts && fabs(muon.pullDyDz(station, 2, arbitrationType, true)) > maxAbsPullY))
541  continue;
542  } else {
543  if (muon.dY(station, 1, arbitrationType) > 999998) // no y-information
544  continue;
545  else
546  existsDTSegY = true;
547 
548  if ((fabs(muon.pullY(station, 1, arbitrationType, true)) > maxAbsPullY &&
549  fabs(muon.dY(station, 1, arbitrationType)) > maxAbsDy) ||
550  (applyAlsoAngularCuts && fabs(muon.pullDyDz(station, 1, arbitrationType, true)) > maxAbsPullY)) {
551  continue;
552  }
553  }
554  }
555 
556  // If we get this far then great this is a good muon
557  return true;
558  }
559 
560  // If we get this far then for sure there are no "good" CSC segments. For
561  // DT, check if there were any segments with y information. If there
562  // were none, but there was a segment with good x, then we're happy. If
563  // there WERE segments with y information, then they must have been shit
564  // since we are here so fail it. Of course, if this is a Loose algorithm
565  // then fail immediately since if we had good x we would already have
566  // returned true
567  if (maxAbsDy < 999999) {
568  if (existsDTSegY)
569  return false;
570  else if (existsGoodDTSegX)
571  return true;
572  } else
573  return false;
574  } // TMOneStation
575 
576  if (type == RPCMu) {
577  if (minNumberOfMatches == 0)
578  return true;
579 
580  int nMatch = 0;
581  for (std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = muon.matches().begin();
582  chamberMatch != muon.matches().end();
583  ++chamberMatch) {
584  if (chamberMatch->detector() != 3)
585  continue;
586 
587  const double trkX = chamberMatch->x;
588  const double errX = chamberMatch->xErr;
589 
590  for (std::vector<reco::MuonRPCHitMatch>::const_iterator rpcMatch = chamberMatch->rpcMatches.begin();
591  rpcMatch != chamberMatch->rpcMatches.end();
592  ++rpcMatch) {
593  const double rpcX = rpcMatch->x;
594 
595  const double dX = fabs(rpcX - trkX);
596  if (dX < maxAbsDx or dX / errX < maxAbsPullX) {
597  ++nMatch;
598  break;
599  }
600  }
601  }
602 
603  if (nMatch >= minNumberOfMatches)
604  return true;
605  else
606  return false;
607  } // RPCMu
608 
609  if (type == ME0Mu) {
610  if (minNumberOfMatches == 0)
611  return true;
612 
613  int nMatch = 0;
614  for (const auto& chamberMatch : muon.matches()) {
615  if (chamberMatch.detector() != MuonSubdetId::ME0)
616  continue;
617 
618  const double trkX = chamberMatch.x;
619  const double errX = chamberMatch.xErr;
620  const double trkY = chamberMatch.y;
621  const double errY = chamberMatch.yErr;
622 
623  for (const auto& segment : chamberMatch.me0Matches) {
624  const double me0X = segment.x;
625  const double me0ErrX = segment.xErr;
626  const double me0Y = segment.y;
627  const double me0ErrY = segment.yErr;
628 
629  const double dX = fabs(me0X - trkX);
630  const double dY = fabs(me0Y - trkY);
631  const double pullX = dX / std::sqrt(errX + me0ErrX);
632  const double pullY = dY / std::sqrt(errY + me0ErrY);
633 
634  if ((dX < maxAbsDx or pullX < maxAbsPullX) and (dY < maxAbsDy or pullY < maxAbsPullY)) {
635  ++nMatch;
636  break;
637  }
638  }
639  }
640 
641  return (nMatch >= minNumberOfMatches);
642  } // ME0Mu
643 
644  if (type == GEMMu) {
645  if (minNumberOfMatches == 0)
646  return true;
647 
648  int nMatch = 0;
649  for (const auto& chamberMatch : muon.matches()) {
650  if (chamberMatch.detector() != MuonSubdetId::GEM)
651  continue;
652 
653  const double trkX = chamberMatch.x;
654  const double errX = chamberMatch.xErr;
655  const double trkY = chamberMatch.y;
656  const double errY = chamberMatch.yErr;
657 
658  for (const auto& segment : chamberMatch.gemMatches) {
659  const double gemX = segment.x;
660  const double gemErrX = segment.xErr;
661  const double gemY = segment.y;
662  const double gemErrY = segment.yErr;
663 
664  const double dX = fabs(gemX - trkX);
665  const double dY = fabs(gemY - trkY);
666  const double pullX = dX / std::sqrt(errX + gemErrX);
667  const double pullY = dY / std::sqrt(errY + gemErrY);
668 
669  if ((dX < maxAbsDx or pullX < maxAbsPullX) and (dY < maxAbsDy or pullY < maxAbsPullY)) {
670  ++nMatch;
671  break;
672  }
673  }
674  }
675 
676  return (nMatch >= minNumberOfMatches);
677  } // GEMMu
678 
679  return goodMuon;
680 }
681 
683  switch (type) {
684  case muon::All:
685  return true;
686  break;
688  return muon.isGlobalMuon();
689  break;
691  return muon.isTrackerMuon();
692  break;
694  return muon.isStandAloneMuon();
695  break;
697  return muon.isTrackerMuon() && muon.numberOfMatches(arbitrationType) > 0;
698  break;
699  case muon::AllArbitrated:
700  return !muon.isTrackerMuon() || muon.numberOfMatches(arbitrationType) > 0;
701  break;
703  return muon.isGlobalMuon() && muon.globalTrack()->normalizedChi2() < 10. &&
704  muon.globalTrack()->hitPattern().numberOfValidMuonHits() > 0;
705  break;
706  // For "Loose" algorithms we choose maximum y quantity cuts of 1E9 instead of
707  // 9999 as before. We do this because the muon methods return 999999 (note
708  // there are six 9's) when the requested information is not available. For
709  // example, if a muon fails to traverse the z measuring superlayer in a station
710  // in the DT, then all methods involving segmentY in this station return
711  // 999999 to demonstrate that the information is missing. In order to not
712  // penalize muons for missing y information in Loose algorithms where we do
713  // not care at all about y information, we raise these limits. In the
714  // TMLastStation and TMOneStation algorithms we actually use this huge number
715  // to determine whether to consider y information at all.
717  return muon.isTrackerMuon() &&
718  isGoodMuon(muon, TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType, true, false);
719  break;
721  return muon.isTrackerMuon() &&
722  isGoodMuon(muon, TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType, true, false);
723  break;
725  return muon.isTrackerMuon() &&
726  isGoodMuon(muon, TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType, false, false);
727  break;
729  return muon.isTrackerMuon() &&
730  isGoodMuon(muon, TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType, false, false);
731  break;
733  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
734  return muon.isTrackerMuon() &&
735  isGoodMuon(muon, TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType, false, false);
736  else
737  return muon.isTrackerMuon() &&
738  isGoodMuon(muon, TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType, false, false);
739  break;
741  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
742  return muon.isTrackerMuon() &&
743  isGoodMuon(muon, TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType, false, false);
744  else
745  return muon.isTrackerMuon() &&
746  isGoodMuon(muon, TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType, false, false);
747  break;
748  //compatibility loose
750  return muon.isTrackerMuon() && isGoodMuon(muon, TM2DCompatibility, 0.7, arbitrationType);
751  break;
752  //compatibility tight
754  return muon.isTrackerMuon() && isGoodMuon(muon, TM2DCompatibility, 1.0, arbitrationType);
755  break;
757  return muon.isGlobalMuon() && muon.isQualityValid() &&
758  fabs(muon.combinedQuality().trkRelChi2 - muon.innerTrack()->normalizedChi2()) < 2.0;
759  break;
761  return muon.isGlobalMuon() && muon.isQualityValid() &&
762  fabs(muon.combinedQuality().staRelChi2 - muon.outerTrack()->normalizedChi2()) < 2.0;
763  break;
764  case muon::GMTkKinkTight:
765  return muon.isGlobalMuon() && muon.isQualityValid() && muon.combinedQuality().trkKink < 100.0;
766  break;
768  return muon.isTrackerMuon() &&
769  isGoodMuon(muon, TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType, false, true);
770  break;
772  return muon.isTrackerMuon() &&
773  isGoodMuon(muon, TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType, false, true);
774  break;
776  return muon.isTrackerMuon() &&
777  isGoodMuon(muon, TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType, false, true);
778  break;
780  return muon.isTrackerMuon() &&
781  isGoodMuon(muon, TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType, false, true);
782  break;
784  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
785  return muon.isTrackerMuon() &&
786  isGoodMuon(muon, TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType, false, false);
787  else
788  return muon.isTrackerMuon() &&
789  isGoodMuon(muon, TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType, true, false);
790  break;
792  if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
793  return muon.isTrackerMuon() &&
794  isGoodMuon(muon, TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType, false, false);
795  else
796  return muon.isTrackerMuon() &&
797  isGoodMuon(muon, TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType, true, false);
798  break;
799  case muon::RPCMuLoose:
800  return muon.isRPCMuon() && isGoodMuon(muon, RPCMu, 2, 20, 4, 1e9, 1e9, 1e9, 1e9, arbitrationType, false, false);
801  break;
802  case muon::AllME0Muons:
803  return muon.isME0Muon();
804  break;
806  return muon.isME0Muon() &&
807  isGoodMuon(muon, ME0Mu, 1, 1e9, 1e9, 1e9, 1e9, 1e9, 1e9, arbitrationType, false, false);
808  break;
809  case muon::AllGEMMuons:
810  return muon.isGEMMuon();
811  break;
813  return muon.isGEMMuon() &&
814  isGoodMuon(muon, GEMMu, 1, 1e9, 1e9, 1e9, 1e9, 1e9, 1e9, arbitrationType, false, false);
815  break;
817  return isLooseTriggerMuon(muon);
818  break;
819  default:
820  return false;
821  }
822 }
823 
825  const reco::Muon& muon1, const reco::Muon& muon2, double pullX, double pullY, bool checkAdjacentChambers) {
826  unsigned int nMatches1 = muon1.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
827  unsigned int nMatches2 = muon2.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
828  unsigned int betterMuon = (muon1.pt() > muon2.pt() ? 1 : 2);
829  for (std::vector<reco::MuonChamberMatch>::const_iterator chamber1 = muon1.matches().begin();
830  chamber1 != muon1.matches().end();
831  ++chamber1)
832  for (std::vector<reco::MuonChamberMatch>::const_iterator chamber2 = muon2.matches().begin();
833  chamber2 != muon2.matches().end();
834  ++chamber2) {
835  // if ( (chamber1->segmentMatches.empty() || chamber2->segmentMatches.empty()) ) continue;
836 
837  // handle case where both muons have information about the same chamber
838  // here we know how close they are
839  if (chamber1->id == chamber2->id) {
840  // found the same chamber
841  if (fabs(chamber1->x - chamber2->x) <
842  pullX * sqrt(chamber1->xErr * chamber1->xErr + chamber2->xErr * chamber2->xErr)) {
843  if (betterMuon == 1)
844  nMatches2--;
845  else
846  nMatches1--;
847  if (nMatches1 == 0 || nMatches2 == 0)
848  return true;
849  continue;
850  }
851  if (fabs(chamber1->y - chamber2->y) <
852  pullY * sqrt(chamber1->yErr * chamber1->yErr + chamber2->yErr * chamber2->yErr)) {
853  if (betterMuon == 1)
854  nMatches2--;
855  else
856  nMatches1--;
857  if (nMatches1 == 0 || nMatches2 == 0)
858  return true;
859  }
860  } else {
861  if (!checkAdjacentChambers)
862  continue;
863  // check if tracks are pointing into overlaping region of the CSC detector
864  if (chamber1->id.subdetId() != MuonSubdetId::CSC || chamber2->id.subdetId() != MuonSubdetId::CSC)
865  continue;
866  CSCDetId id1(chamber1->id);
867  CSCDetId id2(chamber2->id);
868  if (id1.endcap() != id2.endcap())
869  continue;
870  if (id1.station() != id2.station())
871  continue;
872  if (id1.ring() != id2.ring())
873  continue;
874  if (abs(id1.chamber() - id2.chamber()) > 1)
875  continue;
876  // FIXME: we don't handle 18->1; 36->1 transitions since
877  // I don't know how to check for sure how many chambers
878  // are there. Probably need to hard code some checks.
879 
880  // Now we have to make sure that both tracks are close to an edge
881  // FIXME: ignored Y coordinate for now
882  if (fabs(chamber1->edgeX) > chamber1->xErr * pullX)
883  continue;
884  if (fabs(chamber2->edgeX) > chamber2->xErr * pullX)
885  continue;
886  if (chamber1->x * chamber2->x < 0) { // check if the same edge
887  if (betterMuon == 1)
888  nMatches2--;
889  else
890  nMatches1--;
891  if (nMatches1 == 0 || nMatches2 == 0)
892  return true;
893  }
894  }
895  }
896  return false;
897 }
898 
900  // Requirements:
901  // - no depencence on information not availabe in the muon object
902  // - use only robust inputs
903  bool tk_id = muon::isGoodMuon(muon, TMOneStationTight);
904  if (not tk_id)
905  return false;
906  bool layer_requirements = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
907  muon.innerTrack()->hitPattern().pixelLayersWithMeasurement() > 0;
908  bool match_requirements =
909  (muon.expectedNnumberOfMatchedStations() < 2) or (muon.numberOfMatchedStations() > 1) or (muon.pt() < 8);
910  return layer_requirements and match_requirements;
911 }
912 
914  if (!muon.isPFMuon() || !muon.isGlobalMuon())
915  return false;
916 
917  bool muID = isGoodMuon(muon, GlobalMuonPromptTight) && (muon.numberOfMatchedStations() > 1);
918 
919  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
920  muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
921 
922  bool ip =
923  fabs(muon.muonBestTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.muonBestTrack()->dz(vtx.position())) < 0.5;
924 
925  return muID && hits && ip;
926 }
927 
929  return muon.isPFMuon() && (muon.isGlobalMuon() || muon.isTrackerMuon());
930 }
931 
932 bool muon::isMediumMuon(const reco::Muon& muon, bool run2016_hip_mitigation) {
933  if (not isLooseMuon(muon))
934  return false;
935  if (run2016_hip_mitigation) {
936  if (muon.innerTrack()->validFraction() < 0.49)
937  return false;
938  } else {
939  if (muon.innerTrack()->validFraction() < 0.8)
940  return false;
941  }
942 
943  bool goodGlb = muon.isGlobalMuon() && muon.globalTrack()->normalizedChi2() < 3. &&
944  muon.combinedQuality().chi2LocalPosition < 12. && muon.combinedQuality().trkKink < 20.;
945 
946  return (segmentCompatibility(muon) > (goodGlb ? 0.303 : 0.451));
947 }
948 
949 bool muon::isSoftMuon(const reco::Muon& muon, const reco::Vertex& vtx, bool run2016_hip_mitigation) {
950  bool muID = muon::isGoodMuon(muon, TMOneStationTight);
951 
952  if (!muID)
953  return false;
954 
955  bool layers = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
956  muon.innerTrack()->hitPattern().pixelLayersWithMeasurement() > 0;
957 
958  bool ishighq = muon.innerTrack()->quality(reco::Track::highPurity);
959 
960  bool ip = fabs(muon.innerTrack()->dxy(vtx.position())) < 0.3 && fabs(muon.innerTrack()->dz(vtx.position())) < 20.;
961 
962  return layers && ip && (ishighq | run2016_hip_mitigation);
963 }
964 
966  if (!muon.isGlobalMuon())
967  return false;
968 
969  bool muValHits = (muon.globalTrack()->hitPattern().numberOfValidMuonHits() > 0 ||
970  muon.tunePMuonBestTrack()->hitPattern().numberOfValidMuonHits() > 0);
971 
972  bool muMatchedSt = muon.numberOfMatchedStations() > 1;
973  if (!muMatchedSt) {
974  if (muon.isTrackerMuon() && muon.numberOfMatchedStations() == 1) {
975  if (muon.expectedNnumberOfMatchedStations() < 2 || !(muon.stationMask() == 1 || muon.stationMask() == 16) ||
976  muon.numberOfMatchedRPCLayers() > 2)
977  muMatchedSt = true;
978  }
979  }
980 
981  bool muID = muValHits && muMatchedSt;
982 
983  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
984  muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
985 
986  bool momQuality = muon.tunePMuonBestTrack()->ptError() / muon.tunePMuonBestTrack()->pt() < 0.3;
987 
988  bool ip = fabs(muon.innerTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.innerTrack()->dz(vtx.position())) < 0.5;
989 
990  return muID && hits && momQuality && ip;
991 }
992 
994  bool muID = muon.isTrackerMuon() && muon.track().isNonnull() && (muon.numberOfMatchedStations() > 1);
995  if (!muID)
996  return false;
997 
998  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
999  muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
1000 
1001  bool momQuality = muon.tunePMuonBestTrack()->ptError() / muon.tunePMuonBestTrack()->pt() < 0.3;
1002 
1003  bool ip = fabs(muon.innerTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.innerTrack()->dz(vtx.position())) < 0.5;
1004 
1005  return muID && hits && momQuality && ip;
1006 }
1007 
1008 int muon::sharedSegments(const reco::Muon& mu, const reco::Muon& mu2, unsigned int segmentArbitrationMask) {
1009  int ret = 0;
1010 
1011  // Will do with a stupid double loop, since creating and filling a map is probably _more_ inefficient for a single lookup.
1012  for (std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.matches().begin();
1013  chamberMatch != mu.matches().end();
1014  ++chamberMatch) {
1015  if (chamberMatch->segmentMatches.empty())
1016  continue;
1017  for (std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch2 = mu2.matches().begin();
1018  chamberMatch2 != mu2.matches().end();
1019  ++chamberMatch2) {
1020  if (chamberMatch2->segmentMatches.empty())
1021  continue;
1022  if (chamberMatch2->id() != chamberMatch->id())
1023  continue;
1024  for (std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
1025  segmentMatch != chamberMatch->segmentMatches.end();
1026  ++segmentMatch) {
1027  if (!segmentMatch->isMask(segmentArbitrationMask))
1028  continue;
1029  for (std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch2 = chamberMatch2->segmentMatches.begin();
1030  segmentMatch2 != chamberMatch2->segmentMatches.end();
1031  ++segmentMatch2) {
1032  if (!segmentMatch2->isMask(segmentArbitrationMask))
1033  continue;
1034  if ((segmentMatch->cscSegmentRef.isNonnull() &&
1035  segmentMatch->cscSegmentRef == segmentMatch2->cscSegmentRef) ||
1036  (segmentMatch->dtSegmentRef.isNonnull() && segmentMatch->dtSegmentRef == segmentMatch2->dtSegmentRef)) {
1037  ++ret;
1038  } // is the same
1039  } // segment of mu2 in chamber
1040  } // segment of mu1 in chamber
1041  } // chamber of mu2
1042  } // chamber of mu1
1043 
1044  return ret;
1045 }
1046 
1048  const auto& combinedTime = muon.time();
1049  const auto& rpcTime = muon.rpcTime();
1050  bool combinedTimeIsOk = (combinedTime.nDof > 7);
1051  bool rpcTimeIsOk = (rpcTime.nDof > 1 && fabs(rpcTime.timeAtIpInOutErr) < 0.001);
1052  bool outOfTime = false;
1053  if (rpcTimeIsOk) {
1054  if ((fabs(rpcTime.timeAtIpInOut) > 10) && !(combinedTimeIsOk && fabs(combinedTime.timeAtIpInOut) < 10))
1055  outOfTime = true;
1056  } else {
1057  if (combinedTimeIsOk && (combinedTime.timeAtIpInOut > 20 || combinedTime.timeAtIpInOut < -45))
1058  outOfTime = true;
1059  }
1060  return outOfTime;
1061 }
1062 
1064  reco::Vertex const* vertex,
1065  bool run2016_hip_mitigation) {
1066  // https://twiki.cern.ch/twiki/bin/viewauth/CMS/SWGuideMuonIdRun2
1067  unsigned int selectors = muon.selectors();
1068  // Compute Id and Isolation variables
1069  double chIso = muon.pfIsolationR04().sumChargedHadronPt;
1070  double nIso = muon.pfIsolationR04().sumNeutralHadronEt;
1071  double phoIso = muon.pfIsolationR04().sumPhotonEt;
1072  double puIso = muon.pfIsolationR04().sumPUPt;
1073  double dbCorrectedIsolation = chIso + std::max(nIso + phoIso - .5 * puIso, 0.);
1074  double dbCorrectedRelIso = dbCorrectedIsolation / muon.pt();
1075  double tkRelIso = muon.isolationR03().sumPt / muon.pt();
1076 
1077  // Base selectors
1078  if (muon::isLooseMuon(muon))
1079  selectors |= reco::Muon::CutBasedIdLoose;
1080  if (vertex) {
1081  if (muon::isTightMuon(muon, *vertex))
1082  selectors |= reco::Muon::CutBasedIdTight;
1083  if (muon::isSoftMuon(muon, *vertex, run2016_hip_mitigation))
1084  selectors |= reco::Muon::SoftCutBasedId;
1085  if (muon::isHighPtMuon(muon, *vertex))
1087  if (muon::isTrackerHighPtMuon(muon, *vertex))
1088  selectors |= reco::Muon::CutBasedIdTrkHighPt;
1089  }
1090  if (muon::isMediumMuon(muon, run2016_hip_mitigation)) {
1091  selectors |= reco::Muon::CutBasedIdMedium;
1092  if (vertex and fabs(muon.muonBestTrack()->dz(vertex->position())) < 0.1 and
1093  fabs(muon.muonBestTrack()->dxy(vertex->position())) < 0.02)
1095  }
1096 
1097  // PF isolation
1098  if (dbCorrectedRelIso < 0.40)
1099  selectors |= reco::Muon::PFIsoVeryLoose;
1100  if (dbCorrectedRelIso < 0.25)
1101  selectors |= reco::Muon::PFIsoLoose;
1102  if (dbCorrectedRelIso < 0.20)
1103  selectors |= reco::Muon::PFIsoMedium;
1104  if (dbCorrectedRelIso < 0.15)
1105  selectors |= reco::Muon::PFIsoTight;
1106  if (dbCorrectedRelIso < 0.10)
1107  selectors |= reco::Muon::PFIsoVeryTight;
1108  if (dbCorrectedRelIso < 0.05)
1109  selectors |= reco::Muon::PFIsoVeryVeryTight;
1110 
1111  // Tracker isolation
1112  if (tkRelIso < 0.10)
1113  selectors |= reco::Muon::TkIsoLoose;
1114  if (tkRelIso < 0.05)
1115  selectors |= reco::Muon::TkIsoTight;
1116 
1117  // Trigger selectors
1118  if (isLooseTriggerMuon(muon))
1119  selectors |= reco::Muon::TriggerIdLoose;
1120 
1121  // Timing
1122  if (!outOfTimeMuon(muon))
1123  selectors |= reco::Muon::InTimeMuon;
1124 
1125  return static_cast<reco::Muon::Selector>(selectors);
1126 }
float chi2LocalPosition
chi2 value for the STA-TK matching of local position
Definition: MuonQuality.h:19
int chamber() const
Definition: CSCDetId.h:62
type
Definition: HCALResponse.h:21
float segmentX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
int numberOfMatchedRPCLayers(ArbitrationType type=RPCHitAndTrackArbitration) const
TrackRef track() const override
reference to a Track
Definition: Muon.h:46
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
double eta() const final
momentum pseudorapidity
static constexpr int GEM
Definition: MuonSubdetId.h:14
bool isStandAloneMuon() const override
Definition: Muon.h:300
float sumPt
sum-pt of tracks
Definition: MuonIsolation.h:6
AlgorithmType
Definition: MuonSelectors.h:76
bool isMediumMuon(const reco::Muon &, bool run2016_hip_mitigation=false)
virtual TrackRef innerTrack() const
Definition: Muon.h:45
ret
prodAgent to be discontinued
double pt() const final
transverse momentum
float trkKink
value of the kink algorithm applied to the inner track stub
Definition: MuonQuality.h:11
bool isMatchesValid() const
Definition: Muon.h:143
float sumPhotonEt
sum pt of PF photons
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:113
MuonTime time() const
get DT/CSC combined timing information
Definition: Muon.h:132
bool isME0Muon() const
Definition: Muon.h:305
float caloCompatibility(const reco::Muon &muon)
bool isLooseMuon(const reco::Muon &)
float sumNeutralHadronEt
sum pt of neutral hadrons
SelectionType
Selector type.
Definition: MuonSelectors.h:19
bool isTrackerMuon() const override
Definition: Muon.h:299
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
unsigned int expectedNnumberOfMatchedStations(float minDistanceFromEdge=10.0) const
int endcap() const
Definition: CSCDetId.h:85
bool isRPCMuon() const
Definition: Muon.h:303
char const * label
float trackDistErr(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
float dY(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
virtual TrackRef muonBestTrack() const
Definition: Muon.h:60
bool isGlobalMuon() const override
Definition: Muon.h:298
ArbitrationType
define arbitration schemes
Definition: Muon.h:187
float pullDxDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
T sqrt(T t)
Definition: SSEVec.h:19
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
MuonQuality combinedQuality() const
get energy deposition information
Definition: Muon.h:119
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float segmentCompatibility(const reco::Muon &muon, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
static constexpr int ME0
Definition: MuonSubdetId.h:15
Definition: value.py:1
bool isQualityValid() const
Definition: Muon.h:117
bool isSoftMuon(const reco::Muon &, const reco::Vertex &, bool run2016_hip_mitigation=false)
float dX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
float pullY(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
virtual TrackRef outerTrack() const
reference to Track reconstructed in the muon detector only
Definition: Muon.h:48
bool isGoodMuon(const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
main GoodMuon wrapper call
unsigned int stationMask(ArbitrationType type=SegmentAndTrackArbitration) const
int numberOfMatchedStations(ArbitrationType type=SegmentAndTrackArbitration) const
reco::Muon::Selector makeSelectorBitset(reco::Muon const &muon, reco::Vertex const *vertex=0, bool run2016_hip_mitigation=false)
T Max(T a, T b)
Definition: MathUtil.h:44
int ring() const
Definition: CSCDetId.h:68
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)
bool isHighPtMuon(const reco::Muon &, const reco::Vertex &)
SelectionType selectionTypeFromString(const std::string &label)
Definition: MuonSelectors.cc:9
const MuonPFIsolation & pfIsolationR04() const
Definition: Muon.h:172
bool isGEMMuon() const
Definition: Muon.h:304
a lightweight "map" for selection type string label and enum value
Definition: MuonSelectors.h:62
uint64_t selectors() const
Definition: Muon.h:239
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:145
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
MuonTime rpcTime() const
get RPC timing information
Definition: Muon.h:136
float pullDyDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
float trackDist(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
bool isPFMuon() const
Definition: Muon.h:302
bool outOfTimeMuon(const reco::Muon &muon)
int station() const
Definition: CSCDetId.h:79
virtual TrackRef tunePMuonBestTrack() const
Definition: Muon.h:62
bool isTrackerHighPtMuon(const reco::Muon &, const reco::Vertex &)
bool isTightMuon(const reco::Muon &, const reco::Vertex &)
static constexpr int CSC
Definition: MuonSubdetId.h:12
bool isLooseTriggerMuon(const reco::Muon &)
Selector
Definition: Muon.h:200
float caloCompatibility() const
Definition: Muon.h:158
const MuonIsolation & isolationR03() const
Definition: Muon.h:166
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:51
float pullX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
float sumChargedHadronPt
sum-pt of charged Hadron