CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Functions
IPTools Namespace Reference

Functions

std::pair< bool, Measurement1DabsoluteImpactParameter (const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
 Impact parameter without direction (internally used) More...
 
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 
)

Impact parameter without direction (internally used)

Definition at line 23 of file IPTools.cc.

References TrajectoryStateOnSurface::cartesianError(), RecoVertex::convertError(), RecoVertex::convertPos(), VertexDistance::distance(), reco::Vertex::error(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), CartesianTrajectoryError::position(), and reco::Vertex::position().

Referenced by absoluteImpactParameter3D(), absoluteTransverseImpactParameter(), signedImpactParameter3D(), and signedTransverseImpactParameter().

25  {
26  if (!tsos.isValid()) {
27  return pair<bool, Measurement1D>(false, Measurement1D(0., 0.));
28  }
29  GlobalPoint refPoint = tsos.globalPosition();
30  GlobalError refPointErr = tsos.cartesianError().position();
31  GlobalPoint vertexPosition = RecoVertex::convertPos(vertex.position());
32  GlobalError vertexPositionErr = RecoVertex::convertError(vertex.error());
33  return pair<bool, Measurement1D>(
34  true,
35  distanceComputer.distance(VertexState(vertexPosition, vertexPositionErr), VertexState(refPoint, refPointErr)));
36  }
reco::Vertex::Point convertPos(const GlobalPoint &p)
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
const CartesianTrajectoryError cartesianError() const
GlobalPoint globalPosition() const
const Point & position() const
position
Definition: Vertex.h:127
Measurement1D distance(const reco::Vertex &, const reco::Vertex &) const
const GlobalError position() const
Position error submatrix.
Error error() const
return SMatrix
Definition: Vertex.h:163
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 38 of file IPTools.cc.

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

Referenced by PrimaryVertexValidation::analyze(), TracksClusteringFromDisplacedSeed::clusters(), PFCand_AssoMapAlgos::createMappings(), PF_PU_AssoMapAlgos::createMappings(), PF_PU_AssoMapAlgos::FindClosest3D(), EGammaMvaEleEstimator::IDIsoCombinedMvaValue(), EGammaMvaEleEstimatorCSA14::mvaValue(), EGammaMvaEleEstimator::mvaValue(), TracksClusteringFromDisplacedSeed::nearTracks(), Onia2MuMuPAT::produce(), pat::PATElectronProducer::produce(), btagbtvdeep::seedingTracksToFeatures(), and TrackVertexArbitration< VTX >::trackVertexArbitrator().

