Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:32

0001 #ifndef TrackingTools_IPTools_h
0002 #define TrackingTools_IPTools_h
0003 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0004 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
0005 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0006 #include "DataFormats/VertexReco/interface/Vertex.h"
0007 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0008 #include <utility>
0009 #include "DataFormats/CLHEP/interface/Migration.h"
0010 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
0011 #include "RecoVertex/VertexTools/interface/VertexDistance.h"
0012 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0013 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0014 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
0015 
0016 namespace IPTools {
0017   /**
0018    *  Returns the unsigned transverse impact parameter
0019    *  The track is extrapolated to the closest point to the primary vertex in transverse plane
0020    *  then the impact parameter and its error are computed
0021    */
0022   std::pair<bool, Measurement1D> absoluteImpactParameter3D(const reco::TransientTrack& transientTrack,
0023                                                            const reco::Vertex& vertex);
0024 
0025   /**
0026    *  Returns the unsigned 3D impact parameter
0027    *  The track is extrapolated to the closest point to the primary vertex in 3d space
0028    *  then the impact parameter and its error are computed
0029    */
0030 
0031   std::pair<bool, Measurement1D> absoluteTransverseImpactParameter(const reco::TransientTrack& transientTrack,
0032                                                                    const reco::Vertex& vertex);
0033 
0034   /**
0035    *  Returns life time signed transverse impact parameter
0036    *  The track is extrapolated to the closest point to the primary vertex in transverse plane
0037    *  then the impact parameter and its error are computed
0038    */
0039   std::pair<bool, Measurement1D> signedTransverseImpactParameter(const reco::TransientTrack& track,
0040                                                                  const GlobalVector& direction,
0041                                                                  const reco::Vertex& vertex);
0042 
0043   /**
0044    *  Returns life time signed 3D impact parameter
0045    *  The track is extrapolated to the closest point to the primary vertex in 3d space
0046    *  then the impact parameter and its error are computed
0047    */
0048 
0049   std::pair<bool, Measurement1D> signedImpactParameter3D(const reco::TransientTrack& track,
0050                                                          const GlobalVector& direction,
0051                                                          const reco::Vertex& vertex);
0052 
0053   /// Impact parameter without direction (internally used)
0054   std::pair<bool, Measurement1D> absoluteImpactParameter(const TrajectoryStateOnSurface& tsos,
0055                                                          const reco::Vertex& vertex,
0056                                                          VertexDistance& distanceComputer);
0057 
0058   inline TrajectoryStateOnSurface transverseExtrapolate(const TrajectoryStateOnSurface& track,
0059                                                         const GlobalPoint& vertexPosition,
0060                                                         const MagneticField* field) {
0061     TransverseImpactPointExtrapolator extrapolator(field);
0062     return extrapolator.extrapolate(track, vertexPosition);
0063   }
0064 
0065   TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface& state,
0066                                                 const reco::Vertex& vertex,
0067                                                 const GlobalVector& aJetDirection,
0068                                                 const MagneticField* field);
0069 
0070   GlobalVector linearImpactParameter(const TrajectoryStateOnSurface& aTSOS, const GlobalPoint& point);
0071 
0072   std::pair<bool, Measurement1D> linearizedSignedImpactParameter3D(const TrajectoryStateOnSurface& state,
0073                                                                    const GlobalVector& direction,
0074                                                                    const reco::Vertex& vertex);
0075 
0076   inline std::pair<bool, Measurement1D> linearizedSignedImpactParameter3D(const reco::TransientTrack& transientTrack,
0077                                                                           const GlobalVector& direction,
0078                                                                           const reco::Vertex& vertex) {
0079     // extrapolate to the point of closest approach to the jet axis
0080     TrajectoryStateOnSurface closestToJetState =
0081         closestApproachToJet(transientTrack.impactPointState(), vertex, direction, transientTrack.field());
0082     return linearizedSignedImpactParameter3D(closestToJetState, direction, vertex);
0083   }
0084 
0085   std::pair<bool, Measurement1D> signedDecayLength3D(const TrajectoryStateOnSurface& state,
0086                                                      const GlobalVector& direction,
0087                                                      const reco::Vertex& vertex);
0088 
0089   inline std::pair<bool, Measurement1D> signedDecayLength3D(const reco::TransientTrack& transientTrack,
0090                                                             const GlobalVector& direction,
0091                                                             const reco::Vertex& vertex) {
0092     // extrapolate to the point of closest approach to the jet axis
0093     TrajectoryStateOnSurface closestToJetState =
0094         closestApproachToJet(transientTrack.impactPointState(), vertex, direction, transientTrack.field());
0095     return signedDecayLength3D(closestToJetState, direction, vertex);
0096   }
0097 
0098   std::pair<double, Measurement1D> jetTrackDistance(const reco::TransientTrack& track,
0099                                                     const GlobalVector& direction,
0100                                                     const reco::Vertex& vertex);
0101 }  // namespace IPTools
0102 
0103 #endif