CMS 3D CMS Logo

TrackInfoBuilder.cc
Go to the documentation of this file.
10 #include "TVector3.h"
11 
12 namespace btagbtvdeep {
13 
14  // adapted from DeepNtuples
16  : builder_(build),
17  trackMomentum_(0),
18  trackEta_(0),
19  trackEtaRel_(0),
20  trackPtRel_(0),
21  trackPPar_(0),
22  trackDeltaR_(0),
23  trackPtRatio_(0),
24  trackPParRatio_(0),
25  trackSip2dVal_(0),
26  trackSip2dSig_(0),
27  trackSip3dVal_(0),
28  trackSip3dSig_(0),
29  trackJetDistVal_(0),
30  trackJetDistSig_(0) {}
31 
33  const math::XYZVector &jetDir,
34  GlobalVector refjetdirection,
35  const reco::Vertex &pv) {
36  TVector3 jetDir3(jetDir.x(), jetDir.y(), jetDir.z());
37 
38  // deal with PAT/AOD polymorphism to get track
39  const reco::Track *track_ptr = nullptr;
40  auto packed_candidate = dynamic_cast<const pat::PackedCandidate *>(candidate);
41  auto pf_candidate = dynamic_cast<const reco::PFCandidate *>(candidate);
42  if (pf_candidate) {
43  track_ptr = pf_candidate->bestTrack(); // trackRef was sometimes null
44  } else if (packed_candidate && packed_candidate->hasTrackDetails()) {
45  // if PackedCandidate does not have TrackDetails this gives an Exception
46  // because unpackCovariance might be called for pseudoTrack/bestTrack
47  track_ptr = &(packed_candidate->pseudoTrack());
48  }
49 
50  if (!track_ptr) {
51  TVector3 trackMom3(candidate->momentum().x(), candidate->momentum().y(), candidate->momentum().z());
52  trackMomentum_ = candidate->p();
53  trackEta_ = candidate->eta();
54  trackEtaRel_ = reco::btau::etaRel(jetDir, candidate->momentum());
55  trackPtRel_ = trackMom3.Perp(jetDir3);
56  trackPPar_ = jetDir.Dot(candidate->momentum());
57  trackDeltaR_ = reco::deltaR(candidate->momentum(), jetDir);
58  trackPtRatio_ = trackMom3.Perp(jetDir3) / candidate->p();
59  trackPParRatio_ = jetDir.Dot(candidate->momentum()) / candidate->p();
60  trackSip2dVal_ = 0.;
61  trackSip2dSig_ = 0.;
62  trackSip3dVal_ = 0.;
63  trackSip3dSig_ = 0.;
64  trackJetDistVal_ = 0.;
65  trackJetDistSig_ = 0.;
66  return;
67  }
68 
69  math::XYZVector trackMom = track_ptr->momentum();
70  double trackMag = std::sqrt(trackMom.Mag2());
71  TVector3 trackMom3(trackMom.x(), trackMom.y(), trackMom.z());
72 
73  trackMomentum_ = std::sqrt(trackMom.Mag2());
74  trackEta_ = trackMom.Eta();
75  trackEtaRel_ = reco::btau::etaRel(jetDir, trackMom);
76  trackPtRel_ = trackMom3.Perp(jetDir3);
77  trackPPar_ = jetDir.Dot(trackMom);
78  trackDeltaR_ = reco::deltaR(trackMom, jetDir);
79  trackPtRatio_ = trackMom3.Perp(jetDir3) / trackMag;
80  trackPParRatio_ = jetDir.Dot(trackMom) / trackMag;
81 
82  reco::TransientTrack transientTrack;
83  transientTrack = builder_->build(*track_ptr);
84  Measurement1D meas_ip2d = IPTools::signedTransverseImpactParameter(transientTrack, refjetdirection, pv).second;
85  Measurement1D meas_ip3d = IPTools::signedImpactParameter3D(transientTrack, refjetdirection, pv).second;
86  Measurement1D jetdist = IPTools::jetTrackDistance(transientTrack, refjetdirection, pv).second;
87  trackSip2dVal_ = static_cast<float>(meas_ip2d.value());
88  trackSip2dSig_ = static_cast<float>(meas_ip2d.significance());
89  trackSip3dVal_ = static_cast<float>(meas_ip3d.value());
90  trackSip3dSig_ = static_cast<float>(meas_ip3d.significance());
91  trackJetDistVal_ = static_cast<float>(jetdist.value());
92  trackJetDistSig_ = static_cast<float>(jetdist.significance());
93  }
94 
95 } // namespace btagbtvdeep
Vector3DBase
Definition: Vector3DBase.h:8
btagbtvdeep
Definition: BoostedDoubleSVTagInfoFeatures.h:4
Measurement1D
Definition: Measurement1D.h:11
btagbtvdeep::TrackInfoBuilder::TrackInfoBuilder
TrackInfoBuilder(edm::ESHandle< TransientTrackBuilder > &build)
Definition: TrackInfoBuilder.cc:15
PFCandidate.h
reco::Candidate::eta
virtual double eta() const =0
momentum pseudorapidity
Measurement1D::value
double value() const
Definition: Measurement1D.h:25
btagbtvdeep::TrackInfoBuilder::trackSip3dVal_
float trackSip3dVal_
Definition: TrackInfoBuilder.h:48
btagbtvdeep::TrackInfoBuilder::trackMomentum_
float trackMomentum_
Definition: TrackInfoBuilder.h:38
btagbtvdeep::TrackInfoBuilder::trackPtRatio_
float trackPtRatio_
Definition: TrackInfoBuilder.h:44
VertexDistance3D.h
IPTools::signedTransverseImpactParameter
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
newFWLiteAna.build
build
Definition: newFWLiteAna.py:126
IPTools::signedImpactParameter3D
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:81
btagbtvdeep::TrackInfoBuilder::buildTrackInfo
void buildTrackInfo(const reco::Candidate *candidate, const math::XYZVector &jetDir, GlobalVector refjetdirection, const reco::Vertex &pv)
Definition: TrackInfoBuilder.cc:32
btagbtvdeep::TrackInfoBuilder::trackSip3dSig_
float trackSip3dSig_
Definition: TrackInfoBuilder.h:49
btagbtvdeep::TrackInfoBuilder::trackPPar_
float trackPPar_
Definition: TrackInfoBuilder.h:42
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
reco::Track
Definition: Track.h:27
edm::ESHandle< TransientTrackBuilder >
reco::Candidate::p
virtual double p() const =0
magnitude of momentum vector
btagbtvdeep::TrackInfoBuilder::trackSip2dVal_
float trackSip2dVal_
Definition: TrackInfoBuilder.h:46
btagbtvdeep::TrackInfoBuilder::builder_
edm::ESHandle< TransientTrackBuilder > builder_
Definition: TrackInfoBuilder.h:36
IPTools::jetTrackDistance
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:206
Measurement1D::significance
double significance() const
Definition: Measurement1D.h:29
TransientTrackBuilder.h
reco::Candidate::momentum
virtual Vector momentum() const =0
spatial momentum vector
btagbtvdeep::TrackInfoBuilder::trackPParRatio_
float trackPParRatio_
Definition: TrackInfoBuilder.h:45
btagbtvdeep::TrackInfoBuilder::trackJetDistSig_
float trackJetDistSig_
Definition: TrackInfoBuilder.h:52
IPTagInfo.h
PackedCandidate.h
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
TrackInfoBuilder.h
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
btagbtvdeep::TrackInfoBuilder::trackPtRel_
float trackPtRel_
Definition: TrackInfoBuilder.h:41
reco::btau::etaRel
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
Definition: TaggingVariable.h:22
btagbtvdeep::TrackInfoBuilder::trackDeltaR_
float trackDeltaR_
Definition: TrackInfoBuilder.h:43
btagbtvdeep::TrackInfoBuilder::trackSip2dSig_
float trackSip2dSig_
Definition: TrackInfoBuilder.h:47
btagbtvdeep::TrackInfoBuilder::trackEtaRel_
float trackEtaRel_
Definition: TrackInfoBuilder.h:40
TransientTrackRecord.h
reco::Candidate
Definition: Candidate.h:27
IPTools.h
reco::TransientTrack
Definition: TransientTrack.h:19
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
TransientTrackBuilder::build
reco::TransientTrack build(const reco::Track *p) const
Definition: TransientTrackBuilder.cc:20
btagbtvdeep::TrackInfoBuilder::trackJetDistVal_
float trackJetDistVal_
Definition: TrackInfoBuilder.h:51
btagbtvdeep::TrackInfoBuilder::trackEta_
float trackEta_
Definition: TrackInfoBuilder.h:39
JetTagInfo.h
reco::Vertex
Definition: Vertex.h:35