CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
reco::TrackSelector Class Reference

#include <TrackSelector.h>

Public Member Functions

bool operator() (const reco::Track &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
 
bool operator() (const reco::CandidatePtr &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
 
bool operator() (const reco::TrackRef &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
 
 TrackSelector (const edm::ParameterSet &params)
 
 ~TrackSelector ()
 

Private Member Functions

bool trackSelection (const reco::Track &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
 

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::btag::variableJTAParameters varJTApars
 

Detailed Description

Definition at line 14 of file TrackSelector.h.

Constructor & Destructor Documentation

TrackSelector::TrackSelector ( const edm::ParameterSet params)

Definition at line 11 of file TrackSelector.cc.

References edm::ParameterSet::getParameter(), quality, reco::TrackBase::qualityByName(), trackPseudoSelection_cff::qualityClass, selectQuality, AlCaHLTBitMon_QueryRunRegistry::string, reco::TrackBase::undefQuality, useVariableJTA_, and varJTApars.

11  :
12  minPixelHits(params.getParameter<unsigned int>("pixelHitsMin")),
13  minTotalHits(params.getParameter<unsigned int>("totalHitsMin")),
14  minPt(params.getParameter<double>("ptMin")),
15  maxNormChi2(params.getParameter<double>("normChi2Max")),
16  maxJetDeltaR(params.getParameter<double>("jetDeltaRMax")),
17  maxDistToAxis(params.getParameter<double>("maxDistToAxis")),
18  maxDecayLen(params.getParameter<double>("maxDecayLen")),
19  sip2dValMin(params.getParameter<double>("sip2dValMin")),
20  sip2dValMax(params.getParameter<double>("sip2dValMax")),
21  sip2dSigMin(params.getParameter<double>("sip2dSigMin")),
22  sip2dSigMax(params.getParameter<double>("sip2dSigMax")),
23  sip3dValMin(params.getParameter<double>("sip3dValMin")),
24  sip3dValMax(params.getParameter<double>("sip3dValMax")),
25  sip3dSigMin(params.getParameter<double>("sip3dSigMin")),
26  sip3dSigMax(params.getParameter<double>("sip3dSigMax")),
27  useVariableJTA_(params.existsAs<bool>("useVariableJTA") ? params.getParameter<bool>("useVariableJTA") : false)
28 {
30  params.getParameter<std::string>("qualityClass");
31  if (qualityClass == "any" || qualityClass == "Any" ||
32  qualityClass == "ANY" || qualityClass == "") {
33  selectQuality = false;
35  } else {
36  selectQuality = true;
38  }
39  if (useVariableJTA_){
40  varJTApars = {
41  params.getParameter<double>("a_dR"),
42  params.getParameter<double>("b_dR"),
43  params.getParameter<double>("a_pT"),
44  params.getParameter<double>("b_pT"),
45  params.getParameter<double>("min_pT"),
46  params.getParameter<double>("max_pT"),
47  params.getParameter<double>("min_pT_dRcut"),
48  params.getParameter<double>("max_pT_dRcut"),
49  params.getParameter<double>("max_pT_trackPTcut") };
50  }
51 }
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
reco::TrackBase::TrackQuality quality
Definition: TrackSelector.h:43
reco::btag::variableJTAParameters varJTApars
Definition: TrackSelector.h:60
unsigned int minPixelHits
Definition: TrackSelector.h:44
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
unsigned int minTotalHits
Definition: TrackSelector.h:45
reco::TrackSelector::~TrackSelector ( )
inline

Member Function Documentation

bool TrackSelector::operator() ( const reco::Track track,
const reco::btag::TrackIPData ipData,
const reco::Jet jet,
const GlobalPoint pv 
) const

Definition at line 54 of file TrackSelector.cc.

References reco::deltaR(), reco::deltaR2(), maxJetDeltaR, minPt, reco::LeafCandidate::momentum(), reco::TrackBase::momentum(), reco::IPTagInfo< Container, Base >::passVariableJTA(), reco::LeafCandidate::pt(), reco::TrackBase::pt(), trackSelection(), useVariableJTA_, and varJTApars.

Referenced by ~TrackSelector().

58 {
59 
60 
61  bool jtaPassed = false;
62  if (useVariableJTA_) {
64  jet.pt(),track.pt(),
65  reco::deltaR(jet.momentum(),track.momentum()));
66  }
67  else jtaPassed = true;
68 
69  return track.pt() >= minPt &&
70  reco::deltaR2(jet.momentum(),
71  track.momentum()) < maxJetDeltaR*maxJetDeltaR &&
72  jtaPassed &&
73  trackSelection(track, ipData, jet, pv);
74 }
virtual double pt() const final
transverse momentum
bool trackSelection(const reco::Track &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
virtual Vector momentum() const final
spatial momentum vector
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:670
reco::btag::variableJTAParameters varJTApars
Definition: TrackSelector.h:60
static bool passVariableJTA(const btag::variableJTAParameters &params, double jetpt, double trackpt, double jettrackdr)
Definition: IPTagInfo.h:303
double pt() const
track transverse momentum
Definition: TrackBase.h:616
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
bool TrackSelector::operator() ( const reco::CandidatePtr track,
const reco::btag::TrackIPData ipData,
const reco::Jet jet,
const GlobalPoint pv 
) const

Definition at line 77 of file TrackSelector.cc.

References reco::deltaR(), reco::deltaR2(), metsig::jet, maxJetDeltaR, minPt, reco::Candidate::momentum(), reco::LeafCandidate::momentum(), reco::IPTagInfo< Container, Base >::passVariableJTA(), reco::Candidate::pt(), reco::LeafCandidate::pt(), MetAnalyzer::pv(), reco::btag::toTrack(), trackSelection(), useVariableJTA_, and varJTApars.

81 {
82 
83 
84  bool jtaPassed = false;
85  if (useVariableJTA_) {
87  jet.pt(),track->pt(),
88  reco::deltaR(jet.momentum(),track->momentum()));
89  }
90  else jtaPassed = true;
91 
92  return track->pt() >= minPt &&
93  reco::deltaR2(jet.momentum(),
94  track->momentum()) < maxJetDeltaR*maxJetDeltaR &&
95  jtaPassed &&
96  trackSelection(*reco::btag::toTrack(track), ipData, jet, pv);
97 }
virtual double pt() const final
transverse momentum
bool trackSelection(const reco::Track &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
virtual Vector momentum() const final
spatial momentum vector
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
reco::btag::variableJTAParameters varJTApars
Definition: TrackSelector.h:60
static bool passVariableJTA(const btag::variableJTAParameters &params, double jetpt, double trackpt, double jettrackdr)
Definition: IPTagInfo.h:303
def pv(vc)
Definition: MetAnalyzer.py:6
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
virtual double pt() const =0
transverse momentum
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
virtual Vector momentum() const =0
spatial momentum vector
bool reco::TrackSelector::operator() ( const reco::TrackRef track,
const reco::btag::TrackIPData ipData,
const reco::Jet jet,
const GlobalPoint pv 
) const
inline
bool TrackSelector::trackSelection ( const reco::Track track,
const reco::btag::TrackIPData ipData,
const reco::Jet jet,
const GlobalPoint pv 
) const
private

Definition at line 100 of file TrackSelector.cc.

References funct::abs(), reco::btag::TrackIPData::closestToJetAxis, reco::btag::TrackIPData::distanceToJetAxis, reco::TrackBase::hitPattern(), createfilelist::int, reco::btag::TrackIPData::ip2d, reco::btag::TrackIPData::ip3d, maxDecayLen, maxDistToAxis, maxNormChi2, minPixelHits, minTotalHits, reco::TrackBase::normalizedChi2(), reco::HitPattern::numberOfValidHits(), reco::HitPattern::numberOfValidPixelHits(), quality, reco::TrackBase::quality(), selectQuality, Measurement1D::significance(), sip2dSigMax, sip2dSigMin, sip2dValMax, sip2dValMin, sip3dSigMax, sip3dSigMin, sip3dValMax, sip3dValMin, and Measurement1D::value().

Referenced by operator()().

104 {
105 
106  return (!selectQuality || track.quality(quality)) &&
107  (minPixelHits <= 0 ||
109  (minTotalHits <= 0 ||
110  track.hitPattern().numberOfValidHits() >= (int)minTotalHits) &&
111  track.normalizedChi2() < maxNormChi2 &&
113  (ipData.closestToJetAxis - pv).mag() <= maxDecayLen &&
114  ipData.ip2d.value() >= sip2dValMin &&
115  ipData.ip2d.value() <= sip2dValMax &&
116  ipData.ip2d.significance() >= sip2dSigMin &&
117  ipData.ip2d.significance() <= sip2dSigMax &&
118  ipData.ip3d.value() >= sip3dValMin &&
119  ipData.ip3d.value() <= sip3dValMax &&
120  ipData.ip3d.significance() >= sip3dSigMin &&
121  ipData.ip3d.significance() <= sip3dSigMax;
122 }
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:556
int numberOfValidHits() const
Definition: HitPattern.h:823
Measurement1D ip2d
Definition: IPTagInfo.h:31
reco::TrackBase::TrackQuality quality
Definition: TrackSelector.h:43
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int minPixelHits
Definition: TrackSelector.h:44
Measurement1D ip3d
Definition: IPTagInfo.h:32
double significance() const
Definition: Measurement1D.h:32
GlobalPoint closestToJetAxis
Definition: IPTagInfo.h:29
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
double value() const
Definition: Measurement1D.h:28
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:505
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:33
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
unsigned int minTotalHits
Definition: TrackSelector.h:45

Member Data Documentation

double reco::TrackSelector::maxDecayLen
private

Definition at line 50 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::maxDistToAxis
private

Definition at line 49 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::maxJetDeltaR
private

Definition at line 48 of file TrackSelector.h.

Referenced by operator()().

double reco::TrackSelector::maxNormChi2
private

Definition at line 47 of file TrackSelector.h.

Referenced by trackSelection().

unsigned int reco::TrackSelector::minPixelHits
private

Definition at line 44 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::minPt
private

Definition at line 46 of file TrackSelector.h.

Referenced by operator()().

unsigned int reco::TrackSelector::minTotalHits
private

Definition at line 45 of file TrackSelector.h.

Referenced by trackSelection().

reco::TrackBase::TrackQuality reco::TrackSelector::quality
private

Definition at line 43 of file TrackSelector.h.

Referenced by trackSelection(), and TrackSelector().

bool reco::TrackSelector::selectQuality
private

Definition at line 42 of file TrackSelector.h.

Referenced by trackSelection(), and TrackSelector().

double reco::TrackSelector::sip2dSigMax
private

Definition at line 54 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip2dSigMin
private

Definition at line 53 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip2dValMax
private

Definition at line 52 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip2dValMin
private

Definition at line 51 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip3dSigMax
private

Definition at line 58 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip3dSigMin
private

Definition at line 57 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip3dValMax
private

Definition at line 56 of file TrackSelector.h.

Referenced by trackSelection().

double reco::TrackSelector::sip3dValMin
private

Definition at line 55 of file TrackSelector.h.

Referenced by trackSelection().

bool reco::TrackSelector::useVariableJTA_
private

Definition at line 59 of file TrackSelector.h.

Referenced by operator()(), and TrackSelector().

reco::btag::variableJTAParameters reco::TrackSelector::varJTApars
private

Definition at line 60 of file TrackSelector.h.

Referenced by operator()(), and TrackSelector().