CMS 3D CMS Logo

SeedingTrackInfoBuilder.cc
Go to the documentation of this file.
8 
9 namespace btagbtvdeep {
10 
12  : pt_(0),
13  eta_(0),
14  phi_(0),
15  mass_(0),
16  dz_(0),
17  dxy_(0),
18  ip3D_(0),
19  sip3D_(0),
20  ip2D_(0),
21  sip2D_(0),
22  ip3D_signed_(0),
23  sip3D_signed_(0),
24  ip2D_signed_(0),
25  sip2D_signed_(0),
26  chi2reduced_(0),
27  nPixelHits_(0),
28  nHits_(0),
29  jetAxisDistance_(0),
30  jetAxisDlength_(0),
31  trackProbability3D_(0),
32  trackProbability2D_(0) {}
33 
35  const reco::Vertex& pv,
36  const reco::Jet& jet, /*GlobalVector jetdirection,*/
37  float mass,
38  const std::pair<bool, Measurement1D>& ip,
39  const std::pair<bool, Measurement1D>& ip2d,
40  float jet_distance,
41  float jaxis_dlength,
42  HistogramProbabilityEstimator* m_probabilityEstimator,
43  bool m_computeProbabilities = false) {
44  GlobalPoint pvp(pv.x(), pv.y(), pv.z());
45  GlobalVector jetdirection(jet.px(), jet.py(), jet.pz());
46 
47  auto const& aTrack = it->track();
48 
49  pt_ = aTrack.pt();
50  eta_ = aTrack.eta();
51  phi_ = aTrack.phi();
52  dz_ = aTrack.dz(pv.position());
53  dxy_ = aTrack.dxy(pv.position());
54  mass_ = mass;
55 
56  std::pair<bool, Measurement1D> ipSigned = IPTools::signedImpactParameter3D(*it, jetdirection, pv);
57  std::pair<bool, Measurement1D> ip2dSigned = IPTools::signedTransverseImpactParameter(*it, jetdirection, pv);
58 
59  ip3D_ = ip.second.value();
60  sip3D_ = ip.second.significance();
61  ip2D_ = ip2d.second.value();
62  sip2D_ = ip2d.second.significance();
63  ip3D_signed_ = ipSigned.second.value();
64  sip3D_signed_ = ipSigned.second.significance();
65  ip2D_signed_ = ip2dSigned.second.value();
66  sip2D_signed_ = ip2dSigned.second.significance();
67 
68  chi2reduced_ = aTrack.normalizedChi2();
69  nPixelHits_ = aTrack.hitPattern().numberOfValidPixelHits();
70  nHits_ = aTrack.hitPattern().numberOfValidHits();
71 
72  jetAxisDistance_ = std::fabs(jet_distance);
73  jetAxisDlength_ = jaxis_dlength;
74 
75  trackProbability3D_ = 0.5;
76  trackProbability2D_ = 0.5;
77 
78  if (m_computeProbabilities) {
79  //probability with 3D ip
80  std::pair<bool, double> probability =
81  m_probabilityEstimator->probability(false, 0, ip.second.significance(), aTrack, jet, pv);
82  double prob3D = (probability.first ? probability.second : -1.);
83 
84  //probability with 2D ip
85  probability = m_probabilityEstimator->probability(false, 1, ip2d.second.significance(), aTrack, jet, pv);
86  double prob2D = (probability.first ? probability.second : -1.);
87 
88  trackProbability3D_ = prob3D;
89  trackProbability2D_ = prob2D;
90  }
91 
93  trackProbability3D_ = 0.5;
95  trackProbability2D_ = 0.5;
96  }
97 
98 } // namespace btagbtvdeep
Base class for all types of Jets.
Definition: Jet.h:20
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:81
constexpr bool isFinite(T x)
def pv(vc)
Definition: MetAnalyzer.py:7
std::pair< bool, double > probability(bool quality, int ipType, float significance, const reco::Track &, const reco::Jet &, const reco::Vertex &) const
void buildSeedingTrackInfo(const reco::TransientTrack *it, const reco::Vertex &pv, const reco::Jet &jet, float mass, const std::pair< bool, Measurement1D > &ip, const std::pair< bool, Measurement1D > &ip2d, float jet_distance, float jaxis_dlength, HistogramProbabilityEstimator *m_probabilityEstimator, bool m_computeProbabilities)