39  {
40  AnalyticalImpactPointExtrapolator extrapolator(transientTrack.field());
41  VertexDistance3D dist;
43  extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())),
44  vertex,
45  dist);
46  }
reco::Vertex::Point convertPos(const GlobalPoint &p)
const MagneticField * field() const
const Point & position() const
position
Definition: Vertex.h:127
TrajectoryStateOnSurface impactPointState() const
std::pair< bool, Measurement1D > absoluteImpactParameter(const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
Impact parameter without direction (internally used)
Definition: IPTools.cc:23
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 47 of file IPTools.cc.

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

Referenced by btagbtvdeep::seedingTracksToFeatures().

48  {
49  TransverseImpactPointExtrapolator extrapolator(transientTrack.field());
50  VertexDistanceXY dist;
52  extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())),
53  vertex,
54  dist);
55  }
reco::Vertex::Point convertPos(const GlobalPoint &p)
const MagneticField * field() const
const Point & position() const
position
Definition: Vertex.h:127
TrajectoryStateOnSurface impactPointState() const
std::pair< bool, Measurement1D > absoluteImpactParameter(const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
Impact parameter without direction (internally used)
Definition: IPTools.cc:23
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 206 of file IPTools.cc.

References Line::distance(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), data-class-funcs::H, reco::TransientTrack::impactPointState(), TrajectoryStateOnSurface::isValid(), mag(), Line::position(), unit(), cms::cuda::V, histoStyle::weight, reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by BDHadronTrackMonitoringAnalyzer::analyze(), btagbtvdeep::TrackInfoBuilder::buildTrackInfo(), PrimaryVertexAssignment::chargedHadronVertex(), DeepBoostedJetTagInfoProducer::fillParticleFeatures(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), PixelJetPuId::produce(), BoostedDoubleSVProducer::produce(), IPProducer< Container, Base, Helper >::produce(), and btagbtvdeep::seedingTracksToFeatures().

208  {
209  double theLDist_err(0.);
210 
211  //FIXME
212  float weight = 0.; //vertex.trackWeight(track);
213 
214  TrajectoryStateOnSurface stateAtOrigin = track.impactPointState();
215  if (!stateAtOrigin.isValid()) {
216  //TODO: throw instead?
217  return pair<bool, Measurement1D>(false, Measurement1D(0., 0.));
218  }
219 
220  //get the Track line at origin
221  Line::PositionType posTrack(stateAtOrigin.globalPosition());
222  Line::DirectionType dirTrack((stateAtOrigin.globalMomentum()).unit());
223  Line trackLine(posTrack, dirTrack);
224  // get the Jet line
225  // Vertex vertex(vertex);
226  GlobalVector jetVector = direction.unit();
227  Line::PositionType posJet(GlobalPoint(vertex.x(), vertex.y(), vertex.z()));
228  Line::DirectionType dirJet(jetVector);
229  Line jetLine(posJet, dirJet);
230 
231  // now compute the distance between the two lines
232  // If the track has been used to refit the Primary vertex then sign it positively, otherwise negative
233  double theDistanceToJetAxis = (jetLine.distance(trackLine)).mag();
234  if (weight < 1)
235  theDistanceToJetAxis = -theDistanceToJetAxis;
236 
237  // ... and the flight distance along the Jet axis.
238  GlobalPoint V = jetLine.position();
239  GlobalVector Q = dirTrack - jetVector.dot(dirTrack) * jetVector;
240  GlobalVector P = jetVector - jetVector.dot(dirTrack) * dirTrack;
241  double theDistanceAlongJetAxis = P.dot(V - posTrack) / Q.dot(dirTrack);
242 
243  //
244  // get the covariance matrix of the vertex and compute the error on theDistanceToJetAxis
245  //
246 
248 
249  // build the vector of closest approach between lines
250 
251  //FIXME: error not computed.
252  GlobalVector H((jetVector.cross(dirTrack).unit()));
253  CLHEP::HepVector Hh(3);
254  Hh[0] = H.x();
255  Hh[1] = H.y();
256  Hh[2] = H.z();
257 
258  // theLDist_err = sqrt(vertexError.similarity(Hh));
259 
260  // cout << "distance to jet axis : "<< theDistanceToJetAxis <<" and error : "<< theLDist_err<<endl;
261  // Now the impact parameter ...
262 
263  /* GlobalPoint T0 = track.position();
264  GlobalVector D = (T0-V)- (T0-V).dot(dir) * dir;
265  double IP = D.mag();
266  GlobalVector Dold = distance(aTSOS, aJet.vertex(), jetDirection);
267  double IPold = Dold.mag();
268  */
269 
270  Measurement1D DTJA(theDistanceToJetAxis, theLDist_err);
271 
272  return pair<double, Measurement1D>(theDistanceAlongJetAxis, DTJA);
273  }
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:110
Definition: Line.h:10
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
double y() const
y coordinate
Definition: Vertex.h:131
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
GlobalPoint globalPosition() const
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t V
double z() const
z coordinate
Definition: Vertex.h:133
Vector3DBase unit() const
Definition: Vector3DBase.h:54
double x() const
x coordinate
Definition: Vertex.h:129
std::pair< OmniClusterRef, TrackingParticleRef > P
GlobalVector globalMomentum() const
TrajectoryStateOnSurface impactPointState() const
int weight
Definition: histoStyle.py:51
Basic3DVector unit() const
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 198 of file IPTools.cc.

References DeadROC_duringRun::dir, TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), point, createJobs::tmp, and unit().

Referenced by linearizedSignedImpactParameter3D().

198  {
199  Line::PositionType pos(state.globalPosition());
201  Line trackLine(pos, dir);
202  const GlobalPoint& tmp = point;
203  return trackLine.distance(tmp);
204  }
Definition: Line.h:10
GlobalPoint globalPosition() const
GlobalVector globalMomentum() const
tmp
align.sh
Definition: createJobs.py:716
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
Basic3DVector unit() const
pair< bool, Measurement1D > IPTools::linearizedSignedImpactParameter3D ( const TrajectoryStateOnSurface state,
const GlobalVector direction,
const reco::Vertex vertex 
)

