42 for(
int i = 0; selectionTypeStringToEnumMap[
i].
label && (!
found); ++
i)
43 if (! strcmp(label.c_str(), selectionTypeStringToEnumMap[
i].
label)) {
45 value = selectionTypeStringToEnumMap[
i].
value;
49 if (! found)
throw cms::Exception(
"MuonSelectorError") << label <<
" is not a recognized SelectionType";
55 double maxChamberDist,
56 double maxChamberDistPull,
59 unsigned int theMask = 0;
61 for(
int stationIdx = 1; stationIdx < 5; ++stationIdx)
62 for(
int detectorIdx = 1; detectorIdx < 3; ++detectorIdx)
63 if(muon.
trackDist(stationIdx,detectorIdx,arbitrationType) < maxChamberDist &&
64 muon.
trackDist(stationIdx,detectorIdx,arbitrationType)/muon.
trackDistErr(stationIdx,detectorIdx,arbitrationType) < maxChamberDistPull)
65 theMask += 1<<((stationIdx-1)+4*(detectorIdx-1));
77 bool use_weight_regain_at_chamber_boundary =
true;
78 bool use_match_dist_penalty =
true;
80 int nr_of_stations_crossed = 0;
81 int nr_of_stations_with_segment = 0;
82 std::vector<int> stations_w_track(8);
83 std::vector<int> station_has_segmentmatch(8);
84 std::vector<int> station_was_crossed(8);
85 std::vector<float> stations_w_track_at_boundary(8);
86 std::vector<float> station_weight(8);
87 int position_in_stations = 0;
88 float full_weight = 0.;
90 for(
int i = 1;
i<=8; ++
i) {
95 if( muon.
trackDist(
i,1,arbitrationType) < 999999 ) {
96 ++nr_of_stations_crossed;
97 station_was_crossed[
i-1] = 1;
98 if(muon.
trackDist(
i,1,arbitrationType) > -10. ) stations_w_track_at_boundary[
i-1] = muon.
trackDist(
i,1,arbitrationType);
99 else stations_w_track_at_boundary[
i-1] = 0.;
101 if( muon.
segmentX(
i,1,arbitrationType) < 999999 ) {
102 ++nr_of_stations_with_segment;
103 station_has_segmentmatch[
i-1] = 1;
107 if( muon.
trackDist(
i-4,2,arbitrationType) < 999999 ) {
108 ++nr_of_stations_crossed;
109 station_was_crossed[
i-1] = 1;
110 if(muon.
trackDist(
i-4,2,arbitrationType) > -10. ) stations_w_track_at_boundary[
i-1] = muon.
trackDist(
i-4,2,arbitrationType);
111 else stations_w_track_at_boundary[
i-1] = 0.;
113 if( muon.
segmentX(
i-4,2,arbitrationType) < 999999 ) {
114 ++nr_of_stations_with_segment;
115 station_has_segmentmatch[
i-1] = 1;
134 const float attenuate_weight_regain = 0.5;
136 for(
int i = 1;
i<=8; ++
i) {
142 if( station_was_crossed[
i-1] > 0 ) {
148 ++position_in_stations;
150 switch ( nr_of_stations_crossed ) {
152 station_weight[
i-1] = 1.;
155 if ( position_in_stations == 1 ) station_weight[
i-1] = 0.33;
156 else station_weight[
i-1] = 0.67;
159 if ( position_in_stations == 1 ) station_weight[
i-1] = 0.23;
160 else if( position_in_stations == 2 ) station_weight[
i-1] = 0.33;
161 else station_weight[
i-1] = 0.44;
164 if ( position_in_stations == 1 ) station_weight[
i-1] = 0.10;
165 else if( position_in_stations == 2 ) station_weight[
i-1] = 0.20;
166 else if( position_in_stations == 3 ) station_weight[
i-1] = 0.30;
167 else station_weight[
i-1] = 0.40;
174 station_weight[
i-1] = 1./nr_of_stations_crossed;
177 if( use_weight_regain_at_chamber_boundary ) {
178 if(station_has_segmentmatch[
i-1] <= 0 && stations_w_track_at_boundary[
i-1] != 0. ) {
182 station_weight[
i-1] = station_weight[
i-1]*attenuate_weight_regain*0.5*(TMath::Erf(stations_w_track_at_boundary[
i-1]/6.)+1.);
184 else if(station_has_segmentmatch[
i-1] <= 0 && stations_w_track_at_boundary[
i-1] == 0.) {
186 station_weight[
i-1] = 0.;
190 if(station_has_segmentmatch[
i-1] <= 0) station_weight[
i-1] = 0.;
193 if( station_has_segmentmatch[
i-1] > 0 && 42 == 42 ) {
195 if( muon.
dY(
i,1,arbitrationType) < 999999 && muon.
dX(
i,1,arbitrationType) < 999999) {
197 TMath::Sqrt(TMath::Power(muon.
pullX(
i,1,arbitrationType),2.)+TMath::Power(muon.
pullY(
i,1,arbitrationType),2.))> 1. ) {
199 if(use_match_dist_penalty) {
201 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. ) {
202 station_weight[
i-1] *= 1./TMath::Power(
203 TMath::Max((
double)TMath::Sqrt(TMath::Power(muon.
dX(
i,1,arbitrationType),2.)+TMath::Power(muon.
dY(
i,1,arbitrationType),2.)),(
double)1.),.25);
206 station_weight[
i-1] *= 1./TMath::Power(
207 TMath::Sqrt(TMath::Power(muon.
pullX(
i,1,arbitrationType),2.)+TMath::Power(muon.
pullY(
i,1,arbitrationType),2.)),.25);
212 else if (muon.
dY(
i,1,arbitrationType) >= 999999) {
213 if( muon.
pullX(
i,1,arbitrationType) > 1. ) {
215 if(use_match_dist_penalty) {
217 if( muon.
dX(
i,1,arbitrationType) < 3. && muon.
pullX(
i,1,arbitrationType) > 3. ) {
218 station_weight[
i-1] *= 1./TMath::Power(
TMath::Max((
double)muon.
dX(
i,1,arbitrationType),(double)1.),.25);
221 station_weight[
i-1] *= 1./TMath::Power(muon.
pullX(
i,1,arbitrationType),.25);
227 if( muon.
pullY(
i,1,arbitrationType) > 1. ) {
229 if(use_match_dist_penalty) {
231 if( muon.
dY(
i,1,arbitrationType) < 3. && muon.
pullY(
i,1,arbitrationType) > 3. ) {
232 station_weight[
i-1] *= 1./TMath::Power(
TMath::Max((
double)muon.
dY(
i,1,arbitrationType),(double)1.),.25);
235 station_weight[
i-1] *= 1./TMath::Power(muon.
pullY(
i,1,arbitrationType),.25);
243 TMath::Sqrt(TMath::Power(muon.
pullX(
i-4,2,arbitrationType),2.)+TMath::Power(muon.
pullY(
i-4,2,arbitrationType),2.)) > 1. ) {
245 if(use_match_dist_penalty) {
247 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. ) {
248 station_weight[
i-1] *= 1./TMath::Power(
249 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);
252 station_weight[
i-1] *= 1./TMath::Power(
253 TMath::Sqrt(TMath::Power(muon.
pullX(
i-4,2,arbitrationType),2.)+TMath::Power(muon.
pullY(
i-4,2,arbitrationType),2.)),.25);
266 station_weight[
i-1] = 0.;
270 full_weight += station_weight[
i-1];
276 if( nr_of_stations_crossed == 0 ) {
291 double minCompatibility,
294 bool goodMuon =
false;
300 else goodMuon =
false;
318 double maxChamberDist,
319 double maxChamberDistPull,
321 bool syncMinNMatchesNRequiredStationsInBarrelOnly,
322 bool applyAlsoAngularCuts)
325 bool goodMuon =
false;
330 if(minNumberOfMatches == 0)
return true;
332 unsigned int theStationMask = muon.
stationMask(arbitrationType);
333 unsigned int theRequiredStationMask =
RequiredStationMask(muon, maxChamberDist, maxChamberDistPull, arbitrationType);
337 int numRequiredStations = 0;
338 for(
int it = 0; it < 8; ++it) {
339 if(theStationMask & 1<<it) ++numSegs;
340 if(theRequiredStationMask & 1<<it) ++numRequiredStations;
345 if (syncMinNMatchesNRequiredStationsInBarrelOnly) {
347 if (fabs(muon.
eta()) < 1.2) {
348 if(minNumberOfMatches > numRequiredStations)
349 minNumberOfMatches = numRequiredStations;
350 if(minNumberOfMatches < 1)
351 minNumberOfMatches = 1;
354 if(minNumberOfMatches > numRequiredStations)
355 minNumberOfMatches = numRequiredStations;
356 if(minNumberOfMatches < 1)
357 minNumberOfMatches = 1;
360 if(numSegs >= minNumberOfMatches) goodMuon = 1;
367 if(theRequiredStationMask) {
368 for(
int stationIdx = 7; stationIdx >= 0; --stationIdx)
369 if(theRequiredStationMask & 1<<stationIdx){
370 if(theStationMask & 1<<stationIdx) {
371 lastSegBit = stationIdx;
380 for(
int stationIdx = 7; stationIdx >= 0; --stationIdx)
381 if(theStationMask & 1<<stationIdx) {
382 lastSegBit = stationIdx;
387 if(!goodMuon)
return false;
391 station = lastSegBit < 4 ? lastSegBit+1 : lastSegBit-3;
392 detector = lastSegBit < 4 ? 1 : 2;
395 if(fabs(muon.
pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
396 fabs(muon.
dX(station,detector,arbitrationType)) > maxAbsDx)
399 if(applyAlsoAngularCuts && fabs(muon.
pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX)
403 if (maxAbsDy < 999999) {
407 if(fabs(muon.
pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
408 fabs(muon.
dY(station,2,arbitrationType)) > maxAbsDy)
411 if(applyAlsoAngularCuts && fabs(muon.
pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY)
429 for (
int stationIdx = station; stationIdx > 0; --stationIdx) {
430 if(! (theStationMask & 1<<(stationIdx-1)))
433 if(muon.
dY(stationIdx,1,arbitrationType) > 999998)
436 if(fabs(muon.
pullY(stationIdx,1,arbitrationType,1)) > maxAbsPullY &&
437 fabs(muon.
dY(stationIdx,1,arbitrationType)) > maxAbsDy) {
441 if(applyAlsoAngularCuts && fabs(muon.
pullDyDz(stationIdx,1,arbitrationType,1)) > maxAbsPullY)
464 unsigned int theStationMask = muon.
stationMask(arbitrationType);
467 if (! theStationMask)
return false;
474 bool existsGoodDTSegX =
false;
475 bool existsDTSegY =
false;
479 for(
int stationIdx = 0; stationIdx <= 7; ++stationIdx)
480 if(theStationMask & 1<<stationIdx) {
481 station = stationIdx < 4 ? stationIdx+1 : stationIdx-3;
482 detector = stationIdx < 4 ? 1 : 2;
484 if((fabs(muon.
pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
485 fabs(muon.
dX(station,detector,arbitrationType)) > maxAbsDx) ||
486 (applyAlsoAngularCuts && fabs(muon.
pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX))
488 else if (detector == 1)
489 existsGoodDTSegX =
true;
492 if (maxAbsDy < 999999) {
494 if((fabs(muon.
pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
495 fabs(muon.
dY(station,2,arbitrationType)) > maxAbsDy) ||
496 (applyAlsoAngularCuts && fabs(muon.
pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY))
500 if(muon.
dY(station,1,arbitrationType) > 999998)
505 if((fabs(muon.
pullY(station,1,arbitrationType,1)) > maxAbsPullY &&
506 fabs(muon.
dY(station,1,arbitrationType)) > maxAbsDy) ||
507 (applyAlsoAngularCuts && fabs(muon.
pullDyDz(station,1,arbitrationType,1)) > maxAbsPullY)) {
524 if (maxAbsDy < 999999) {
527 else if (existsGoodDTSegX)
535 if ( minNumberOfMatches == 0 )
return true;
538 for ( std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = muon.
matches().begin();
539 chamberMatch != muon.
matches().end(); ++chamberMatch )
541 if ( chamberMatch->detector() != 3 )
continue;
543 const double trkX = chamberMatch->x;
544 const double errX = chamberMatch->xErr;
546 for ( std::vector<reco::MuonRPCHitMatch>::const_iterator rpcMatch = chamberMatch->rpcMatches.begin();
547 rpcMatch != chamberMatch->rpcMatches.end(); ++rpcMatch )
549 const double rpcX = rpcMatch->x;
551 const double dX = fabs(rpcX-trkX);
552 if ( dX < maxAbsDx
or dX/errX < maxAbsPullX )
560 if ( nMatch >= minNumberOfMatches )
return true;
604 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,
true,
false);
607 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,
true,
false);
610 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,
false,
false);
613 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,
false,
false);
616 if (muon.
pt() < 8. && fabs(muon.
eta()) < 1.2)
617 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,
false,
false);
619 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,
false,
false);
622 if (muon.
pt() < 8. && fabs(muon.
eta()) < 1.2)
623 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,
false,
false);
625 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,
false,
false);
645 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,
false,
true);
648 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,
false,
true);
651 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,
false,
true);
654 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,
false,
true);
657 if (muon.
pt() < 8. && fabs(muon.
eta()) < 1.2)
658 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,
false,
false);
660 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,
true,
false);
663 if (muon.
pt() < 8. && fabs(muon.
eta()) < 1.2)
664 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,
false,
false);
666 return muon.
isTrackerMuon() &&
isGoodMuon(muon,
TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,
true,
false);
669 return muon.
isRPCMuon() &&
isGoodMuon(muon,
RPCMu, 2, 20, 4, 1e9, 1e9, 1e9, 1e9, arbitrationType,
false,
false);
677 double pullX,
double pullY,
bool checkAdjacentChambers)
681 unsigned int betterMuon = ( muon1.
pt() > muon2.
pt() ? 1 : 2 );
682 for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber1 = muon1.
matches().begin();
683 chamber1 != muon1.
matches().end(); ++chamber1 )
684 for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber2 = muon2.
matches().begin();
685 chamber2 != muon2.
matches().end(); ++chamber2 )
692 if ( chamber1->id == chamber2->id ){
694 if ( fabs(chamber1->x-chamber2->x) <
695 pullX *
sqrt(chamber1->xErr*chamber1->xErr+chamber2->xErr*chamber2->xErr) )
697 if ( betterMuon == 1 )
701 if ( nMatches1==0 || nMatches2==0 )
return true;
704 if ( fabs(chamber1->y-chamber2->y) <
705 pullY *
sqrt(chamber1->yErr*chamber1->yErr+chamber2->yErr*chamber2->yErr) )
707 if ( betterMuon == 1 )
711 if ( nMatches1==0 || nMatches2==0 )
return true;
714 if ( ! checkAdjacentChambers )
continue;
722 if ( id1.
ring() != id2.
ring() )
continue;
730 if ( fabs(chamber1->edgeX) > chamber1->xErr*pullX )
continue;
731 if ( fabs(chamber2->edgeX) > chamber2->xErr*pullX )
continue;
732 if ( chamber1->x * chamber2->x < 0 ) {
733 if ( betterMuon == 1 )
737 if ( nMatches1==0 || nMatches2==0 )
return true;
752 bool hits = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
753 muon.
innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
758 return muID && hits && ip;
783 if(!muID)
return false;
785 bool layers = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
786 muon.
innerTrack()->hitPattern().pixelLayersWithMeasurement() > 0;
792 return layers && ip && ishighq;
799 if(!muID)
return false;
801 bool hits = muon.
innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
802 muon.
innerTrack()->hitPattern().numberOfValidPixelHits() > 0;
808 return muID && hits && momQuality && ip;
816 for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.
matches().begin();
817 chamberMatch != mu.
matches().end(); ++chamberMatch) {
818 if (chamberMatch->segmentMatches.empty())
continue;
819 for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch2 = mu2.
matches().begin();
820 chamberMatch2 != mu2.
matches().end(); ++chamberMatch2) {
821 if (chamberMatch2->segmentMatches.empty())
continue;
822 if (chamberMatch2->id() != chamberMatch->id())
continue;
823 for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
824 segmentMatch != chamberMatch->segmentMatches.end(); ++segmentMatch) {
825 if (!segmentMatch->isMask(segmentArbitrationMask))
continue;
826 for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch2 = chamberMatch2->segmentMatches.begin();
827 segmentMatch2 != chamberMatch2->segmentMatches.end(); ++segmentMatch2) {
828 if (!segmentMatch2->isMask(segmentArbitrationMask))
continue;
829 if ((segmentMatch->cscSegmentRef.isNonnull() && segmentMatch->cscSegmentRef == segmentMatch2->cscSegmentRef) ||
830 (segmentMatch-> dtSegmentRef.isNonnull() && segmentMatch-> dtSegmentRef == segmentMatch2-> dtSegmentRef) ) {
float chi2LocalPosition
chi2 value for the STA-TK matching of local position
float segmentX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
tuple ret
prodAgent to be discontinued
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
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
virtual TrackRef innerTrack() const
bool isTrackerMuon() const
bool isGlobalMuon() const
float trkKink
value of the kink algorithm applied to the inner track stub
bool isMatchesValid() const
bool isStandAloneMuon() const
float trkRelChi2
chi2 value for the inner track stub with respect to the global track
bool isMediumMuon(const reco::Muon &)
const Point & position() const
position
float caloCompatibility(const reco::Muon &muon)
bool isLooseMuon(const reco::Muon &)
SelectionType
Selector type.
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
bool isSoftMuon(const reco::Muon &, const reco::Vertex &)
virtual const Track * bestTrack() const
best track pointer
float trackDistErr(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
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
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
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
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)
a lightweight "map" for selection type string label and enum value
unsigned int RequiredStationMask(const reco::Muon &muon, double maxChamberDist, double maxChamberDistPull, reco::Muon::ArbitrationType arbitrationType)
std::vector< MuonChamberMatch > & matches()
get muon matching information
int numberOfMatches(ArbitrationType type=SegmentAndTrackArbitration) const
get number of chambers with matched segments
float pullDyDz(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const
float trackDist(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration) const
virtual double eta() const final
momentum pseudorapidity
bool isTightMuon(const reco::Muon &, const reco::Vertex &)
float caloCompatibility() const
virtual double pt() const final
transverse momentum
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
float pullX(int station, int muonSubdetId, ArbitrationType type=SegmentAndTrackArbitration, bool includeSegmentError=true) const