CMS 3D CMS Logo

TrackPairInfoBuilder.cc
Go to the documentation of this file.
6 
7 namespace btagbtvdeep {
8 
10  :
11 
12  track_pt_(0),
13  track_eta_(0),
14  track_phi_(0),
15  track_dz_(0),
16  track_dxy_(0),
17 
18  pca_distance_(0),
19  pca_significance_(0),
20 
21  pcaSeed_x_(0),
22  pcaSeed_y_(0),
23  pcaSeed_z_(0),
24  pcaSeed_xerr_(0),
25  pcaSeed_yerr_(0),
26  pcaSeed_zerr_(0),
27  pcaTrack_x_(0),
28  pcaTrack_y_(0),
29  pcaTrack_z_(0),
30  pcaTrack_xerr_(0),
31  pcaTrack_yerr_(0),
32  pcaTrack_zerr_(0),
33 
34  dotprodTrack_(0),
35  dotprodSeed_(0),
36  pcaSeed_dist_(0),
37  pcaTrack_dist_(0),
38 
39  track_candMass_(0),
40  track_ip2d_(0),
41  track_ip2dSig_(0),
42  track_ip3d_(0),
43  track_ip3dSig_(0),
44 
45  dotprodTrackSeed2D_(0),
46  dotprodTrackSeed2DV_(0),
47  dotprodTrackSeed3D_(0),
48  dotprodTrackSeed3DV_(0),
49 
50  pca_jetAxis_dist_(0),
51  pca_jetAxis_dotprod_(0),
52  pca_jetAxis_dEta_(0),
53  pca_jetAxis_dPhi_(0)
54 
55  {}
56 
58  const reco::TransientTrack* tt,
59  const reco::Vertex& pv,
60  float mass,
61  GlobalVector jetdirection,
62  const std::pair<bool, Measurement1D>& t_ip,
63  const std::pair<bool, Measurement1D>& t_ip2d) {
64  GlobalPoint pvp(pv.x(), pv.y(), pv.z());
65 
66  VertexDistance3D distanceComputer;
68 
69  auto const& iImpactState = it->impactPointState();
70  auto const& tImpactState = tt->impactPointState();
71 
72  if (dist.calculate(tImpactState, iImpactState)) {
73  GlobalPoint ttPoint = dist.points().first;
74  GlobalError ttPointErr = tImpactState.cartesianError().position();
75  GlobalPoint seedPosition = dist.points().second;
76  GlobalError seedPositionErr = iImpactState.cartesianError().position();
77 
79  distanceComputer.distance(VertexState(seedPosition, seedPositionErr), VertexState(ttPoint, ttPointErr));
80 
82 
83  GlobalVector pairMomentum((Basic3DVector<float>)(it->track().momentum() + tt->track().momentum()));
84  GlobalVector pvToPCA(cp - pvp);
85 
86  float pvToPCAseed = (seedPosition - pvp).mag();
87  float pvToPCAtrack = (ttPoint - pvp).mag();
88  float distance = dist.distance();
89 
90  GlobalVector trackDir2D(tImpactState.globalDirection().x(), tImpactState.globalDirection().y(), 0.);
91  GlobalVector seedDir2D(iImpactState.globalDirection().x(), iImpactState.globalDirection().y(), 0.);
92  GlobalVector trackPCADir2D(ttPoint.x() - pvp.x(), ttPoint.y() - pvp.y(), 0.);
93  GlobalVector seedPCADir2D(seedPosition.x() - pvp.x(), seedPosition.y() - pvp.y(), 0.);
94 
95  float dotprodTrack = (ttPoint - pvp).unit().dot(tImpactState.globalDirection().unit());
96  float dotprodSeed = (seedPosition - pvp).unit().dot(iImpactState.globalDirection().unit());
97 
99  Line::DirectionType dir(jetdirection);
100  Line::DirectionType pairMomentumDir(pairMomentum);
101  Line jetLine(pos, dir);
102  Line PCAMomentumLine(cp, pairMomentumDir);
103 
104  track_pt_ = tt->track().pt();
105  track_eta_ = tt->track().eta();
106  track_phi_ = tt->track().phi();
107  track_dz_ = tt->track().dz(pv.position());
108  track_dxy_ = tt->track().dxy(pv.position());
109 
111  pca_significance_ = m.significance();
112 
113  pcaSeed_x_ = seedPosition.x();
114  pcaSeed_y_ = seedPosition.y();
115  pcaSeed_z_ = seedPosition.z();
116  pcaSeed_xerr_ = seedPositionErr.cxx();
117  pcaSeed_yerr_ = seedPositionErr.cyy();
118  pcaSeed_zerr_ = seedPositionErr.czz();
119  pcaTrack_x_ = ttPoint.x();
120  pcaTrack_y_ = ttPoint.y();
121  pcaTrack_z_ = ttPoint.z();
122  pcaTrack_xerr_ = ttPointErr.cxx();
123  pcaTrack_yerr_ = ttPointErr.cyy();
124  pcaTrack_zerr_ = ttPointErr.czz();
125 
128  pcaSeed_dist_ = pvToPCAseed;
129  pcaTrack_dist_ = pvToPCAtrack;
130 
132  track_ip2d_ = t_ip2d.second.value();
133  track_ip2dSig_ = t_ip2d.second.significance();
134  track_ip3d_ = t_ip.second.value();
135  track_ip3dSig_ = t_ip.second.significance();
136 
137  dotprodTrackSeed2D_ = trackDir2D.unit().dot(seedDir2D.unit());
138  dotprodTrackSeed3D_ = iImpactState.globalDirection().unit().dot(tImpactState.globalDirection().unit());
139  dotprodTrackSeed2DV_ = trackPCADir2D.unit().dot(seedPCADir2D.unit());
140  dotprodTrackSeed3DV_ = (seedPosition - pvp).unit().dot((ttPoint - pvp).unit());
141 
142  pca_jetAxis_dist_ = jetLine.distance(cp).mag();
143  pca_jetAxis_dotprod_ = pairMomentum.unit().dot(jetdirection.unit());
144  pca_jetAxis_dEta_ = std::fabs(pvToPCA.eta() - jetdirection.eta());
145  pca_jetAxis_dPhi_ = std::fabs(pvToPCA.phi() - jetdirection.phi());
146  }
147  }
148 
149 } // namespace btagbtvdeep
Definition: Line.h:10
T z() const
Definition: PV3DBase.h:61
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T eta() const
Definition: PV3DBase.h:73
GlobalPoint crossingPoint() const override
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Definition: TTTypes.h:54
T mag() const
Definition: PV3DBase.h:64
std::pair< GlobalPoint, GlobalPoint > points() const override
Basic3DVector unit() const
GlobalVector distance(const Line &aLine) const
Definition: Line.h:35
void buildTrackPairInfo(const reco::TransientTrack *it, const reco::TransientTrack *tt, const reco::Vertex &pv, float mass, GlobalVector jetdirection, const std::pair< bool, Measurement1D > &t_ip, const std::pair< bool, Measurement1D > &t_ip2d)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Vector3DBase unit() const
Definition: Vector3DBase.h:54
float distance() const override
T dot(const Basic3DVector &rh) const
Scalar product, or "dot" product, with a vector of same type.
def cp(fromDir, toDir, listOfFiles, overwrite=False, smallList=False)