CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackSelector.cc
Go to the documentation of this file.
1 #include <Math/GenVector/VectorUtil.h>
2 
4 
10 
12 
13 using namespace reco;
14 using namespace ROOT::Math;
15 
17  minPixelHits(params.getParameter<unsigned int>("pixelHitsMin")),
18  minTotalHits(params.getParameter<unsigned int>("totalHitsMin")),
19  minPt(params.getParameter<double>("ptMin")),
20  maxNormChi2(params.getParameter<double>("normChi2Max")),
21  maxJetDeltaR(params.getParameter<double>("jetDeltaRMax")),
22  maxDistToAxis(params.getParameter<double>("maxDistToAxis")),
23  maxDecayLen(params.getParameter<double>("maxDecayLen")),
24  sip2dValMin(params.getParameter<double>("sip2dValMin")),
25  sip2dValMax(params.getParameter<double>("sip2dValMax")),
26  sip2dSigMin(params.getParameter<double>("sip2dSigMin")),
27  sip2dSigMax(params.getParameter<double>("sip2dSigMax")),
28  sip3dValMin(params.getParameter<double>("sip3dValMin")),
29  sip3dValMax(params.getParameter<double>("sip3dValMax")),
30  sip3dSigMin(params.getParameter<double>("sip3dSigMin")),
31  sip3dSigMax(params.getParameter<double>("sip3dSigMax"))
32 {
33  std::string qualityClass =
34  params.getParameter<std::string>("qualityClass");
35  if (qualityClass == "any" || qualityClass == "Any" ||
36  qualityClass == "ANY" || qualityClass == "") {
37  selectQuality = false;
39  } else {
40  selectQuality = true;
42  }
43 }
44 
45 bool
47  const TrackIPTagInfo::TrackIPData &ipData,
48  const Jet &jet,
49  const GlobalPoint &pv) const
50 {
51  return (!selectQuality || track.quality(quality)) &&
52  (minPixelHits <= 0 ||
53  track.hitPattern().numberOfValidPixelHits() >= (int)minPixelHits) &&
54  (minTotalHits <= 0 ||
55  track.hitPattern().numberOfValidHits() >= (int)minTotalHits) &&
56  track.pt() >= minPt &&
57  track.normalizedChi2() < maxNormChi2 &&
58  VectorUtil::DeltaR(jet.momentum(),
59  track.momentum()) < maxJetDeltaR &&
61  (ipData.closestToJetAxis - pv).mag() <= maxDecayLen &&
62  ipData.ip2d.value() >= sip2dValMin &&
63  ipData.ip2d.value() <= sip2dValMax &&
64  ipData.ip2d.significance() >= sip2dSigMin &&
65  ipData.ip2d.significance() <= sip2dSigMax &&
66  ipData.ip3d.value() >= sip3dValMin &&
67  ipData.ip3d.value() <= sip3dValMax &&
68  ipData.ip3d.significance() >= sip3dSigMin &&
69  ipData.ip3d.significance() <= sip3dSigMax;
70 }
T getParameter(std::string const &) const
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:150
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:111
int numberOfValidHits() const
Definition: HitPattern.h:554
Base class for all types of Jets.
Definition: Jet.h:21
virtual Vector momentum() const
spatial momentum vector
#define abs(x)
Definition: mlp_lapack.h:159
reco::TrackBase::TrackQuality quality
Definition: TrackSelector.h:25
TrackSelector(const edm::ParameterSet &params)
double pt() const
track transverse momentum
Definition: TrackBase.h:131
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
unsigned int minPixelHits
Definition: TrackSelector.h:26
double significance() const
Definition: Measurement1D.h:32
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
double value() const
Definition: Measurement1D.h:28
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:377
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
bool operator()(const reco::Track &track, const reco::TrackIPTagInfo::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
int numberOfValidPixelHits() const
Definition: HitPattern.h:566
unsigned int minTotalHits
Definition: TrackSelector.h:27