#include <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 |
bool | useVariableJTA_ |
reco::TrackIPTagInfo::variableJTAParameters | varJTApars |
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, AlCaHLTBitMon_QueryRunRegistry::string, reco::TrackBase::undefQuality, useVariableJTA_, and varJTApars.
: minPixelHits(params.getParameter<unsigned int>("pixelHitsMin")), minTotalHits(params.getParameter<unsigned int>("totalHitsMin")), minPt(params.getParameter<double>("ptMin")), maxNormChi2(params.getParameter<double>("normChi2Max")), maxJetDeltaR(params.getParameter<double>("jetDeltaRMax")), maxDistToAxis(params.getParameter<double>("maxDistToAxis")), maxDecayLen(params.getParameter<double>("maxDecayLen")), sip2dValMin(params.getParameter<double>("sip2dValMin")), sip2dValMax(params.getParameter<double>("sip2dValMax")), sip2dSigMin(params.getParameter<double>("sip2dSigMin")), sip2dSigMax(params.getParameter<double>("sip2dSigMax")), sip3dValMin(params.getParameter<double>("sip3dValMin")), sip3dValMax(params.getParameter<double>("sip3dValMax")), sip3dSigMin(params.getParameter<double>("sip3dSigMin")), sip3dSigMax(params.getParameter<double>("sip3dSigMax")), useVariableJTA_(params.existsAs<bool>("useVariableJTA") ? params.getParameter<bool>("useVariableJTA") : false) { std::string qualityClass = params.getParameter<std::string>("qualityClass"); if (qualityClass == "any" || qualityClass == "Any" || qualityClass == "ANY" || qualityClass == "") { selectQuality = false; quality = reco::TrackBase::undefQuality; } else { selectQuality = true; quality = reco::TrackBase::qualityByName(qualityClass); } if (useVariableJTA_){ varJTApars = { params.getParameter<double>("a_dR"), params.getParameter<double>("b_dR"), params.getParameter<double>("a_pT"), params.getParameter<double>("b_pT"), params.getParameter<double>("min_pT"), params.getParameter<double>("max_pT"), params.getParameter<double>("min_pT_dRcut"), params.getParameter<double>("max_pT_dRcut"), params.getParameter<double>("max_pT_trackPTcut") }; } }
reco::TrackSelector::~TrackSelector | ( | ) | [inline] |
Definition at line 16 of file TrackSelector.h.
{}
bool TrackSelector::operator() | ( | const reco::Track & | track, |
const reco::TrackIPTagInfo::TrackIPData & | ipData, | ||
const reco::Jet & | jet, | ||
const GlobalPoint & | pv | ||
) | const |
Definition at line 59 of file TrackSelector.cc.
References abs, reco::TrackIPTagInfo::TrackIPData::closestToJetAxis, reco::TrackIPTagInfo::TrackIPData::distanceToJetAxis, reco::TrackBase::hitPattern(), reco::TrackIPTagInfo::TrackIPData::ip2d, reco::TrackIPTagInfo::TrackIPData::ip3d, maxDecayLen, maxDistToAxis, maxJetDeltaR, maxNormChi2, minPixelHits, minPt, minTotalHits, reco::LeafCandidate::momentum(), reco::TrackBase::momentum(), reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), reco::TrackIPTagInfo::passVariableJTA(), reco::TrackBase::pt(), reco::LeafCandidate::pt(), quality, reco::TrackBase::quality(), selectQuality, Measurement1D::significance(), sip2dSigMax, sip2dSigMin, sip2dValMax, sip2dValMin, sip3dSigMax, sip3dSigMin, sip3dValMax, sip3dValMin, useVariableJTA_, Measurement1D::value(), and varJTApars.
{ bool jtaPassed = false; if (useVariableJTA_) { jtaPassed = TrackIPTagInfo::passVariableJTA( varJTApars, jet.pt(),track.pt(), VectorUtil::DeltaR(jet.momentum(),track.momentum())); } else jtaPassed = true; return (!selectQuality || track.quality(quality)) && (minPixelHits <= 0 || track.hitPattern().numberOfValidPixelHits() >= (int)minPixelHits) && (minTotalHits <= 0 || track.hitPattern().numberOfValidHits() >= (int)minTotalHits) && track.pt() >= minPt && track.normalizedChi2() < maxNormChi2 && VectorUtil::DeltaR(jet.momentum(), track.momentum()) < maxJetDeltaR && jtaPassed && std::abs(ipData.distanceToJetAxis.value()) <= maxDistToAxis && (ipData.closestToJetAxis - pv).mag() <= maxDecayLen && ipData.ip2d.value() >= sip2dValMin && ipData.ip2d.value() <= sip2dValMax && ipData.ip2d.significance() >= sip2dSigMin && ipData.ip2d.significance() <= sip2dSigMax && ipData.ip3d.value() >= sip3dValMin && ipData.ip3d.value() <= sip3dValMax && ipData.ip3d.significance() >= sip3dSigMin && ipData.ip3d.significance() <= sip3dSigMax; }
double reco::TrackSelector::maxDecayLen [private] |
Definition at line 32 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::maxDistToAxis [private] |
Definition at line 31 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::maxJetDeltaR [private] |
Definition at line 30 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::maxNormChi2 [private] |
Definition at line 29 of file TrackSelector.h.
Referenced by operator()().
unsigned int reco::TrackSelector::minPixelHits [private] |
Definition at line 26 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::minPt [private] |
Definition at line 28 of file TrackSelector.h.
Referenced by operator()().
unsigned int reco::TrackSelector::minTotalHits [private] |
Definition at line 27 of file TrackSelector.h.
Referenced by operator()().
Definition at line 25 of file TrackSelector.h.
Referenced by operator()(), and TrackSelector().
bool reco::TrackSelector::selectQuality [private] |
Definition at line 24 of file TrackSelector.h.
Referenced by operator()(), and TrackSelector().
double reco::TrackSelector::sip2dSigMax [private] |
Definition at line 36 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip2dSigMin [private] |
Definition at line 35 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip2dValMax [private] |
Definition at line 34 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip2dValMin [private] |
Definition at line 33 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip3dSigMax [private] |
Definition at line 40 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip3dSigMin [private] |
Definition at line 39 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip3dValMax [private] |
Definition at line 38 of file TrackSelector.h.
Referenced by operator()().
double reco::TrackSelector::sip3dValMin [private] |
Definition at line 37 of file TrackSelector.h.
Referenced by operator()().
bool reco::TrackSelector::useVariableJTA_ [private] |
Definition at line 41 of file TrackSelector.h.
Referenced by operator()(), and TrackSelector().
Definition at line 42 of file TrackSelector.h.
Referenced by operator()(), and TrackSelector().