Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "RecoBTag/FeatureTools/interface/SeedingTrackInfoBuilder.h"
0002 #include "DataFormats/GeometrySurface/interface/Line.h"
0003 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0004 #include "TrackingTools/IPTools/interface/IPTools.h"
0005 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
0006 #include "FWCore/Utilities/interface/isFinite.h"
0007 #include "RecoBTag/TrackProbability/interface/HistogramProbabilityEstimator.h"
0008 
0009 namespace btagbtvdeep {
0010 
0011   SeedingTrackInfoBuilder::SeedingTrackInfoBuilder()
0012       : pt_(0),
0013         eta_(0),
0014         phi_(0),
0015         mass_(0),
0016         dz_(0),
0017         dxy_(0),
0018         ip3D_(0),
0019         sip3D_(0),
0020         ip2D_(0),
0021         sip2D_(0),
0022         ip3D_signed_(0),
0023         sip3D_signed_(0),
0024         ip2D_signed_(0),
0025         sip2D_signed_(0),
0026         chi2reduced_(0),
0027         nPixelHits_(0),
0028         nHits_(0),
0029         jetAxisDistance_(0),
0030         jetAxisDlength_(0),
0031         trackProbability3D_(0),
0032         trackProbability2D_(0) {}
0033 
0034   void SeedingTrackInfoBuilder::buildSeedingTrackInfo(const reco::TransientTrack* it,
0035                                                       const reco::Vertex& pv,
0036                                                       const reco::Jet& jet, /*GlobalVector jetdirection,*/
0037                                                       float mass,
0038                                                       const std::pair<bool, Measurement1D>& ip,
0039                                                       const std::pair<bool, Measurement1D>& ip2d,
0040                                                       float jet_distance,
0041                                                       float jaxis_dlength,
0042                                                       HistogramProbabilityEstimator* m_probabilityEstimator,
0043                                                       bool m_computeProbabilities = false) {
0044     GlobalPoint pvp(pv.x(), pv.y(), pv.z());
0045     GlobalVector jetdirection(jet.px(), jet.py(), jet.pz());
0046 
0047     auto const& aTrack = it->track();
0048 
0049     pt_ = aTrack.pt();
0050     eta_ = aTrack.eta();
0051     phi_ = aTrack.phi();
0052     dz_ = aTrack.dz(pv.position());
0053     dxy_ = aTrack.dxy(pv.position());
0054     mass_ = mass;
0055 
0056     std::pair<bool, Measurement1D> ipSigned = IPTools::signedImpactParameter3D(*it, jetdirection, pv);
0057     std::pair<bool, Measurement1D> ip2dSigned = IPTools::signedTransverseImpactParameter(*it, jetdirection, pv);
0058 
0059     ip3D_ = ip.second.value();
0060     sip3D_ = ip.second.significance();
0061     ip2D_ = ip2d.second.value();
0062     sip2D_ = ip2d.second.significance();
0063     ip3D_signed_ = ipSigned.second.value();
0064     sip3D_signed_ = ipSigned.second.significance();
0065     ip2D_signed_ = ip2dSigned.second.value();
0066     sip2D_signed_ = ip2dSigned.second.significance();
0067 
0068     chi2reduced_ = aTrack.normalizedChi2();
0069     nPixelHits_ = aTrack.hitPattern().numberOfValidPixelHits();
0070     nHits_ = aTrack.hitPattern().numberOfValidHits();
0071 
0072     jetAxisDistance_ = std::fabs(jet_distance);
0073     jetAxisDlength_ = jaxis_dlength;
0074 
0075     trackProbability3D_ = 0.5;
0076     trackProbability2D_ = 0.5;
0077 
0078     if (m_computeProbabilities) {
0079       //probability with 3D ip
0080       std::pair<bool, double> probability =
0081           m_probabilityEstimator->probability(false, 0, ip.second.significance(), aTrack, jet, pv);
0082       double prob3D = (probability.first ? probability.second : -1.);
0083 
0084       //probability with 2D ip
0085       probability = m_probabilityEstimator->probability(false, 1, ip2d.second.significance(), aTrack, jet, pv);
0086       double prob2D = (probability.first ? probability.second : -1.);
0087 
0088       trackProbability3D_ = prob3D;
0089       trackProbability2D_ = prob2D;
0090     }
0091 
0092     if (!edm::isFinite(trackProbability3D_))
0093       trackProbability3D_ = 0.5;
0094     if (!edm::isFinite(trackProbability2D_))
0095       trackProbability2D_ = 0.5;
0096   }
0097 
0098 }  // namespace btagbtvdeep