Definition at line 141 of file IPTools.cc.

References funct::abs(), TrajectoryStateOnSurface::cartesianError(), reco::Vertex::covariance(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), linearImpactParameter(), CartesianTrajectoryError::matrix(), mathSSE::sqrt(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by linearizedSignedImpactParameter3D().

143  {
144  //Check if extrapolation has been successfull
145  if (!closestToJetState.isValid()) {
146  return pair<bool, Measurement1D>(false, Measurement1D(0., 0.));
147  }
148 
149  GlobalPoint vertexPosition(vertex.x(), vertex.y(), vertex.z());
150  GlobalVector impactParameter = linearImpactParameter(closestToJetState, vertexPosition);
151  GlobalVector jetDir = direction.unit();
152  GlobalVector flightDistance(closestToJetState.globalPosition() - vertexPosition);
153  double theDistanceAlongJetAxis = jetDir.dot(flightDistance);
154  double signedIP = impactParameter.mag() *
155  ((theDistanceAlongJetAxis != 0) ? theDistanceAlongJetAxis / abs(theDistanceAlongJetAxis) : 1.);
156 
157  GlobalVector ipDirection = impactParameter.unit();
158  //GlobalPoint closestPoint = closestToJetState.globalPosition();
159  GlobalVector momentumAtClosestPoint = closestToJetState.globalMomentum();
160  GlobalVector momentumDir = momentumAtClosestPoint.unit();
161 
162  AlgebraicVector3 deriv_v;
163  deriv_v[0] = -ipDirection.x();
164  deriv_v[1] = -ipDirection.y();
165  deriv_v[2] = -ipDirection.z();
166 
167  AlgebraicVector6 deriv;
168  deriv[0] = ipDirection.x();
169  deriv[1] = ipDirection.y();
170  deriv[2] = ipDirection.z();
171  deriv[3] = -(momentumDir.dot(flightDistance) * ipDirection.x()) / momentumAtClosestPoint.mag();
172  deriv[4] = -(momentumDir.dot(flightDistance) * ipDirection.y()) / momentumAtClosestPoint.mag();
173  deriv[5] = -(momentumDir.dot(flightDistance) * ipDirection.z()) / momentumAtClosestPoint.mag();
174 
175  double trackError2 = ROOT::Math::Similarity(deriv, closestToJetState.cartesianError().matrix());
176  double vertexError2 = ROOT::Math::Similarity(deriv_v, vertex.covariance());
177  double ipError = sqrt(trackError2 + vertexError2);
178 
179  return pair<bool, Measurement1D>(true, Measurement1D(signedIP, ipError));
180  }
GlobalVector linearImpactParameter(const TrajectoryStateOnSurface &aTSOS, const GlobalPoint &point)
Definition: IPTools.cc:198
double y() const
y coordinate
Definition: Vertex.h:131
T y() const
Definition: PV3DBase.h:60
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:148
T mag() const
Definition: PV3DBase.h:64
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
z coordinate
Definition: Vertex.h:133
Vector3DBase unit() const
Definition: Vector3DBase.h:54
double x() const
x coordinate
Definition: Vertex.h:129
ROOT::Math::SVector< double, 3 > AlgebraicVector3
ROOT::Math::SVector< double, 6 > AlgebraicVector6
T x() const
Definition: PV3DBase.h:59
std::pair<bool, Measurement1D> IPTools::linearizedSignedImpactParameter3D ( const reco::TransientTrack transientTrack,
const GlobalVector direction,
const reco::Vertex vertex 
)
inline

Definition at line 76 of file IPTools.h.

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

78  {
79  // extrapolate to the point of closest approach to the jet axis
80  TrajectoryStateOnSurface closestToJetState =
81  closestApproachToJet(transientTrack.impactPointState(), vertex, direction, transientTrack.field());
82  return linearizedSignedImpactParameter3D(closestToJetState, direction, vertex);
83  }
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:182
const MagneticField * field() const
std::pair< bool, Measurement1D > linearizedSignedImpactParameter3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:141
TrajectoryStateOnSurface impactPointState() const
pair< bool, Measurement1D > IPTools::signedDecayLength3D ( const TrajectoryStateOnSurface state,
const GlobalVector direction,
const reco::Vertex vertex 
)

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

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

Definition at line 105 of file IPTools.cc.

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

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

107  {
108  //Check if extrapolation has been successfull
109  if (!closestToJetState.isValid()) {
110  return pair<bool, Measurement1D>(false, Measurement1D(0., 0.));
111  }
112 
113  GlobalVector jetDirection = direction.unit();
114  GlobalPoint vertexPosition(vertex.x(), vertex.y(), vertex.z());
115 
116  double decayLen = jetDirection.dot(closestToJetState.globalPosition() - vertexPosition);
117 
118  //error calculation
119 
121  j[0] = jetDirection.x();
122  j[1] = jetDirection.y();
123  j[2] = jetDirection.z();
125  jj[0] = jetDirection.x();
126  jj[1] = jetDirection.y();
127  jj[2] = jetDirection.z();
128  jj[3] = 0.;
129  jj[4] = 0.;
130  jj[5] = 0.;
131 
132  //TODO: FIXME: the extrapolation uncertainty is very relevant here should be taken into account!!
133  double trackError2 = ROOT::Math::Similarity(jj, closestToJetState.cartesianError().matrix());
134  double vertexError2 = ROOT::Math::Similarity(j, vertex.covariance());
135 
136  double decayLenError = sqrt(trackError2 + vertexError2);
137 
138  return pair<bool, Measurement1D>(true, Measurement1D(decayLen, decayLenError));
139  }
double y() const
y coordinate
Definition: Vertex.h:131
T y() const
Definition: PV3DBase.h:60
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:148
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
double z() const
z coordinate
Definition: Vertex.h:133
Vector3DBase unit() const
Definition: Vector3DBase.h:54
double x() const
x coordinate
Definition: Vertex.h:129
ROOT::Math::SVector< double, 3 > AlgebraicVector3
ROOT::Math::SVector< double, 6 > AlgebraicVector6
T x() const
Definition: PV3DBase.h:59
std::pair<bool, Measurement1D> IPTools::signedDecayLength3D ( const reco::TransientTrack transientTrack,
const GlobalVector direction,
const reco::Vertex vertex 
)
inline

Definition at line 89 of file IPTools.h.

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

91  {
92  // extrapolate to the point of closest approach to the jet axis
93  TrajectoryStateOnSurface closestToJetState =
94  closestApproachToJet(transientTrack.impactPointState(), vertex, direction, transientTrack.field());
95  return signedDecayLength3D(closestToJetState, direction, vertex);
96  }
std::pair< bool, Measurement1D > signedDecayLength3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:105
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:182
const MagneticField * field() const
TrajectoryStateOnSurface impactPointState() const
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 81 of file IPTools.cc.

References absoluteImpactParameter(), RecoVertex::convertPos(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), reco::Vertex::position(), mps_fire::result, jetcorrextractor::sign(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by tadqm::TrackAnalyzer::analyze(), PrimaryVertexValidation::analyze(), btagbtvdeep::SeedingTrackInfoBuilder::buildSeedingTrackInfo(), btagbtvdeep::TrackInfoBuilder::buildTrackInfo(), pat::PATMuonProducer::embedHighLevel(), pat::PATElectronProducer::embedHighLevel(), DeepBoostedJetTagInfoProducer::fillParticleFeatures(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), reco::tau::RecoTauImpactParameterSignificancePlugin::operator()(), PFTauTransverseImpactParameters::produce(), SoftPFElectronTagInfoProducer::produce(), SoftPFMuonTagInfoProducer::produce(), IPProducer< Container, Base, Helper >::produce(), and SoftLepton::tag().

83  {
84  //Extrapolate to closest point on transverse plane
85  AnalyticalImpactPointExtrapolator extrapolator(track.field());
86  TrajectoryStateOnSurface closestIn3DSpaceState =
87  extrapolator.extrapolate(track.impactPointState(), RecoVertex::convertPos(vertex.position()));
88 
89  //Compute absolute value
90  VertexDistance3D dist;
91  pair<bool, Measurement1D> result = absoluteImpactParameter(closestIn3DSpaceState, vertex, dist);
92  if (!result.first)
93  return result;
94 
95  //Compute Sign
96  GlobalPoint impactPoint = closestIn3DSpaceState.globalPosition();
97  GlobalVector IPVec(impactPoint.x() - vertex.x(), impactPoint.y() - vertex.y(), impactPoint.z() - vertex.z());
98  double prod = IPVec.dot(direction);
99  double sign = (prod >= 0) ? 1. : -1.;
100 
101  //Apply sign to the result
102  return pair<bool, Measurement1D>(result.first, Measurement1D(sign * result.second.value(), result.second.error()));
103  }
reco::Vertex::Point convertPos(const GlobalPoint &p)
double y() const
y coordinate
Definition: Vertex.h:131
double sign(double x)
T y() const
Definition: PV3DBase.h:60
const MagneticField * field() const
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
const Point & position() const
position
Definition: Vertex.h:127
tuple result
Definition: mps_fire.py:311
T z() const
Definition: PV3DBase.h:61
double z() const
z coordinate
Definition: Vertex.h:133
double x() const
x coordinate
Definition: Vertex.h:129
TrajectoryStateOnSurface impactPointState() const
std::pair< bool, Measurement1D > absoluteImpactParameter(const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
Impact parameter without direction (internally used)
Definition: IPTools.cc:23
T x() const
Definition: PV3DBase.h:59
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 57 of file IPTools.cc.

References absoluteImpactParameter(), RecoVertex::convertPos(), reco::TransientTrack::field(), reco::TransientTrack::impactPointState(), reco::Vertex::position(), mps_fire::result, jetcorrextractor::sign(), reco::Vertex::x(), and reco::Vertex::y().

Referenced by tadqm::TrackAnalyzer::analyze(), PrimaryVertexValidation::analyze(), btagbtvdeep::SeedingTrackInfoBuilder::buildSeedingTrackInfo(), btagbtvdeep::TrackInfoBuilder::buildTrackInfo(), DeepBoostedJetTagInfoProducer::fillParticleFeatures(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), PFConversionProducer::produce(), PFTrackProducer::produce(), PFTauTransverseImpactParameters::produce(), SoftPFElectronTagInfoProducer::produce(), SoftPFMuonTagInfoProducer::produce(), IPProducer< Container, Base, Helper >::produce(), ConvBremPFTrackFinder::runConvBremFinder(), and SoftLepton::tag().

59  {
60  //Extrapolate to closest point on transverse plane
61  TransverseImpactPointExtrapolator extrapolator(track.field());
62  TrajectoryStateOnSurface closestOnTransversePlaneState =
63  extrapolator.extrapolate(track.impactPointState(), RecoVertex::convertPos(vertex.position()));
64 
65  //Compute absolute value
66  VertexDistanceXY dist;
67  pair<bool, Measurement1D> result = absoluteImpactParameter(closestOnTransversePlaneState, vertex, dist);
68  if (!result.first)
69  return result;
70 
71  //Compute Sign
72  GlobalPoint impactPoint = closestOnTransversePlaneState.globalPosition();
73  GlobalVector IPVec(impactPoint.x() - vertex.x(), impactPoint.y() - vertex.y(), 0.);
74  double prod = IPVec.dot(direction);
75  double sign = (prod >= 0) ? 1. : -1.;
76 
77  //Apply sign to the result
78  return pair<bool, Measurement1D>(result.first, Measurement1D(sign * result.second.value(), result.second.error()));
79  }
reco::Vertex::Point convertPos(const GlobalPoint &p)
double y() const
y coordinate
Definition: Vertex.h:131
double sign(double x)
T y() const
Definition: PV3DBase.h:60
const MagneticField * field() const
const Point & position() const
position
Definition: Vertex.h:127
tuple result
Definition: mps_fire.py:311
double x() const
x coordinate
Definition: Vertex.h:129
TrajectoryStateOnSurface impactPointState() const
std::pair< bool, Measurement1D > absoluteImpactParameter(const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
Impact parameter without direction (internally used)
Definition: IPTools.cc:23
T x() const
Definition: PV3DBase.h:59
TrajectoryStateOnSurface IPTools::transverseExtrapolate ( const TrajectoryStateOnSurface track,
const GlobalPoint vertexPosition,
const MagneticField field 
)
inline

Definition at line 58 of file IPTools.h.

References TransverseImpactPointExtrapolator::extrapolate().

60  {
61  TransverseImpactPointExtrapolator extrapolator(field);
62  return extrapolator.extrapolate(track, vertexPosition);
63  }