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  useVariableJTA_(params.existsAs<bool>("useVariableJTA") ? params.getParameter<bool>("useVariableJTA") : false)
33 {
34  std::string qualityClass =
35  params.getParameter<std::string>("qualityClass");
36  if (qualityClass == "any" || qualityClass == "Any" ||
37  qualityClass == "ANY" || qualityClass == "") {
38  selectQuality = false;
40  } else {
41  selectQuality = true;
43  }
44  if (useVariableJTA_){
45  varJTApars = {
46  params.getParameter<double>("a_dR"),
47  params.getParameter<double>("b_dR"),
48  params.getParameter<double>("a_pT"),
49  params.getParameter<double>("b_pT"),
50  params.getParameter<double>("min_pT"),
51  params.getParameter<double>("max_pT"),
52  params.getParameter<double>("min_pT_dRcut"),
53  params.getParameter<double>("max_pT_dRcut"),
54  params.getParameter<double>("max_pT_trackPTcut") };
55  }
56 }
57 
58 bool
60  const btag::TrackIPData &ipData,
61  const Jet &jet,
62  const GlobalPoint &pv) const
63 {
64 
65 
66  bool jtaPassed = false;
67  if (useVariableJTA_) {
69  jet.pt(),track.pt(),
70  VectorUtil::DeltaR(jet.momentum(),track.momentum()));
71  }
72  else jtaPassed = true;
73 
74  return (!selectQuality || track.quality(quality)) &&
75  (minPixelHits <= 0 ||
76  track.hitPattern().numberOfValidPixelHits() >= (int)minPixelHits) &&
77  (minTotalHits <= 0 ||
78  track.hitPattern().numberOfValidHits() >= (int)minTotalHits) &&
79  track.pt() >= minPt &&
80  track.normalizedChi2() < maxNormChi2 &&
82  track.momentum()) < maxJetDeltaR &&
83  jtaPassed &&
85  (ipData.closestToJetAxis - pv).mag() <= maxDecayLen &&
86  ipData.ip2d.value() >= sip2dValMin &&
87  ipData.ip2d.value() <= sip2dValMax &&
88  ipData.ip2d.significance() >= sip2dSigMin &&
89  ipData.ip2d.significance() <= sip2dSigMax &&
90  ipData.ip3d.value() >= sip3dValMin &&
91  ipData.ip3d.value() <= sip3dValMax &&
92  ipData.ip3d.significance() >= sip3dSigMin &&
93  ipData.ip3d.significance() <= sip3dSigMax;
94 }
T getParameter(std::string const &) const
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:548
int numberOfValidHits() const
Definition: HitPattern.h:801
Base class for all types of Jets.
Definition: Jet.h:20
Measurement1D ip2d
Definition: IPTagInfo.h:31
virtual Vector momentum() const
spatial momentum vector
reco::TrackBase::TrackQuality quality
Definition: TrackSelector.h:25
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:662
bool operator()(const reco::Track &track, const reco::btag::TrackIPData &ipData, const reco::Jet &jet, const GlobalPoint &pv) const
virtual double pt() const
transverse momentum
TrackSelector(const edm::ParameterSet &params)
reco::btag::variableJTAParameters varJTApars
Definition: TrackSelector.h:42
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:608
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int minPixelHits
Definition: TrackSelector.h:26
Measurement1D ip3d
Definition: IPTagInfo.h:32
double significance() const
Definition: Measurement1D.h:32
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:123
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:437
double value() const
Definition: Measurement1D.h:28
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:497
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:33
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:17
int numberOfValidPixelHits() const
Definition: HitPattern.h:816
volatile std::atomic< bool > shutdown_flag false
unsigned int minTotalHits
Definition: TrackSelector.h:27