20 throw cms::Exception(
"MuonSelectorError") << label <<
" is not a recognized SelectionType";
35 throw cms::Exception(
"MuonSelectorError") << label <<
" is not a recognized reco::Muon::Selector";
41 double maxChamberDist,
42 double maxChamberDistPull,
44 unsigned int theMask = 0;
46 for (
int stationIdx = 1; stationIdx < 5; ++stationIdx) {
47 for (
int detectorIdx = 1; detectorIdx < 3; ++detectorIdx) {
48 float dist = muon.
trackDist(stationIdx, detectorIdx, arbitrationType);
49 if (dist < maxChamberDist &&
50 dist < maxChamberDistPull * muon.
trackDistErr(stationIdx, detectorIdx, arbitrationType))
51 theMask += 1 << ((stationIdx - 1) + 4 * (detectorIdx - 1));
62 bool use_weight_regain_at_chamber_boundary =
true;
63 bool use_match_dist_penalty =
true;
65 int nr_of_stations_crossed = 0;
66 int nr_of_stations_with_segment = 0;
67 std::vector<int> stations_w_track(8);
68 std::vector<int> station_has_segmentmatch(8);
69 std::vector<int> station_was_crossed(8);
70 std::vector<float> stations_w_track_at_boundary(8);
71 std::vector<float> station_weight(8);
72 int position_in_stations = 0;
73 float full_weight = 0.;
75 for (
int i = 1;
i <= 8; ++
i) {
80 float thisTrackDist = muon.
trackDist(
i, 1, arbitrationType);
81 if (thisTrackDist < 999999) {
82 ++nr_of_stations_crossed;
83 station_was_crossed[
i - 1] = 1;
84 if (thisTrackDist > -10.)
85 stations_w_track_at_boundary[
i - 1] = thisTrackDist;
87 stations_w_track_at_boundary[
i - 1] = 0.;
90 if (muon.
segmentX(
i, 1, arbitrationType) < 999999) {
91 ++nr_of_stations_with_segment;
92 station_has_segmentmatch[
i - 1] = 1;
95 float thisTrackDist = muon.
trackDist(
i - 4, 2, arbitrationType);
96 if (thisTrackDist < 999999) {
97 ++nr_of_stations_crossed;
98 station_was_crossed[
i - 1] = 1;
99 if (thisTrackDist > -10.)
100 stations_w_track_at_boundary[
i - 1] = thisTrackDist;
102 stations_w_track_at_boundary[
i - 1] = 0.;
105 if (muon.
segmentX(
i - 4, 2, arbitrationType) < 999999) {
106 ++nr_of_stations_with_segment;
107 station_has_segmentmatch[
i - 1] = 1;
126 const float attenuate_weight_regain = 0.5;
128 for (
int i = 1;
i <= 8; ++
i) {
134 if (station_was_crossed[
i - 1] > 0) {
140 ++position_in_stations;
142 switch (nr_of_stations_crossed) {
144 station_weight[
i - 1] = 1.f;
147 if (position_in_stations == 1)
148 station_weight[
i - 1] = 0.33f;
150 station_weight[
i - 1] = 0.67f;
153 if (position_in_stations == 1)
154 station_weight[
i - 1] = 0.23f;
155 else if (position_in_stations == 2)
156 station_weight[
i - 1] = 0.33f;
158 station_weight[
i - 1] = 0.44f;
161 if (position_in_stations == 1)
162 station_weight[
i - 1] = 0.10f;
163 else if (position_in_stations == 2)
164 station_weight[
i - 1] = 0.20f;
165 else if (position_in_stations == 3)
166 station_weight[
i - 1] = 0.30f;
168 station_weight[
i - 1] = 0.40f;
175 station_weight[
i - 1] = 1.f / nr_of_stations_crossed;
178 if (use_weight_regain_at_chamber_boundary) {
179 if (station_has_segmentmatch[
i - 1] <= 0 && stations_w_track_at_boundary[
i - 1] != 0.) {
184 station_weight[
i - 1] = station_weight[
i - 1] * attenuate_weight_regain * 0.5f *
185 (std::erf(stations_w_track_at_boundary[
i - 1] / 6.
f) + 1.f);
186 }
else if (station_has_segmentmatch[
i - 1] <= 0 &&
187 stations_w_track_at_boundary[
i - 1] == 0.
f) {
189 station_weight[
i - 1] = 0.f;
192 if (station_has_segmentmatch[
i - 1] <= 0)
193 station_weight[
i - 1] = 0.f;
197 if (station_has_segmentmatch[
i - 1] > 0) {
199 if (muon.
dY(
i, 1, arbitrationType) < 999999.f &&
200 muon.
dX(
i, 1, arbitrationType) < 999999.f) {
201 const float pullTot2 =
203 if (pullTot2 > 1.
f) {
207 if (use_match_dist_penalty) {
209 if (dxy2 < 9.f && pullTot2 > 9.
f) {
211 station_weight[
i - 1] *= 1.f /
std::pow(dxy2, .125);
213 station_weight[
i - 1] *= 1.f /
std::pow(pullTot2, .125);
217 }
else if (muon.
dY(
i, 1, arbitrationType) >= 999999.f) {
219 if (muon.
pullX(
i, 1, arbitrationType) > 1.f) {
221 if (use_match_dist_penalty) {
223 if (muon.
dX(
i, 1, arbitrationType) < 3.f && muon.
pullX(
i, 1, arbitrationType) > 3.f) {
224 if (muon.
dX(
i, 1, arbitrationType) > 1.f)
225 station_weight[
i - 1] *= 1.
f /
std::pow(muon.
dX(
i, 1, arbitrationType), .25);
227 station_weight[
i - 1] *= 1.f /
std::pow(muon.
pullX(
i, 1, arbitrationType), .25);
233 if (muon.
pullY(
i, 1, arbitrationType) > 1.f) {
235 if (use_match_dist_penalty) {
237 if (muon.
dY(
i, 1, arbitrationType) < 3. && muon.
pullY(
i, 1, arbitrationType) > 3.) {
238 if (muon.
dY(
i, 1, arbitrationType) > 1.f)
239 station_weight[
i - 1] *= 1.
f /
std::pow(muon.
dY(
i, 1, arbitrationType), .25);
241 station_weight[
i - 1] *= 1.f /
std::pow(muon.
pullY(
i, 1, arbitrationType), .25);
247 const float pullTot2 =
249 if (pullTot2 > 1.
f) {
251 if (use_match_dist_penalty) {
255 if (dxy2 < 9.f && pullTot2 > 9.
f) {
257 station_weight[
i - 1] *= 1.f /
std::pow(dxy2, .125);
259 station_weight[
i - 1] *= 1.f /
std::pow(pullTot2, .125);
271 station_weight[
i - 1] = 0.f;
275 full_weight += station_weight[
i - 1];
281 if (nr_of_stations_crossed == 0) {
295 double minCompatibility,
299 bool goodMuon =
false;
325 double maxChamberDist,
326 double maxChamberDistPull,
328 bool syncMinNMatchesNRequiredStationsInBarrelOnly,
329 bool applyAlsoAngularCuts) {
332 bool goodMuon =
false;
337 if (minNumberOfMatches == 0)
340 unsigned int theStationMask = muon.
stationMask(arbitrationType);
341 unsigned int theRequiredStationMask =
346 int numRequiredStations = 0;
347 for (
int it = 0; it < 8; ++it) {
348 if (theStationMask & 1 << it)
350 if (theRequiredStationMask & 1 << it)
351 ++numRequiredStations;
356 if (syncMinNMatchesNRequiredStationsInBarrelOnly) {
359 if (minNumberOfMatches > numRequiredStations)
360 minNumberOfMatches = numRequiredStations;
361 if (minNumberOfMatches < 1)
362 minNumberOfMatches = 1;
365 if (minNumberOfMatches > numRequiredStations)
366 minNumberOfMatches = numRequiredStations;
367 if (minNumberOfMatches < 1)
368 minNumberOfMatches = 1;
371 if (numSegs >= minNumberOfMatches)
379 if (theRequiredStationMask) {
380 for (
int stationIdx = 7; stationIdx >= 0; --stationIdx)
381 if (theRequiredStationMask & 1 << stationIdx) {
382 if (theStationMask & 1 << stationIdx) {
383 lastSegBit = stationIdx;
392 for (
int stationIdx = 7; stationIdx >= 0; --stationIdx)
393 if (theStationMask & 1 << stationIdx) {
394 lastSegBit = stationIdx;
404 station = lastSegBit < 4 ? lastSegBit + 1 : lastSegBit - 3;
405 detector = lastSegBit < 4 ? 1 : 2;
408 if (
std::abs(muon.
pullX(station, detector, arbitrationType,
true)) > maxAbsPullX &&
409 std::abs(muon.
dX(station, detector, arbitrationType)) > maxAbsDx)
412 if (applyAlsoAngularCuts &&
std::abs(muon.
pullDxDz(station, detector, arbitrationType,
true)) > maxAbsPullX)
416 if (maxAbsDy < 999999) {
420 if (
std::abs(muon.
pullY(station, 2, arbitrationType,
true)) > maxAbsPullY &&
421 std::abs(muon.
dY(station, 2, arbitrationType)) > maxAbsDy)
424 if (applyAlsoAngularCuts &&
std::abs(muon.
pullDyDz(station, 2, arbitrationType,
true)) > maxAbsPullY)
442 for (
int stationIdx = station; stationIdx > 0; --stationIdx) {
443 if (!(theStationMask & 1 << (stationIdx - 1)))
446 if (muon.
dY(stationIdx, 1, arbitrationType) > 999998)
449 if (
std::abs(muon.
pullY(stationIdx, 1, arbitrationType,
true)) > maxAbsPullY &&
450 std::abs(muon.
dY(stationIdx, 1, arbitrationType)) > maxAbsDy) {
454 if (applyAlsoAngularCuts &&
std::abs(muon.
pullDyDz(stationIdx, 1, arbitrationType,
true)) > maxAbsPullY)
477 unsigned int theStationMask = muon.
stationMask(arbitrationType);
488 bool existsGoodDTSegX =
false;
489 bool existsDTSegY =
false;
493 for (
int stationIdx = 0; stationIdx <= 7; ++stationIdx)
494 if (theStationMask & 1 << stationIdx) {
495 station = stationIdx < 4 ? stationIdx + 1 : stationIdx - 3;
496 detector = stationIdx < 4 ? 1 : 2;
498 if ((
std::abs(muon.
pullX(station, detector, arbitrationType,
true)) > maxAbsPullX &&
499 std::abs(muon.
dX(station, detector, arbitrationType)) > maxAbsDx) ||
500 (applyAlsoAngularCuts &&
std::abs(muon.
pullDxDz(station, detector, arbitrationType,
true)) > maxAbsPullX))
502 else if (detector == 1)
503 existsGoodDTSegX =
true;
506 if (maxAbsDy < 999999) {
508 if ((
std::abs(muon.
pullY(station, 2, arbitrationType,
true)) > maxAbsPullY &&
509 std::abs(muon.
dY(station, 2, arbitrationType)) > maxAbsDy) ||
510 (applyAlsoAngularCuts &&
std::abs(muon.
pullDyDz(station, 2, arbitrationType,
true)) > maxAbsPullY))
513 if (muon.
dY(station, 1, arbitrationType) > 999998)
518 if ((
std::abs(muon.
pullY(station, 1, arbitrationType,
true)) > maxAbsPullY &&
519 std::abs(muon.
dY(station, 1, arbitrationType)) > maxAbsDy) ||
520 (applyAlsoAngularCuts &&
std::abs(muon.
pullDyDz(station, 1, arbitrationType,
true)) > maxAbsPullY)) {
537 if (maxAbsDy < 999999) {
540 else if (existsGoodDTSegX)
547 if (minNumberOfMatches == 0)
551 for (
const auto& chamberMatch : muon.
matches()) {
555 const float trkX = chamberMatch.x;
556 const float errX = chamberMatch.xErr;
558 for (
const auto& rpcMatch : chamberMatch.rpcMatches) {
559 const float rpcX = rpcMatch.x;
560 const float dX =
std::abs(rpcX - trkX);
561 if (dX < maxAbsDx
or dX < maxAbsPullX * errX) {
568 if (nMatch >= minNumberOfMatches)
575 if (minNumberOfMatches == 0)
579 for (
const auto& chamberMatch : muon.
matches()) {
583 const float trkX = chamberMatch.x;
584 const float errX2 = chamberMatch.xErr * chamberMatch.xErr;
585 const float trkY = chamberMatch.y;
586 const float errY2 = chamberMatch.yErr * chamberMatch.yErr;
588 for (
const auto& segment : chamberMatch.me0Matches) {
589 const float me0X = segment.x;
590 const float me0ErrX2 = segment.xErr * segment.xErr;
591 const float me0Y = segment.y;
592 const float me0ErrY2 = segment.yErr * segment.yErr;
594 const float dX =
std::abs(me0X - trkX);
595 const float dY =
std::abs(me0Y - trkY);
596 const float invPullX2 = errX2 + me0ErrX2;
597 const float invPullY2 = errY2 + me0ErrY2;
599 if ((dX < maxAbsDx
or dX < maxAbsPullX *
std::sqrt(invPullX2)) and
600 (dY < maxAbsDy
or dY < maxAbsPullY *
std::sqrt(invPullY2))) {
607 return (nMatch >= minNumberOfMatches);
611 if (minNumberOfMatches == 0)
615 for (
const auto& chamberMatch : muon.
matches()) {
619 const float trkX = chamberMatch.x;
620 const float errX2 = chamberMatch.xErr * chamberMatch.xErr;
621 const float trkY = chamberMatch.y;
622 const float errY2 = chamberMatch.yErr * chamberMatch.yErr;
624 for (
const auto& segment : chamberMatch.gemMatches) {
625 const float gemX = segment.x;
626 const float gemErrX2 = segment.xErr * segment.xErr;
627 const float gemY = segment.y;
628 const float gemErrY2 = segment.yErr * segment.yErr;
630 const float dX =
std::abs(gemX - trkX);
631 const float dY =
std::abs(gemY - trkY);
632 const float invPullX2 = errX2 + gemErrX2;
633 const float invPullY2 = errY2 + gemErrY2;
635 if ((dX < maxAbsDx
or dX < maxAbsPullX *
std::sqrt(invPullX2)) and
636 (dY < maxAbsDy
or dY < maxAbsPullY *
std::sqrt(invPullY2))) {
643 return (nMatch >= minNumberOfMatches);
671 muon.
globalTrack()->hitPattern().numberOfValidMuonHits() > 0;
685 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType,
true,
false);
689 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType,
true,
false);
693 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType,
false,
false);
697 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType,
false,
false);
702 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType,
false,
false);
705 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType,
false,
false);
710 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType,
false,
false);
713 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType,
false,
false);
736 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType,
false,
true);
740 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType,
false,
true);
744 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType,
false,
true);
748 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType,
false,
true);
753 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 1E9, 1E9, 1E9, 1E9, arbitrationType,
false,
false);
756 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 1E9, 1E9, -3, -3, arbitrationType,
true,
false);
761 isGoodMuon(muon,
TMOneStation, 1, 3, 3, 3, 3, 1E9, 1E9, arbitrationType,
false,
false);
764 isGoodMuon(muon,
TMLastStation, 2, 3, 3, 3, 3, -3, -3, arbitrationType,
true,
false);
767 return muon.
isRPCMuon() &&
isGoodMuon(muon,
RPCMu, 2, 20, 4, 1e9, 1e9, 1e9, 1e9, arbitrationType,
false,
false);
774 isGoodMuon(muon,
ME0Mu, 1, 1e9, 1e9, 1e9, 1e9, 1e9, 1e9, arbitrationType,
false,
false);
781 isGoodMuon(muon,
GEMMu, 1, 1e9, 1e9, 1e9, 1e9, 1e9, 1e9, arbitrationType,
false,
false);
792 const reco::Muon& muon1,
const reco::Muon& muon2,
double pullX,
double pullY,
bool checkAdjacentChambers) {
795 unsigned int betterMuon = (muon1.
pt() > muon2.
pt() ? 1 : 2);
796 for (std::vector<reco::MuonChamberMatch>::const_iterator chamber1 = muon1.
matches().begin();
797 chamber1 != muon1.
matches().end();
799 for (std::vector<reco::MuonChamberMatch>::const_iterator chamber2 = muon2.
matches().begin();
800 chamber2 != muon2.
matches().end();
806 if (chamber1->id == chamber2->id) {
808 if (
std::abs(chamber1->x - chamber2->x) <
809 pullX *
sqrt(chamber1->xErr * chamber1->xErr + chamber2->xErr * chamber2->xErr)) {
814 if (nMatches1 == 0 || nMatches2 == 0)
818 if (
std::abs(chamber1->y - chamber2->y) <
819 pullY *
sqrt(chamber1->yErr * chamber1->yErr + chamber2->yErr * chamber2->yErr)) {
824 if (nMatches1 == 0 || nMatches2 == 0)
828 if (!checkAdjacentChambers)
849 if (
std::abs(chamber1->edgeX) > chamber1->xErr * pullX)
851 if (
std::abs(chamber2->edgeX) > chamber2->xErr * pullX)
853 if (chamber1->x * chamber2->x < 0) {
858 if (nMatches1 == 0 || nMatches2 == 0)
873 bool layer_requirements = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
874 muon.
innerTrack()->hitPattern().pixelLayersWithMeasurement() > 0;
875 bool match_requirements =
877 return layer_requirements and match_requirements;
886 bool hits = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
887 muon.
innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
892 return muID && hits && ip;
902 if (run2016_hip_mitigation) {
903 if (muon.
innerTrack()->validFraction() < 0.49)
922 bool layers = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
923 muon.
innerTrack()->hitPattern().pixelLayersWithMeasurement() > 0;
930 return layers && ip && (ishighq | run2016_hip_mitigation);
937 bool muValHits = (muon.
globalTrack()->hitPattern().numberOfValidMuonHits() > 0 ||
949 bool muID = muValHits && muMatchedSt;
951 bool hits = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
952 muon.
innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
959 return muID && hits && momQuality && ip;
967 bool hits = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
968 muon.
innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
975 return muID && hits && momQuality && ip;
982 for (std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.
matches().begin();
983 chamberMatch != mu.
matches().end();
985 if (chamberMatch->segmentMatches.empty())
987 for (std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch2 = mu2.
matches().begin();
988 chamberMatch2 != mu2.
matches().end();
990 if (chamberMatch2->segmentMatches.empty())
992 if (chamberMatch2->id() != chamberMatch->id())
994 for (std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
995 segmentMatch != chamberMatch->segmentMatches.end();
997 if (!segmentMatch->isMask(segmentArbitrationMask))
999 for (std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch2 = chamberMatch2->segmentMatches.begin();
1000 segmentMatch2 != chamberMatch2->segmentMatches.end();
1002 if (!segmentMatch2->isMask(segmentArbitrationMask))
1004 if ((segmentMatch->cscSegmentRef.isNonnull() &&
1005 segmentMatch->cscSegmentRef == segmentMatch2->cscSegmentRef) ||
1006 (segmentMatch->dtSegmentRef.isNonnull() && segmentMatch->dtSegmentRef == segmentMatch2->dtSegmentRef)) {
1018 const auto& combinedTime = muon.
time();
1019 const auto& rpcTime = muon.
rpcTime();
1020 bool combinedTimeIsOk = (combinedTime.nDof > 7);
1021 bool rpcTimeIsOk = (rpcTime.nDof > 1 &&
std::abs(rpcTime.timeAtIpInOutErr) < 0.001);
1022 bool outOfTime =
false;
1024 if ((
std::abs(rpcTime.timeAtIpInOut) > 10) && !(combinedTimeIsOk &&
std::abs(combinedTime.timeAtIpInOut) < 10))
1027 if (combinedTimeIsOk && (combinedTime.timeAtIpInOut > 20 || combinedTime.timeAtIpInOut < -45))
1035 bool run2016_hip_mitigation) {
1043 double dbCorrectedIsolation = chIso +
std::max(nIso + phoIso - .5 * puIso, 0.);
1044 double dbCorrectedRelIso = dbCorrectedIsolation / muon.
pt();
1068 if (dbCorrectedRelIso < 0.40)
1070 if (dbCorrectedRelIso < 0.25)
1072 if (dbCorrectedRelIso < 0.20)
1074 if (dbCorrectedRelIso < 0.15)
1076 if (dbCorrectedRelIso < 0.10)
1078 if (dbCorrectedRelIso < 0.05)
1082 if (tkRelIso < 0.10)
1084 if (tkRelIso < 0.05)
float chi2LocalPosition
chi2 value for the STA-TK matching of local position
float segmentX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
int numberOfMatchedRPCLayers(ArbitrationType type=RPCHitAndTrackArbitration) const
tuple ret
prodAgent to be discontinued
bool isNonnull() const
Checks for non-null.
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
float sumPt
sum-pt of tracks
double pt() const final
transverse momentum
reco::Muon::Selector value
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::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
bool isMediumMuon(const reco::Muon &, bool run2016_hip_mitigation=false)
virtual TrackRef innerTrack() const
float trkKink
value of the kink algorithm applied to the inner track stub
bool isMatchesValid() const
float sumPhotonEt
sum pt of PF photons
TrackRef track() const override
reference to a Track
reco::Muon::Selector makeSelectorBitset(reco::Muon const &muon, reco::Vertex const *vertex=nullptr, bool run2016_hip_mitigation=false)
float trkRelChi2
chi2 value for the inner track stub with respect to the global track
const Point & position() const
position
MuonTime time() const
get DT/CSC combined timing information
bool isTrackerMuon() const override
float caloCompatibility(const reco::Muon &muon)
bool isLooseMuon(const reco::Muon &)
float sumNeutralHadronEt
sum pt of neutral hadrons
SelectionType
Selector type.
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
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
ArbitrationType
define arbitration schemes
float pullDxDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
MuonQuality combinedQuality() const
get energy deposition information
Abs< T >::type abs(const T &t)
float segmentCompatibility(const reco::Muon &muon, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
bool isQualityValid() const
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
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
float staRelChi2
chi2 value for the outer track stub with respect to the global track
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)
const reco::Muon::ArbitrationType arbitrationType
const MuonPFIsolation & pfIsolationR04() const
uint64_t selectors() const
unsigned int RequiredStationMask(const reco::Muon &muon, double maxChamberDist, double maxChamberDistPull, reco::Muon::ArbitrationType arbitrationType)
std::vector< MuonChamberMatch > & matches()
get muon matching information
static const SelectorStringToEnum selectorStringToEnumMap[]
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
MuonTime rpcTime() const
get RPC timing information
float pullDyDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
reco::Muon::Selector selectorFromString(const std::string &label)
float trackDist(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
bool outOfTimeMuon(const reco::Muon &muon)
virtual TrackRef tunePMuonBestTrack() const
bool isTrackerHighPtMuon(const reco::Muon &, const reco::Vertex &)
bool isTightMuon(const reco::Muon &, const reco::Vertex &)
bool isGlobalMuon() const override
bool isLooseTriggerMuon(const reco::Muon &)
float caloCompatibility() const
Power< A, B >::type pow(const A &a, const B &b)
const MuonIsolation & isolationR03() const
static const SelectionTypeStringToEnum selectionTypeStringToEnumMap[]
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
bool isStandAloneMuon() const override
float pullX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
float sumChargedHadronPt
sum-pt of charged Hadron
double eta() const final
momentum pseudorapidity