Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:27

0001 #include "RecoBTag/FeatureTools/interface/TrackPairInfoBuilder.h"
0002 #include "DataFormats/Candidate/interface/Candidate.h"
0003 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0004 #include "DataFormats/VertexReco/interface/Vertex.h"
0005 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0006 
0007 namespace btagbtvdeep {
0008 
0009   TrackPairInfoBuilder::TrackPairInfoBuilder()
0010       :
0011 
0012         track_pt_(0),
0013         track_eta_(0),
0014         track_phi_(0),
0015         track_dz_(0),
0016         track_dxy_(0),
0017 
0018         pca_distance_(0),
0019         pca_significance_(0),
0020 
0021         pcaSeed_x_(0),
0022         pcaSeed_y_(0),
0023         pcaSeed_z_(0),
0024         pcaSeed_xerr_(0),
0025         pcaSeed_yerr_(0),
0026         pcaSeed_zerr_(0),
0027         pcaTrack_x_(0),
0028         pcaTrack_y_(0),
0029         pcaTrack_z_(0),
0030         pcaTrack_xerr_(0),
0031         pcaTrack_yerr_(0),
0032         pcaTrack_zerr_(0),
0033 
0034         dotprodTrack_(0),
0035         dotprodSeed_(0),
0036         pcaSeed_dist_(0),
0037         pcaTrack_dist_(0),
0038 
0039         track_candMass_(0),
0040         track_ip2d_(0),
0041         track_ip2dSig_(0),
0042         track_ip3d_(0),
0043         track_ip3dSig_(0),
0044 
0045         dotprodTrackSeed2D_(0),
0046         dotprodTrackSeed2DV_(0),
0047         dotprodTrackSeed3D_(0),
0048         dotprodTrackSeed3DV_(0),
0049 
0050         pca_jetAxis_dist_(0),
0051         pca_jetAxis_dotprod_(0),
0052         pca_jetAxis_dEta_(0),
0053         pca_jetAxis_dPhi_(0)
0054 
0055   {}
0056 
0057   void TrackPairInfoBuilder::buildTrackPairInfo(const reco::TransientTrack* it,
0058                                                 const reco::TransientTrack* tt,
0059                                                 const reco::Vertex& pv,
0060                                                 float mass,
0061                                                 GlobalVector jetdirection,
0062                                                 const std::pair<bool, Measurement1D>& t_ip,
0063                                                 const std::pair<bool, Measurement1D>& t_ip2d) {
0064     GlobalPoint pvp(pv.x(), pv.y(), pv.z());
0065 
0066     VertexDistance3D distanceComputer;
0067     TwoTrackMinimumDistance dist;
0068 
0069     auto const& iImpactState = it->impactPointState();
0070     auto const& tImpactState = tt->impactPointState();
0071 
0072     if (dist.calculate(tImpactState, iImpactState)) {
0073       GlobalPoint ttPoint = dist.points().first;
0074       GlobalError ttPointErr = tImpactState.cartesianError().position();
0075       GlobalPoint seedPosition = dist.points().second;
0076       GlobalError seedPositionErr = iImpactState.cartesianError().position();
0077 
0078       Measurement1D m =
0079           distanceComputer.distance(VertexState(seedPosition, seedPositionErr), VertexState(ttPoint, ttPointErr));
0080 
0081       GlobalPoint cp(dist.crossingPoint());
0082 
0083       GlobalVector pairMomentum((Basic3DVector<float>)(it->track().momentum() + tt->track().momentum()));
0084       GlobalVector pvToPCA(cp - pvp);
0085 
0086       float pvToPCAseed = (seedPosition - pvp).mag();
0087       float pvToPCAtrack = (ttPoint - pvp).mag();
0088       float distance = dist.distance();
0089 
0090       GlobalVector trackDir2D(tImpactState.globalDirection().x(), tImpactState.globalDirection().y(), 0.);
0091       GlobalVector seedDir2D(iImpactState.globalDirection().x(), iImpactState.globalDirection().y(), 0.);
0092       GlobalVector trackPCADir2D(ttPoint.x() - pvp.x(), ttPoint.y() - pvp.y(), 0.);
0093       GlobalVector seedPCADir2D(seedPosition.x() - pvp.x(), seedPosition.y() - pvp.y(), 0.);
0094 
0095       float dotprodTrack = (ttPoint - pvp).unit().dot(tImpactState.globalDirection().unit());
0096       float dotprodSeed = (seedPosition - pvp).unit().dot(iImpactState.globalDirection().unit());
0097 
0098       Line::PositionType pos(pvp);
0099       Line::DirectionType dir(jetdirection);
0100       Line::DirectionType pairMomentumDir(pairMomentum);
0101       Line jetLine(pos, dir);
0102       Line PCAMomentumLine(cp, pairMomentumDir);
0103 
0104       track_pt_ = tt->track().pt();
0105       track_eta_ = tt->track().eta();
0106       track_phi_ = tt->track().phi();
0107       track_dz_ = tt->track().dz(pv.position());
0108       track_dxy_ = tt->track().dxy(pv.position());
0109 
0110       pca_distance_ = distance;
0111       pca_significance_ = m.significance();
0112 
0113       pcaSeed_x_ = seedPosition.x();
0114       pcaSeed_y_ = seedPosition.y();
0115       pcaSeed_z_ = seedPosition.z();
0116       pcaSeed_xerr_ = seedPositionErr.cxx();
0117       pcaSeed_yerr_ = seedPositionErr.cyy();
0118       pcaSeed_zerr_ = seedPositionErr.czz();
0119       pcaTrack_x_ = ttPoint.x();
0120       pcaTrack_y_ = ttPoint.y();
0121       pcaTrack_z_ = ttPoint.z();
0122       pcaTrack_xerr_ = ttPointErr.cxx();
0123       pcaTrack_yerr_ = ttPointErr.cyy();
0124       pcaTrack_zerr_ = ttPointErr.czz();
0125 
0126       dotprodTrack_ = dotprodTrack;
0127       dotprodSeed_ = dotprodSeed;
0128       pcaSeed_dist_ = pvToPCAseed;
0129       pcaTrack_dist_ = pvToPCAtrack;
0130 
0131       track_candMass_ = mass;
0132       track_ip2d_ = t_ip2d.second.value();
0133       track_ip2dSig_ = t_ip2d.second.significance();
0134       track_ip3d_ = t_ip.second.value();
0135       track_ip3dSig_ = t_ip.second.significance();
0136 
0137       dotprodTrackSeed2D_ = trackDir2D.unit().dot(seedDir2D.unit());
0138       dotprodTrackSeed3D_ = iImpactState.globalDirection().unit().dot(tImpactState.globalDirection().unit());
0139       dotprodTrackSeed2DV_ = trackPCADir2D.unit().dot(seedPCADir2D.unit());
0140       dotprodTrackSeed3DV_ = (seedPosition - pvp).unit().dot((ttPoint - pvp).unit());
0141 
0142       pca_jetAxis_dist_ = jetLine.distance(cp).mag();
0143       pca_jetAxis_dotprod_ = pairMomentum.unit().dot(jetdirection.unit());
0144       pca_jetAxis_dEta_ = std::fabs(pvToPCA.eta() - jetdirection.eta());
0145       pca_jetAxis_dPhi_ = std::fabs(pvToPCA.phi() - jetdirection.phi());
0146     }
0147   }
0148 
0149 }  // namespace btagbtvdeep