#include <RecoBTag/SecondaryVertex/interface/TrackSelector.h>
Public Member Functions | |
bool | operator() (const reco::Track &track, const reco::TrackIPTagInfo::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const |
TrackSelector (const edm::ParameterSet ¶ms) | |
~TrackSelector () | |
Private Attributes | |
double | maxDecayLen |
double | maxDistToAxis |
double | maxJetDeltaR |
double | maxNormChi2 |
unsigned int | minPixelHits |
double | minPt |
unsigned int | minTotalHits |
reco::TrackBase::TrackQuality | quality |
bool | selectQuality |
double | sip2dSigMax |
double | sip2dSigMin |
double | sip2dValMax |
double | sip2dValMin |
double | sip3dSigMax |
double | sip3dSigMin |
double | sip3dValMax |
double | sip3dValMin |
Definition at line 13 of file TrackSelector.h.
TrackSelector::TrackSelector | ( | const edm::ParameterSet & | params | ) |
Definition at line 16 of file TrackSelector.cc.
References edm::ParameterSet::getParameter(), quality, reco::TrackBase::qualityByName(), selectQuality, and reco::TrackBase::undefQuality.
00016 : 00017 minPixelHits(params.getParameter<unsigned int>("pixelHitsMin")), 00018 minTotalHits(params.getParameter<unsigned int>("totalHitsMin")), 00019 minPt(params.getParameter<double>("ptMin")), 00020 maxNormChi2(params.getParameter<double>("normChi2Max")), 00021 maxJetDeltaR(params.getParameter<double>("jetDeltaRMax")), 00022 maxDistToAxis(params.getParameter<double>("maxDistToAxis")), 00023 maxDecayLen(params.getParameter<double>("maxDecayLen")), 00024 sip2dValMin(params.getParameter<double>("sip2dValMin")), 00025 sip2dValMax(params.getParameter<double>("sip2dValMax")), 00026 sip2dSigMin(params.getParameter<double>("sip2dSigMin")), 00027 sip2dSigMax(params.getParameter<double>("sip2dSigMax")), 00028 sip3dValMin(params.getParameter<double>("sip3dValMin")), 00029 sip3dValMax(params.getParameter<double>("sip3dValMax")), 00030 sip3dSigMin(params.getParameter<double>("sip3dSigMin")), 00031 sip3dSigMax(params.getParameter<double>("sip3dSigMax")) 00032 { 00033 std::string qualityClass = 00034 params.getParameter<std::string>("qualityClass"); 00035 if (qualityClass == "any" || qualityClass == "Any" || 00036 qualityClass == "ANY" || qualityClass == "") { 00037 selectQuality = false; 00038 quality = reco::TrackBase::undefQuality; 00039 } else { 00040 selectQuality = true; 00041 quality = reco::TrackBase::qualityByName(qualityClass); 00042 } 00043 }
reco::TrackSelector::~TrackSelector | ( | ) | [inline] |
bool TrackSelector::operator() | ( | const reco::Track & | track, | |
const reco::TrackIPTagInfo::TrackIPData & | ipData, | |||
const reco::Jet & | jet, | |||
const GlobalPoint & | pv | |||
) | const |
Definition at line 46 of file TrackSelector.cc.
References funct::abs(), reco::TrackIPTagInfo::TrackIPData::closestToJetAxis, reco::TrackIPTagInfo::TrackIPData::distanceToJetAxis, reco::TrackBase::hitPattern(), int, reco::TrackIPTagInfo::TrackIPData::ip2d, reco::TrackIPTagInfo::TrackIPData::ip3d, muonGeometry::mag(), maxDecayLen, maxDistToAxis, maxJetDeltaR, maxNormChi2, minPixelHits, minPt, minTotalHits, reco::TrackBase::momentum(), reco::Particle::momentum(), reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), reco::TrackBase::pt(), quality, reco::TrackBase::quality(), selectQuality, Measurement1D::significance(), sip2dSigMax, sip2dSigMin, sip2dValMax, sip2dValMin, sip3dSigMax, sip3dSigMin, sip3dValMax, sip3dValMin, and Measurement1D::value().
00050 { 00051 return (!selectQuality || track.quality(quality)) && 00052 (minPixelHits <= 0 || 00053 track.hitPattern().numberOfValidPixelHits() >= (int)minPixelHits) && 00054 (minTotalHits <= 0 || 00055 track.hitPattern().numberOfValidHits() >= (int)minPixelHits) && 00056 track.pt() >= minPt && 00057 track.normalizedChi2() < maxNormChi2 && 00058 VectorUtil::DeltaR(jet.momentum(), 00059 track.momentum()) < maxJetDeltaR && 00060 std::abs(ipData.distanceToJetAxis) <= maxDistToAxis && 00061 (ipData.closestToJetAxis - pv).mag() <= maxDecayLen && 00062 ipData.ip2d.value() >= sip2dValMin && 00063 ipData.ip2d.value() <= sip2dValMax && 00064 ipData.ip2d.significance() >= sip2dSigMin && 00065 ipData.ip2d.significance() <= sip2dSigMax && 00066 ipData.ip3d.value() >= sip3dValMin && 00067 ipData.ip3d.value() <= sip3dValMax && 00068 ipData.ip3d.significance() >= sip3dSigMin && 00069 ipData.ip3d.significance() <= sip3dSigMax; 00070 }
double reco::TrackSelector::maxDecayLen [private] |
double reco::TrackSelector::maxDistToAxis [private] |
double reco::TrackSelector::maxJetDeltaR [private] |
double reco::TrackSelector::maxNormChi2 [private] |
unsigned int reco::TrackSelector::minPixelHits [private] |
double reco::TrackSelector::minPt [private] |
unsigned int reco::TrackSelector::minTotalHits [private] |
bool reco::TrackSelector::selectQuality [private] |
double reco::TrackSelector::sip2dSigMax [private] |
double reco::TrackSelector::sip2dSigMin [private] |
double reco::TrackSelector::sip2dValMax [private] |
double reco::TrackSelector::sip2dValMin [private] |
double reco::TrackSelector::sip3dSigMax [private] |
double reco::TrackSelector::sip3dSigMin [private] |
double reco::TrackSelector::sip3dValMax [private] |
double reco::TrackSelector::sip3dValMin [private] |