Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-13 03:23:59

0001 #include "DataFormats/PatCandidates/interface/PackedCandidate.h"
0002 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0003 
0004 #include "DataFormats/BTauReco/interface/SeedingTrackFeatures.h"
0005 #include "DataFormats/BTauReco/interface/TrackPairFeatures.h"
0006 
0007 #include "RecoBTag/FeatureTools/interface/TrackPairInfoBuilder.h"
0008 #include "RecoBTag/FeatureTools/interface/SeedingTrackInfoBuilder.h"
0009 
0010 #include "DataFormats/PatCandidates/interface/Jet.h"
0011 
0012 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0013 #include "TrackingTools/IPTools/interface/IPTools.h"
0014 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0015 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0016 
0017 #include "RecoBTag/TrackProbability/interface/HistogramProbabilityEstimator.h"
0018 
0019 #include "RecoBTag/FeatureTools/interface/SeedingTracksConverter.h"
0020 
0021 namespace btagbtvdeep {
0022 
0023   void seedingTracksToFeatures(const std::vector<reco::TransientTrack>& selectedTracks,
0024                                const std::vector<float>& masses,
0025                                const reco::Jet& jet,
0026                                const reco::Vertex& pv,
0027                                HistogramProbabilityEstimator* probabilityEstimator,
0028                                bool computeProbabilities,
0029                                std::vector<btagbtvdeep::SeedingTrackFeatures>& seedingT_features_vector)
0030 
0031   {
0032     GlobalVector jetdirection(jet.px(), jet.py(), jet.pz());
0033     GlobalPoint pvp(pv.x(), pv.y(), pv.z());
0034 
0035     std::multimap<double, std::pair<btagbtvdeep::SeedingTrackInfoBuilder, std::vector<btagbtvdeep::TrackPairFeatures>>>
0036         sortedSeedsMap;
0037     std::multimap<double, btagbtvdeep::TrackPairInfoBuilder> sortedNeighboursMap;
0038 
0039     std::vector<btagbtvdeep::TrackPairFeatures> tp_features_vector;
0040 
0041     sortedSeedsMap.clear();
0042     seedingT_features_vector.clear();
0043 
0044     std::vector<std::pair<bool, Measurement1D>> absIP3D(selectedTracks.size());
0045     std::vector<std::pair<bool, Measurement1D>> absIP2D(selectedTracks.size());
0046     std::vector<bool> absIP3D_filled(selectedTracks.size(), false);
0047     std::vector<bool> absIP2D_filled(selectedTracks.size(), false);
0048 
0049     unsigned int selTrackCount = 0;
0050 
0051     for (auto const& it : selectedTracks) {
0052       selTrackCount += 1;
0053       sortedNeighboursMap.clear();
0054       tp_features_vector.clear();
0055 
0056       if (reco::deltaR(it.track(), jet) > 0.4)
0057         continue;
0058 
0059       std::pair<bool, Measurement1D> ip = IPTools::absoluteImpactParameter3D(it, pv);
0060 
0061       absIP3D[selTrackCount - 1] = ip;
0062       absIP3D_filled[selTrackCount - 1] = true;
0063 
0064       std::pair<double, Measurement1D> jet_dist = IPTools::jetTrackDistance(it, jetdirection, pv);
0065       TrajectoryStateOnSurface closest =
0066           IPTools::closestApproachToJet(it.impactPointState(), pv, jetdirection, it.field());
0067       float length = 999;
0068       if (closest.isValid())
0069         length = (closest.globalPosition() - pvp).mag();
0070 
0071       if (!(ip.first && ip.second.value() >= 0.0 && ip.second.significance() >= 1.0 && ip.second.value() <= 9999. &&
0072             ip.second.significance() <= 9999. && it.track().normalizedChi2() < 5. &&
0073             std::fabs(it.track().dxy(pv.position())) < 2 && std::fabs(it.track().dz(pv.position())) < 17 &&
0074             jet_dist.second.value() < 0.07 && length < 5.))
0075         continue;
0076 
0077       std::pair<bool, Measurement1D> ip2d = IPTools::absoluteTransverseImpactParameter(it, pv);
0078 
0079       absIP2D[selTrackCount - 1] = ip2d;
0080       absIP2D_filled[selTrackCount - 1] = true;
0081 
0082       btagbtvdeep::SeedingTrackInfoBuilder seedInfo;
0083       seedInfo.buildSeedingTrackInfo(&(it),
0084                                      pv,
0085                                      jet,
0086                                      masses[selTrackCount - 1],
0087                                      ip,
0088                                      ip2d,
0089                                      jet_dist.second.value(),
0090                                      length,
0091                                      probabilityEstimator,
0092                                      computeProbabilities);
0093 
0094       unsigned int neighbourTrackCount = 0;
0095 
0096       for (auto const& tt : selectedTracks) {
0097         neighbourTrackCount += 1;
0098 
0099         if (neighbourTrackCount == selTrackCount)
0100           continue;
0101         if (std::fabs(pv.z() - tt.track().vz()) > 0.1)
0102           continue;
0103 
0104         //avoid calling IPs twice
0105         if (!absIP2D_filled[neighbourTrackCount - 1]) {
0106           absIP2D[neighbourTrackCount - 1] = IPTools::absoluteTransverseImpactParameter(tt, pv);
0107           absIP2D_filled[neighbourTrackCount - 1] = true;
0108         }
0109 
0110         if (!absIP3D_filled[neighbourTrackCount - 1]) {
0111           absIP3D[neighbourTrackCount - 1] = IPTools::absoluteImpactParameter3D(tt, pv);
0112           absIP3D_filled[neighbourTrackCount - 1] = true;
0113         }
0114 
0115         std::pair<bool, Measurement1D> t_ip = absIP3D[neighbourTrackCount - 1];
0116         std::pair<bool, Measurement1D> t_ip2d = absIP2D[neighbourTrackCount - 1];
0117 
0118         btagbtvdeep::TrackPairInfoBuilder trackPairInfo;
0119         trackPairInfo.buildTrackPairInfo(&(it), &(tt), pv, masses[neighbourTrackCount - 1], jetdirection, t_ip, t_ip2d);
0120         sortedNeighboursMap.insert(std::make_pair(trackPairInfo.pca_distance(), trackPairInfo));
0121       }
0122 
0123       int max_counter = 0;
0124 
0125       for (auto const& im : sortedNeighboursMap) {
0126         if (max_counter >= 20)
0127           break;
0128         btagbtvdeep::TrackPairFeatures tp_features;
0129 
0130         auto const& tp = im.second;
0131 
0132         tp_features.pt = (tp.track_pt() == 0) ? 0 : 1.0 / tp.track_pt();
0133         tp_features.eta = tp.track_eta();
0134         tp_features.phi = tp.track_phi();
0135         tp_features.mass = tp.track_candMass();
0136         tp_features.dz = logWithOffset(tp.track_dz());
0137         tp_features.dxy = logWithOffset(tp.track_dxy());
0138         tp_features.ip3D = log(tp.track_ip3d());
0139         tp_features.sip3D = log(tp.track_ip3dSig());
0140         tp_features.ip2D = log(tp.track_ip2d());
0141         tp_features.sip2D = log(tp.track_ip2dSig());
0142         tp_features.distPCA = log(tp.pca_distance());
0143         tp_features.dsigPCA = log(tp.pca_significance());
0144         tp_features.x_PCAonSeed = tp.pcaSeed_x();
0145         tp_features.y_PCAonSeed = tp.pcaSeed_y();
0146         tp_features.z_PCAonSeed = tp.pcaSeed_z();
0147         tp_features.xerr_PCAonSeed = tp.pcaSeed_xerr();
0148         tp_features.yerr_PCAonSeed = tp.pcaSeed_yerr();
0149         tp_features.zerr_PCAonSeed = tp.pcaSeed_zerr();
0150         tp_features.x_PCAonTrack = tp.pcaTrack_x();
0151         tp_features.y_PCAonTrack = tp.pcaTrack_y();
0152         tp_features.z_PCAonTrack = tp.pcaTrack_z();
0153         tp_features.xerr_PCAonTrack = tp.pcaTrack_xerr();
0154         tp_features.yerr_PCAonTrack = tp.pcaTrack_yerr();
0155         tp_features.zerr_PCAonTrack = tp.pcaTrack_zerr();
0156         tp_features.dotprodTrack = tp.dotprodTrack();
0157         tp_features.dotprodSeed = tp.dotprodSeed();
0158         tp_features.dotprodTrackSeed2D = tp.dotprodTrackSeed2D();
0159         tp_features.dotprodTrackSeed3D = tp.dotprodTrackSeed3D();
0160         tp_features.dotprodTrackSeedVectors2D = tp.dotprodTrackSeed2DV();
0161         tp_features.dotprodTrackSeedVectors3D = tp.dotprodTrackSeed3DV();
0162         tp_features.pvd_PCAonSeed = log(tp.pcaSeed_dist());
0163         tp_features.pvd_PCAonTrack = log(tp.pcaTrack_dist());
0164         tp_features.dist_PCAjetAxis = log(tp.pca_jetAxis_dist());
0165         tp_features.dotprod_PCAjetMomenta = tp.pca_jetAxis_dotprod();
0166         tp_features.deta_PCAjetDirs = log(tp.pca_jetAxis_dEta());
0167         tp_features.dphi_PCAjetDirs = tp.pca_jetAxis_dPhi();
0168 
0169         max_counter = max_counter + 1;
0170         tp_features_vector.push_back(tp_features);
0171       }
0172 
0173       sortedSeedsMap.insert(std::make_pair(-seedInfo.sip3d_Signed(), std::make_pair(seedInfo, tp_features_vector)));
0174     }
0175 
0176     int max_counter_seed = 0;
0177 
0178     for (auto const& im : sortedSeedsMap) {
0179       if (max_counter_seed >= 10)
0180         break;
0181 
0182       btagbtvdeep::SeedingTrackFeatures seed_features;
0183 
0184       auto const& seed = im.second.first;
0185 
0186       seed_features.nearTracks = im.second.second;
0187       seed_features.pt = (seed.pt() == 0) ? 0 : 1.0 / seed.pt();
0188       seed_features.eta = seed.eta();
0189       seed_features.phi = seed.phi();
0190       seed_features.mass = seed.mass();
0191       seed_features.dz = logWithOffset(seed.dz());
0192       seed_features.dxy = logWithOffset(seed.dxy());
0193       seed_features.ip3D = log(seed.ip3d());
0194       seed_features.sip3D = log(seed.sip3d());
0195       seed_features.ip2D = log(seed.ip2d());
0196       seed_features.sip2D = log(seed.sip2d());
0197       seed_features.signedIp3D = logWithOffset(seed.ip3d_Signed());
0198       seed_features.signedSip3D = logWithOffset(seed.sip3d_Signed());
0199       seed_features.signedIp2D = logWithOffset(seed.ip2d_Signed());
0200       seed_features.signedSip2D = logWithOffset(seed.sip2d_Signed());
0201       seed_features.trackProbability3D = seed.trackProbability3D();
0202       seed_features.trackProbability2D = seed.trackProbability2D();
0203       seed_features.chi2reduced = seed.chi2reduced();
0204       seed_features.nPixelHits = seed.nPixelHits();
0205       seed_features.nHits = seed.nHits();
0206       seed_features.jetAxisDistance = log(seed.jetAxisDistance());
0207       seed_features.jetAxisDlength = log(seed.jetAxisDlength());
0208 
0209       max_counter_seed = max_counter_seed + 1;
0210       seedingT_features_vector.push_back(seed_features);
0211     }
0212 
0213     if (sortedSeedsMap.size() < 10) {
0214       for (unsigned int i = sortedSeedsMap.size(); i < 10; i++) {
0215         std::vector<btagbtvdeep::TrackPairFeatures> tp_features_zeropad(20);
0216         btagbtvdeep::SeedingTrackFeatures seed_features_zeropad{};
0217         seed_features_zeropad.nearTracks = tp_features_zeropad;
0218         seedingT_features_vector.push_back(seed_features_zeropad);
0219       }
0220     }
0221   }
0222 
0223 }  // namespace btagbtvdeep