CMS 3D CMS Logo

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