CMS 3D CMS Logo

Functions

IPTools Namespace Reference

Functions

std::pair< bool, Measurement1DabsoluteImpactParameter (const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
 Impact parameter without direction (internally used)
std::pair< bool, Measurement1DabsoluteImpactParameter3D (const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
std::pair< bool, Measurement1DabsoluteTransverseImpactParameter (const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
TrajectoryStateOnSurface closestApproachToJet (const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
std::pair< double, Measurement1DjetTrackDistance (const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
GlobalVector linearImpactParameter (const TrajectoryStateOnSurface &aTSOS, const GlobalPoint &point)
std::pair< bool, Measurement1DlinearizedSignedImpactParameter3D (const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
std::pair< bool, Measurement1DlinearizedSignedImpactParameter3D (const reco::TransientTrack &transientTrack, const GlobalVector &direction, const reco::Vertex &vertex)
std::pair< bool, Measurement1DsignedDecayLength3D (const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
std::pair< bool, Measurement1DsignedDecayLength3D (const reco::TransientTrack &transientTrack, const GlobalVector &direction, const reco::Vertex &vertex)
std::pair< bool, Measurement1DsignedImpactParameter3D (const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
std::pair< bool, Measurement1DsignedTransverseImpactParameter (const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
TrajectoryStateOnSurface transverseExtrapolate (const TrajectoryStateOnSurface &track, const GlobalPoint &vertexPosition, const MagneticField *field)

Function Documentation

std::pair< bool, Measurement1D > IPTools::absoluteImpactParameter ( const TrajectoryStateOnSurface tsos,
const reco::Vertex vertex,
VertexDistance distanceComputer 
)
std::pair< bool, Measurement1D > IPTools::absoluteImpactParameter3D ( const reco::TransientTrack transientTrack,
const reco::Vertex vertex 
)

Returns the unsigned transverse impact parameter The track is extrapolated to the closest point to the primary vertex in transverse plane then the impact parameter and its error are computed

Definition at line 37 of file IPTools.cc.

References absoluteImpactParameter(), RecoVertex::convertPos(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), and reco::Vertex::position().

Referenced by InclusiveVertexFinder::nearTracks(), TrackVertexArbitrator::produce(), and InclusiveVertexFinder::produce().

    {
      AnalyticalImpactPointExtrapolator extrapolator(transientTrack.field());
      VertexDistance3D dist;
      return absoluteImpactParameter(extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())), vertex, dist);
    }
std::pair< bool, Measurement1D > IPTools::absoluteTransverseImpactParameter ( const reco::TransientTrack transientTrack,
const reco::Vertex vertex 
)

Returns the unsigned 3D impact parameter The track is extrapolated to the closest point to the primary vertex in 3d space then the impact parameter and its error are computed

Definition at line 43 of file IPTools.cc.

References absoluteImpactParameter(), RecoVertex::convertPos(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), and reco::Vertex::position().

Referenced by pat::PATElectronProducer::produce(), and pat::PATMuonProducer::produce().

    {
      TransverseImpactPointExtrapolator extrapolator(transientTrack.field());
      VertexDistanceXY dist;
      return absoluteImpactParameter(extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())), vertex,  dist);
    }
TrajectoryStateOnSurface IPTools::closestApproachToJet ( const TrajectoryStateOnSurface state,
const reco::Vertex vertex,
const GlobalVector aJetDirection,
const MagneticField field 
)
pair< double, Measurement1D > IPTools::jetTrackDistance ( const reco::TransientTrack track,
const GlobalVector direction,
const reco::Vertex vertex 
)

Definition at line 200 of file IPTools.cc.

References Line::distance(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), reco::TransientTrack::impactPointState(), TrajectoryStateOnSurface::isValid(), mag(), P, Line::position(), csvLumiCalc::unit, CommonMethods::weight(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by QualityCutsAnalyzer::LoopOverJetTracksAssociation(), and TrackIPProducer::produce().

                                                                                                                                   {
    double  theLDist_err(0.);
  
    //FIXME
    float weight=0.;//vertex.trackWeight(track);

    TrajectoryStateOnSurface stateAtOrigin = track.impactPointState(); 
    if(!stateAtOrigin.isValid())
      {
        //TODO: throw instead?
        return pair<bool,Measurement1D>(false,Measurement1D(0.,0.));
      }
   
    //get the Track line at origin
    Line::PositionType posTrack(stateAtOrigin.globalPosition());
    Line::DirectionType dirTrack((stateAtOrigin.globalMomentum()).unit());
    Line trackLine(posTrack,dirTrack);
    // get the Jet  line 
    // Vertex vertex(vertex);
    GlobalVector jetVector = direction.unit();    
    Line::PositionType posJet(GlobalPoint(vertex.x(),vertex.y(),vertex.z()));
    Line::DirectionType dirJet(jetVector);
    Line jetLine(posJet,dirJet);
  
    // now compute the distance between the two lines
    // If the track has been used to refit the Primary vertex then sign it positively, otherwise negative
    double theDistanceToJetAxis = (jetLine.distance(trackLine)).mag();
    if (weight<1) theDistanceToJetAxis= -theDistanceToJetAxis;

    // ... and the flight distance along the Jet axis.
    GlobalPoint  V = jetLine.position();    
    GlobalVector Q = dirTrack - jetVector.dot(dirTrack) * jetVector;
    GlobalVector P = jetVector - jetVector.dot(dirTrack) * dirTrack;
    double theDistanceAlongJetAxis = P.dot(V-posTrack)/Q.dot(dirTrack);

    //
    // get the covariance matrix of the vertex and compute the error on theDistanceToJetAxis
    //
    

    // build the vector of closest approach between lines


    //FIXME: error not computed.
    GlobalVector H((jetVector.cross(dirTrack).unit()));
    CLHEP::HepVector Hh(3);
    Hh[0] = H.x();
    Hh[1] = H.y();
    Hh[2] = H.z();
    
    //  theLDist_err = sqrt(vertexError.similarity(Hh));

    //    cout << "distance to jet axis : "<< theDistanceToJetAxis <<" and error : "<< theLDist_err<<endl;
    // Now the impact parameter ...

    /*    GlobalPoint T0 = track.position();
          GlobalVector D = (T0-V)- (T0-V).dot(dir) * dir;
          double IP = D.mag();    
          GlobalVector Dold = distance(aTSOS, aJet.vertex(), jetDirection);
          double IPold = Dold.mag();
    */



  
    Measurement1D DTJA(theDistanceToJetAxis,theLDist_err);
  
    return pair<double,Measurement1D> (theDistanceAlongJetAxis,DTJA);
  }
GlobalVector IPTools::linearImpactParameter ( const TrajectoryStateOnSurface state,
const GlobalPoint point 
)

Compute the impact parameter of a track, linearized from the given state, with respect to a given point

Definition at line 191 of file IPTools.cc.

References dir, TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), point, pos, tmp, and csvLumiCalc::unit.

Referenced by linearizedSignedImpactParameter3D().

                                                                                                         {

    Line::PositionType pos(state.globalPosition());
    Line::DirectionType dir((state.globalMomentum()).unit());
    Line trackLine(pos,dir);
    GlobalPoint  tmp=point; 
    return  trackLine.distance(tmp);
  }
pair< bool, Measurement1D > IPTools::linearizedSignedImpactParameter3D ( const TrajectoryStateOnSurface state,
const GlobalVector direction,
const reco::Vertex vertex 
)

Definition at line 133 of file IPTools.cc.

References abs, TrajectoryStateOnSurface::cartesianError(), reco::Vertex::covariance(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), linearImpactParameter(), CartesianTrajectoryError::matrix(), mathSSE::sqrt(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by linearizedSignedImpactParameter3D().

  {
    //Check if extrapolation has been successfull
    if(!closestToJetState.isValid()) {
      return pair<bool,Measurement1D>(false,Measurement1D(0.,0.));
    }

    GlobalPoint vertexPosition(vertex.x(),vertex.y(),vertex.z());
    GlobalVector impactParameter = linearImpactParameter(closestToJetState, vertexPosition);
    GlobalVector jetDir = direction.unit();
    GlobalVector flightDistance(closestToJetState.globalPosition()-vertexPosition);
    double theDistanceAlongJetAxis = jetDir.dot(flightDistance);
    double signedIP = impactParameter.mag()*((theDistanceAlongJetAxis!=0)?theDistanceAlongJetAxis/abs(theDistanceAlongJetAxis):1.);


    GlobalVector ipDirection = impactParameter.unit();
    GlobalPoint closestPoint = closestToJetState.globalPosition();
    GlobalVector momentumAtClosestPoint = closestToJetState.globalMomentum(); 
    GlobalVector momentumDir = momentumAtClosestPoint.unit();
    
    AlgebraicVector3 deriv_v;
    deriv_v[0] = - ipDirection.x();
    deriv_v[1] = - ipDirection.y();
    deriv_v[2] = - ipDirection.z();

 
    AlgebraicVector6 deriv;
    deriv[0] = ipDirection.x();
    deriv[1] = ipDirection.y();
    deriv[2] = ipDirection.z();
    deriv[3] = -  (momentumDir.dot(flightDistance)*ipDirection.x())/momentumAtClosestPoint.mag();
    deriv[4] = -  (momentumDir.dot(flightDistance)*ipDirection.y())/momentumAtClosestPoint.mag();
    deriv[5] = -  (momentumDir.dot(flightDistance)*ipDirection.z())/momentumAtClosestPoint.mag();

    double trackError2 = ROOT::Math::Similarity(deriv , closestToJetState.cartesianError().matrix());
    double vertexError2 = ROOT::Math::Similarity(deriv_v , vertex.covariance());
    double ipError = sqrt(trackError2+vertexError2);

    return pair<bool,Measurement1D>(true,Measurement1D(signedIP,ipError));
  }
std::pair<bool,Measurement1D> IPTools::linearizedSignedImpactParameter3D ( const reco::TransientTrack transientTrack,
const GlobalVector direction,
const reco::Vertex vertex 
) [inline]

Definition at line 73 of file IPTools.h.

References closestApproachToJet(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), and linearizedSignedImpactParameter3D().

        {
          // extrapolate to the point of closest approach to the jet axis
          TrajectoryStateOnSurface closestToJetState = closestApproachToJet(transientTrack.impactPointState(), vertex, direction,transientTrack.field());
         return linearizedSignedImpactParameter3D(closestToJetState,direction,vertex);
        }
pair< bool, Measurement1D > IPTools::signedDecayLength3D ( const TrajectoryStateOnSurface state,
const GlobalVector direction,
const reco::Vertex vertex 
)

chech it!!!!!!!!!!!!!!!!!!!!!!!

chech it!!!!!!!!!!!!!!!!!!!!!!!

Definition at line 94 of file IPTools.cc.

References TrajectoryStateOnSurface::cartesianError(), reco::Vertex::covariance(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), j, CartesianTrajectoryError::matrix(), mathSSE::sqrt(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by QualityCutsAnalyzer::LoopOverJetTracksAssociation(), and signedDecayLength3D().

                                                                                                        {

    //Check if extrapolation has been successfull
    if(!closestToJetState.isValid()) {
      return pair<bool,Measurement1D>(false,Measurement1D(0.,0.));
    }

    GlobalVector jetDirection = direction.unit();
    GlobalPoint vertexPosition(vertex.x(),vertex.y(),vertex.z());
  
    double decayLen = jetDirection.dot(closestToJetState.globalPosition()-vertexPosition);

    //error calculation

    AlgebraicVector3 j;
    j[0] = jetDirection.x();
    j[1] = jetDirection.y();
    j[2] = jetDirection.z();
    AlgebraicVector6 jj;
    jj[0] = jetDirection.x();
    jj[1] = jetDirection.y();
    jj[2] = jetDirection.z();
    jj[3] =0.;
    jj[4] =0.;
    jj[5] =0.;
     
    //TODO: FIXME: the extrapolation uncertainty is very relevant here should be taken into account!!
    double trackError2 = ROOT::Math::Similarity(jj,closestToJetState.cartesianError().matrix());
    double vertexError2 = ROOT::Math::Similarity(j,vertex.covariance());

    double decayLenError = sqrt(trackError2+vertexError2);
    
    return pair<bool,Measurement1D>(true,Measurement1D(decayLen,decayLenError));

  }
std::pair<bool,Measurement1D> IPTools::signedDecayLength3D ( const reco::TransientTrack transientTrack,
const GlobalVector direction,
const reco::Vertex vertex 
) [inline]

Definition at line 84 of file IPTools.h.

References closestApproachToJet(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), and signedDecayLength3D().

        {
          // extrapolate to the point of closest approach to the jet axis
          TrajectoryStateOnSurface closestToJetState = closestApproachToJet(transientTrack.impactPointState(), vertex, direction,transientTrack.field());
         return signedDecayLength3D(closestToJetState,direction,vertex);
        }
pair< bool, Measurement1D > IPTools::signedImpactParameter3D ( const reco::TransientTrack track,
const GlobalVector direction,
const reco::Vertex vertex 
)

Returns life time signed 3D impact parameter The track is extrapolated to the closest point to the primary vertex in 3d space then the impact parameter and its error are computed

Definition at line 71 of file IPTools.cc.

References absoluteImpactParameter(), RecoVertex::convertPos(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), reco::Vertex::position(), CrabTask::prod, query::result, reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by pat::PATMuonProducer::embedHighLevel(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), TrackIPProducer::produce(), and SoftLepton::tag().

                                                                                                                  {
    //Extrapolate to closest point on transverse plane
    AnalyticalImpactPointExtrapolator extrapolator(track.field());
    TrajectoryStateOnSurface closestIn3DSpaceState = extrapolator.extrapolate(track.impactPointState(),RecoVertex::convertPos(vertex.position()));

    //Compute absolute value
    VertexDistance3D dist;
    pair<bool,Measurement1D> result = absoluteImpactParameter(closestIn3DSpaceState, vertex,  dist);
    if(!result.first) return result;

    //Compute Sign
    GlobalPoint impactPoint = closestIn3DSpaceState.globalPosition();
    GlobalVector IPVec(impactPoint.x()-vertex.x(),impactPoint.y()-vertex.y(),impactPoint.z()-vertex.z());
    double prod = IPVec.dot(direction);
    double sign = (prod>=0) ? 1. : -1.;

    //Apply sign to the result
    return pair<bool,Measurement1D>(result.first,Measurement1D(sign*result.second.value(), result.second.error()));
  }
pair< bool, Measurement1D > IPTools::signedTransverseImpactParameter ( const reco::TransientTrack track,
const GlobalVector direction,
const reco::Vertex vertex 
)

Returns life time signed transverse impact parameter The track is extrapolated to the closest point to the primary vertex in transverse plane then the impact parameter and its error are computed

Definition at line 50 of file IPTools.cc.

References absoluteImpactParameter(), RecoVertex::convertPos(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), reco::Vertex::position(), CrabTask::prod, query::result, reco::Vertex::x(), and reco::Vertex::y().

Referenced by CaloRecoTauAlgorithm::buildCaloTau(), PFRecoTauAlgorithm::buildPFTau(), HPSPFRecoTauAlgorithm::buildPFTau(), pat::PATMuonProducer::embedHighLevel(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), reco::tau::RecoTauImpactParameterSignificancePlugin::operator()(), TrackIPProducer::produce(), ConvBremPFTrackFinder::runConvBremFinder(), and SoftLepton::tag().

                                                                                                                  {
    //Extrapolate to closest point on transverse plane
    TransverseImpactPointExtrapolator extrapolator(track.field());
    TrajectoryStateOnSurface closestOnTransversePlaneState = extrapolator.extrapolate(track.impactPointState(),RecoVertex::convertPos(vertex.position()));
    
    //Compute absolute value
    VertexDistanceXY dist;
    pair<bool,Measurement1D> result = absoluteImpactParameter(closestOnTransversePlaneState, vertex,  dist);
    if(!result.first) return result;

    //Compute Sign
    GlobalPoint impactPoint    = closestOnTransversePlaneState.globalPosition();
    GlobalVector IPVec(impactPoint.x()-vertex.x(),impactPoint.y()-vertex.y(),0.);
    double prod = IPVec.dot(direction);
    double sign = (prod>=0) ? 1. : -1.;
    
    //Apply sign to the result
    return pair<bool,Measurement1D>(result.first,Measurement1D(sign*result.second.value(), result.second.error()));
  }
TrajectoryStateOnSurface IPTools::transverseExtrapolate ( const TrajectoryStateOnSurface track,
const GlobalPoint vertexPosition,
const MagneticField field 
) [inline]

Definition at line 56 of file IPTools.h.

References TransverseImpactPointExtrapolator::extrapolate().

        {
                 TransverseImpactPointExtrapolator extrapolator(field);
                 return extrapolator.extrapolate(track, vertexPosition);
        }