CMS 3D CMS Logo

IPTools.cc
Go to the documentation of this file.
2 
9 #include "CLHEP/Vector/ThreeVector.h"
10 #include "CLHEP/Vector/LorentzVector.h"
11 #include "CLHEP/Matrix/Vector.h"
12 #include <string>
13 
18 
19 namespace IPTools {
20  using namespace std;
21  using namespace reco;
22 
23  std::pair<bool, Measurement1D> absoluteImpactParameter(const TrajectoryStateOnSurface& tsos,
24  const reco::Vertex& vertex,
25  VertexDistance& distanceComputer) {
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  }
37 
38  std::pair<bool, Measurement1D> absoluteImpactParameter3D(const reco::TransientTrack& transientTrack,
39  const reco::Vertex& vertex) {
40  AnalyticalImpactPointExtrapolator extrapolator(transientTrack.field());
41  VertexDistance3D dist;
43  extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())),
44  vertex,
45  dist);
46  }
47  std::pair<bool, Measurement1D> absoluteTransverseImpactParameter(const reco::TransientTrack& transientTrack,
48  const reco::Vertex& vertex) {
49  TransverseImpactPointExtrapolator extrapolator(transientTrack.field());
50  VertexDistanceXY dist;
52  extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())),
53  vertex,
54  dist);
55  }
56 
57  pair<bool, Measurement1D> signedTransverseImpactParameter(const TransientTrack& track,
58  const GlobalVector& direction,
59  const Vertex& vertex) {
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  }
80 
81  pair<bool, Measurement1D> signedImpactParameter3D(const TransientTrack& track,
82  const GlobalVector& direction,
83  const Vertex& vertex) {
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  }
104 
105  pair<bool, Measurement1D> signedDecayLength3D(const TrajectoryStateOnSurface& closestToJetState,
106  const GlobalVector& direction,
107  const Vertex& vertex) {
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  }
140 
141  pair<bool, Measurement1D> linearizedSignedImpactParameter3D(const TrajectoryStateOnSurface& closestToJetState,
142  const GlobalVector& direction,
143  const Vertex& vertex) {
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  }
181 
183  const Vertex& vertex,
184  const GlobalVector& direction,
185  const MagneticField* field) {
187  Line::DirectionType dir(direction.unit());
188  Line jetLine(pos, dir);
189 
190  AnalyticalTrajectoryExtrapolatorToLine extrapolator(field);
191 
192  return extrapolator.extrapolate(state, jetLine);
193  }
194 
201  Line trackLine(pos, dir);
202  const GlobalPoint& tmp = point;
203  return trackLine.distance(tmp);
204  }
205 
206  pair<double, Measurement1D> jetTrackDistance(const TransientTrack& track,
207  const GlobalVector& direction,
208  const Vertex& vertex) {
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  }
274 
275 } // namespace IPTools
Vector3DBase
Definition: Vector3DBase.h:8
class-composition.H
H
Definition: class-composition.py:31
cms::cuda::V
cudaStream_t T uint32_t const T *__restrict__ const uint32_t *__restrict__ uint32_t int cudaStream_t V
Definition: HistoContainer.h:99
TrajectoryStateOnSurface.h
AlgebraicVector3
ROOT::Math::SVector< double, 3 > AlgebraicVector3
Definition: AlgebraicROOTObjects.h:12
Measurement1D
Definition: Measurement1D.h:11
MessageLogger.h
IPTools::signedDecayLength3D
std::pair< bool, Measurement1D > signedDecayLength3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:105
AnalyticalImpactPointExtrapolator.h
TrajectoryStateOnSurface::cartesianError
const CartesianTrajectoryError cartesianError() const
Definition: TrajectoryStateOnSurface.h:71
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
IPTools::absoluteImpactParameter3D
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:38
TrajectoryStateOnSurface::globalPosition
GlobalPoint globalPosition() const
Definition: TrajectoryStateOnSurface.h:65
pos
Definition: PixelAliasList.h:18
Measurement1D.h
IPTools::closestApproachToJet
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:182
Validation_hcalonly_cfi.sign
sign
Definition: Validation_hcalonly_cfi.py:32
ConvertToFromReco.h
VertexDistance3D.h
TransverseImpactPointExtrapolator
Definition: TransverseImpactPointExtrapolator.h:26
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
TransientTrack.h
IPTools::signedTransverseImpactParameter
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
class-composition.Q
Q
Definition: class-composition.py:82
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
Vector3DBase::unit
Vector3DBase unit() const
Definition: Vector3DBase.h:54
VertexDistance::distance
Measurement1D distance(const reco::Vertex &, const reco::Vertex &) const
Definition: VertexDistance.cc:26
VertexDistanceXY.h
IPTools::signedImpactParameter3D
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:81
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
dumpMFGeometry_cfg.prod
prod
Definition: dumpMFGeometry_cfg.py:24
VertexDistance3D
Definition: VertexDistance3D.h:13
Line.h
GlobalPoint
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Point3DBase< float, GlobalTag >
IPTools::jetTrackDistance
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:206
VertexDistance
Definition: VertexDistance.h:16
bphysicsOniaDQM_cfi.vertex
vertex
Definition: bphysicsOniaDQM_cfi.py:7
AnalyticalImpactPointExtrapolator
Definition: AnalyticalImpactPointExtrapolator.h:26
RecoVertex::convertError
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
math::GlobalVector
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:28
reco::TransientTrack::impactPointState
TrajectoryStateOnSurface impactPointState() const
Definition: TransientTrack.h:98
IPTools
Definition: IPTools.h:16
Line
Definition: Line.h:10
IPTools::linearizedSignedImpactParameter3D
std::pair< bool, Measurement1D > linearizedSignedImpactParameter3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:141
GlobalErrorBase< double, ErrorMatrixTag >
VertexDistanceXY
Definition: VertexDistanceXY.h:11
AlgebraicVector6
ROOT::Math::SVector< double, 6 > AlgebraicVector6
Definition: AlgebraicROOTObjects.h:15
reco::TransientTrack::field
const MagneticField * field() const
Definition: TransientTrack.h:110
CartesianTrajectoryError::matrix
const AlgebraicSymMatrix66 & matrix() const
Definition: CartesianTrajectoryError.h:28
unit
Basic3DVector unit() const
Definition: Basic3DVectorLD.h:162
Line::distance
GlobalVector distance(const Line &aLine) const
Definition: Line.h:35
TrajectoryStateOnSurface::globalMomentum
GlobalVector globalMomentum() const
Definition: TrajectoryStateOnSurface.h:66
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
std
Definition: JetResolutionObject.h:76
IPTools.h
reco::TransientTrack
Definition: TransientTrack.h:19
VertexState
Definition: VertexState.h:13
findQualityFiles.jj
string jj
Definition: findQualityFiles.py:188
Line::position
PositionType position() const
Definition: Line.h:24
math::GlobalPoint
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
RecoVertex::convertPos
reco::Vertex::Point convertPos(const GlobalPoint &p)
Definition: ConvertToFromReco.h:7
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
mps_fire.result
result
Definition: mps_fire.py:303
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
IPTools::absoluteTransverseImpactParameter
std::pair< bool, Measurement1D > absoluteTransverseImpactParameter(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:47
IPTools::linearImpactParameter
GlobalVector linearImpactParameter(const TrajectoryStateOnSurface &aTSOS, const GlobalPoint &point)
Definition: IPTools.cc:198
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
P
std::pair< OmniClusterRef, TrackingParticleRef > P
Definition: BDHadronTrackMonitoringAnalyzer.cc:202
point
*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
MagneticField
Definition: MagneticField.h:19
AnalyticalTrajectoryExtrapolatorToLine
Definition: AnalyticalTrajectoryExtrapolatorToLine.h:28
CartesianTrajectoryError::position
const GlobalError position() const
Position error submatrix.
Definition: CartesianTrajectoryError.h:41
AnalyticalTrajectoryExtrapolatorToLine::extrapolate
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &fts, const Line &L) const
extrapolation from FreeTrajectoryState
Definition: AnalyticalTrajectoryExtrapolatorToLine.cc:22
AnalyticalTrajectoryExtrapolatorToLine.h
weight
Definition: weight.py:1
reco::Vertex
Definition: Vertex.h:35
TrajectoryStateOnSurface::isValid
bool isValid() const
Definition: TrajectoryStateOnSurface.h:54
IPTools::absoluteImpactParameter
std::pair< bool, Measurement1D > absoluteImpactParameter(const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
Impact parameter without direction (internally used)
Definition: IPTools.cc:23
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23