CMS 3D CMS Logo

TrackPairInfoBuilder.cc
Go to the documentation of this file.
6 
7 namespace btagbtvdeep{
8 
10 
11  track_pt_(0),
12  track_eta_(0),
13  track_phi_(0),
14  track_dz_(0),
15  track_dxy_(0),
16 
17  pca_distance_(0),
18  pca_significance_(0),
19 
20  pcaSeed_x_(0),
21  pcaSeed_y_(0),
22  pcaSeed_z_(0),
23  pcaSeed_xerr_(0),
24  pcaSeed_yerr_(0),
25  pcaSeed_zerr_(0),
26  pcaTrack_x_(0),
27  pcaTrack_y_(0),
28  pcaTrack_z_(0),
29  pcaTrack_xerr_(0),
30  pcaTrack_yerr_(0),
31  pcaTrack_zerr_(0),
32 
33  dotprodTrack_(0),
34  dotprodSeed_(0),
35  pcaSeed_dist_(0),
36  pcaTrack_dist_(0),
37 
38  track_candMass_(0),
39  track_ip2d_(0),
40  track_ip2dSig_(0),
41  track_ip3d_(0),
42  track_ip3dSig_(0),
43 
44  dotprodTrackSeed2D_(0),
45  dotprodTrackSeed2DV_(0),
46  dotprodTrackSeed3D_(0),
47  dotprodTrackSeed3DV_(0),
48 
49  pca_jetAxis_dist_(0),
50  pca_jetAxis_dotprod_(0),
51  pca_jetAxis_dEta_(0),
52  pca_jetAxis_dPhi_(0)
53 
54 {
55 
56 }
57 
58 
60  const std::pair<bool,Measurement1D> & t_ip, const std::pair<bool,Measurement1D> & t_ip2d )
61  {
62 
63  GlobalPoint pvp(pv.x(),pv.y(),pv.z());
64 
65  VertexDistance3D distanceComputer;
67 
68  auto const& iImpactState = it->impactPointState();
69  auto const& tImpactState = tt->impactPointState();
70 
71  if(dist.calculate(tImpactState,iImpactState)) {
72 
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 
78  Measurement1D m = distanceComputer.distance(VertexState(seedPosition,seedPositionErr), VertexState(ttPoint, ttPointErr));
79 
80  GlobalPoint cp(dist.crossingPoint());
81 
82  GlobalVector pairMomentum((Basic3DVector<float>) (it->track().momentum()+tt->track().momentum()));
83  GlobalVector pvToPCA(cp-pvp);
84 
85  float pvToPCAseed = (seedPosition-pvp).mag();
86  float pvToPCAtrack = (ttPoint-pvp).mag();
87  float distance = dist.distance();
88 
89  GlobalVector trackDir2D(tImpactState.globalDirection().x(),tImpactState.globalDirection().y(),0.);
90  GlobalVector seedDir2D(iImpactState.globalDirection().x(),iImpactState.globalDirection().y(),0.);
91  GlobalVector trackPCADir2D(ttPoint.x()-pvp.x(),ttPoint.y()-pvp.y(),0.);
92  GlobalVector seedPCADir2D(seedPosition.x()-pvp.x(),seedPosition.y()-pvp.y(),0.);
93 
94  float dotprodTrack = (ttPoint-pvp).unit().dot(tImpactState.globalDirection().unit());
95  float dotprodSeed = (seedPosition-pvp).unit().dot(iImpactState.globalDirection().unit());
96 
98  Line::DirectionType dir(jetdirection);
99  Line::DirectionType pairMomentumDir(pairMomentum);
100  Line jetLine(pos,dir);
101  Line PCAMomentumLine(cp,pairMomentumDir);
102 
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 
110 
113 
114  pcaSeed_x_=seedPosition.x();
115  pcaSeed_y_=seedPosition.y();
116  pcaSeed_z_=seedPosition.z();
117  pcaSeed_xerr_=seedPositionErr.cxx();
118  pcaSeed_yerr_=seedPositionErr.cyy();
119  pcaSeed_zerr_=seedPositionErr.czz();
120  pcaTrack_x_=ttPoint.x();
121  pcaTrack_y_=ttPoint.y();
122  pcaTrack_z_=ttPoint.z();
123  pcaTrack_xerr_=ttPointErr.cxx();
124  pcaTrack_yerr_=ttPointErr.cyy();
125  pcaTrack_zerr_=ttPointErr.czz();
126 
129  pcaSeed_dist_=pvToPCAseed;
130  pcaTrack_dist_=pvToPCAtrack;
131 
133  track_ip2d_=t_ip2d.second.value();
134  track_ip2dSig_=t_ip2d.second.significance();
135  track_ip3d_=t_ip.second.value();
136  track_ip3dSig_=t_ip.second.significance();
137 
138  dotprodTrackSeed2D_=trackDir2D.unit().dot(seedDir2D.unit());
139  dotprodTrackSeed3D_=iImpactState.globalDirection().unit().dot(tImpactState.globalDirection().unit());
140  dotprodTrackSeed2DV_=trackPCADir2D.unit().dot(seedPCADir2D.unit());
141  dotprodTrackSeed3DV_=(seedPosition-pvp).unit().dot((ttPoint-pvp).unit());
142 
143  pca_jetAxis_dist_=jetLine.distance(cp).mag();
144  pca_jetAxis_dotprod_=pairMomentum.unit().dot(jetdirection.unit());
145  pca_jetAxis_dEta_=std::fabs(pvToPCA.eta()-jetdirection.eta());
146  pca_jetAxis_dPhi_=std::fabs(pvToPCA.phi()-jetdirection.phi());
147 
148 
149  }
150  }
151 
152 
153 }
std::pair< GlobalPoint, GlobalPoint > points() const override
Definition: Line.h:10
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
double y() const
y coordinate
Definition: Vertex.h:113
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:684
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:714
const Point & position() const
position
Definition: Vertex.h:109
GlobalVector distance(const Line &aLine) const
Definition: Line.h:38
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
T mag() const
Definition: PV3DBase.h:67
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:690
double pt() const
track transverse momentum
Definition: TrackBase.h:660
T z() const
Definition: PV3DBase.h:64
def pv(vc)
Definition: MetAnalyzer.py:7
double z() const
z coordinate
Definition: Vertex.h:115
GlobalPoint crossingPoint() const override
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:648
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)
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double x() const
x coordinate
Definition: Vertex.h:111
double significance() const
Definition: Measurement1D.h:29
const Track & track() const
float distance() const override
T eta() const
Definition: PV3DBase.h:76
TrajectoryStateOnSurface impactPointState() const
dbl *** dir
Definition: mlp_gen.cc:35
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:630
T x() const
Definition: PV3DBase.h:62