test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 namespace IPTools
21 {
22  using namespace std;
23  using namespace reco;
24 
25 
26  std::pair<bool,Measurement1D> absoluteImpactParameter(const TrajectoryStateOnSurface & tsos , const reco::Vertex & vertex, VertexDistance & distanceComputer) {
27  if(!tsos.isValid()) {
28  return pair<bool,Measurement1D>(false,Measurement1D(0.,0.)) ;
29  }
30  GlobalPoint refPoint = tsos.globalPosition();
31  GlobalError refPointErr = tsos.cartesianError().position();
32  GlobalPoint vertexPosition = RecoVertex::convertPos(vertex.position());
33  GlobalError vertexPositionErr = RecoVertex::convertError(vertex.error());
34  return pair<bool,Measurement1D>(true,distanceComputer.distance(VertexState(vertexPosition,vertexPositionErr), VertexState(refPoint, refPointErr)));
35  }
36 
37  std::pair<bool,Measurement1D> absoluteImpactParameter3D(const reco::TransientTrack & transientTrack, const reco::Vertex & vertex)
38  {
39  AnalyticalImpactPointExtrapolator extrapolator(transientTrack.field());
40  VertexDistance3D dist;
41  return absoluteImpactParameter(extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())), vertex, dist);
42  }
43  std::pair<bool,Measurement1D> absoluteTransverseImpactParameter(const reco::TransientTrack & transientTrack, const reco::Vertex & vertex)
44  {
45  TransverseImpactPointExtrapolator extrapolator(transientTrack.field());
46  VertexDistanceXY dist;
47  return absoluteImpactParameter(extrapolator.extrapolate(transientTrack.impactPointState(), RecoVertex::convertPos(vertex.position())), vertex, dist);
48  }
49 
50  pair<bool,Measurement1D> signedTransverseImpactParameter(const TransientTrack & track,
51  const GlobalVector & direction, const Vertex & vertex){
52  //Extrapolate to closest point on transverse plane
53  TransverseImpactPointExtrapolator extrapolator(track.field());
54  TrajectoryStateOnSurface closestOnTransversePlaneState = extrapolator.extrapolate(track.impactPointState(),RecoVertex::convertPos(vertex.position()));
55 
56  //Compute absolute value
57  VertexDistanceXY dist;
58  pair<bool,Measurement1D> result = absoluteImpactParameter(closestOnTransversePlaneState, vertex, dist);
59  if(!result.first) return result;
60 
61  //Compute Sign
62  GlobalPoint impactPoint = closestOnTransversePlaneState.globalPosition();
63  GlobalVector IPVec(impactPoint.x()-vertex.x(),impactPoint.y()-vertex.y(),0.);
64  double prod = IPVec.dot(direction);
65  double sign = (prod>=0) ? 1. : -1.;
66 
67  //Apply sign to the result
68  return pair<bool,Measurement1D>(result.first,Measurement1D(sign*result.second.value(), result.second.error()));
69  }
70 
71  pair<bool,Measurement1D> signedImpactParameter3D(const TransientTrack & track,
72  const GlobalVector & direction, const Vertex & vertex){
73  //Extrapolate to closest point on transverse plane
74  AnalyticalImpactPointExtrapolator extrapolator(track.field());
75  TrajectoryStateOnSurface closestIn3DSpaceState = extrapolator.extrapolate(track.impactPointState(),RecoVertex::convertPos(vertex.position()));
76 
77  //Compute absolute value
78  VertexDistance3D dist;
79  pair<bool,Measurement1D> result = absoluteImpactParameter(closestIn3DSpaceState, vertex, dist);
80  if(!result.first) return result;
81 
82  //Compute Sign
83  GlobalPoint impactPoint = closestIn3DSpaceState.globalPosition();
84  GlobalVector IPVec(impactPoint.x()-vertex.x(),impactPoint.y()-vertex.y(),impactPoint.z()-vertex.z());
85  double prod = IPVec.dot(direction);
86  double sign = (prod>=0) ? 1. : -1.;
87 
88  //Apply sign to the result
89  return pair<bool,Measurement1D>(result.first,Measurement1D(sign*result.second.value(), result.second.error()));
90  }
91 
92 
93 
94  pair<bool,Measurement1D> signedDecayLength3D(const TrajectoryStateOnSurface & closestToJetState,
95  const GlobalVector & direction, const Vertex & vertex) {
96 
97  //Check if extrapolation has been successfull
98  if(!closestToJetState.isValid()) {
99  return pair<bool,Measurement1D>(false,Measurement1D(0.,0.));
100  }
101 
102  GlobalVector jetDirection = direction.unit();
103  GlobalPoint vertexPosition(vertex.x(),vertex.y(),vertex.z());
104 
105  double decayLen = jetDirection.dot(closestToJetState.globalPosition()-vertexPosition);
106 
107  //error calculation
108 
110  j[0] = jetDirection.x();
111  j[1] = jetDirection.y();
112  j[2] = jetDirection.z();
114  jj[0] = jetDirection.x();
115  jj[1] = jetDirection.y();
116  jj[2] = jetDirection.z();
117  jj[3] =0.;
118  jj[4] =0.;
119  jj[5] =0.;
120 
121  //TODO: FIXME: the extrapolation uncertainty is very relevant here should be taken into account!!
122  double trackError2 = ROOT::Math::Similarity(jj,closestToJetState.cartesianError().matrix());
123  double vertexError2 = ROOT::Math::Similarity(j,vertex.covariance());
124 
125  double decayLenError = sqrt(trackError2+vertexError2);
126 
127  return pair<bool,Measurement1D>(true,Measurement1D(decayLen,decayLenError));
128 
129  }
130 
131 
132 
133  pair<bool,Measurement1D> linearizedSignedImpactParameter3D(const TrajectoryStateOnSurface & closestToJetState ,
134  const GlobalVector & direction, const Vertex & vertex)
135  {
136  //Check if extrapolation has been successfull
137  if(!closestToJetState.isValid()) {
138  return pair<bool,Measurement1D>(false,Measurement1D(0.,0.));
139  }
140 
141  GlobalPoint vertexPosition(vertex.x(),vertex.y(),vertex.z());
142  GlobalVector impactParameter = linearImpactParameter(closestToJetState, vertexPosition);
143  GlobalVector jetDir = direction.unit();
144  GlobalVector flightDistance(closestToJetState.globalPosition()-vertexPosition);
145  double theDistanceAlongJetAxis = jetDir.dot(flightDistance);
146  double signedIP = impactParameter.mag()*((theDistanceAlongJetAxis!=0)?theDistanceAlongJetAxis/abs(theDistanceAlongJetAxis):1.);
147 
148 
149  GlobalVector ipDirection = impactParameter.unit();
150  //GlobalPoint closestPoint = closestToJetState.globalPosition();
151  GlobalVector momentumAtClosestPoint = closestToJetState.globalMomentum();
152  GlobalVector momentumDir = momentumAtClosestPoint.unit();
153 
154  AlgebraicVector3 deriv_v;
155  deriv_v[0] = - ipDirection.x();
156  deriv_v[1] = - ipDirection.y();
157  deriv_v[2] = - ipDirection.z();
158 
159 
160  AlgebraicVector6 deriv;
161  deriv[0] = ipDirection.x();
162  deriv[1] = ipDirection.y();
163  deriv[2] = ipDirection.z();
164  deriv[3] = - (momentumDir.dot(flightDistance)*ipDirection.x())/momentumAtClosestPoint.mag();
165  deriv[4] = - (momentumDir.dot(flightDistance)*ipDirection.y())/momentumAtClosestPoint.mag();
166  deriv[5] = - (momentumDir.dot(flightDistance)*ipDirection.z())/momentumAtClosestPoint.mag();
167 
168  double trackError2 = ROOT::Math::Similarity(deriv , closestToJetState.cartesianError().matrix());
169  double vertexError2 = ROOT::Math::Similarity(deriv_v , vertex.covariance());
170  double ipError = sqrt(trackError2+vertexError2);
171 
172  return pair<bool,Measurement1D>(true,Measurement1D(signedIP,ipError));
173  }
174 
175 
176 
177  TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface & state,const Vertex & vertex, const GlobalVector& direction,const MagneticField * field) {
178 
179  Line::PositionType pos(GlobalPoint(vertex.x(),vertex.y(),vertex.z()));
180  Line::DirectionType dir(direction.unit());
181  Line jetLine(pos,dir);
182 
183  AnalyticalTrajectoryExtrapolatorToLine extrapolator(field);
184 
185  return extrapolator.extrapolate(state, jetLine);
186  }
187 
192 
193  Line::PositionType pos(state.globalPosition());
195  Line trackLine(pos,dir);
197  return trackLine.distance(tmp);
198  }
199 
200  pair<double,Measurement1D> jetTrackDistance(const TransientTrack & track, const GlobalVector & direction, const Vertex & vertex) {
201  double theLDist_err(0.);
202 
203  //FIXME
204  float weight=0.;//vertex.trackWeight(track);
205 
206  TrajectoryStateOnSurface stateAtOrigin = track.impactPointState();
207  if(!stateAtOrigin.isValid())
208  {
209  //TODO: throw instead?
210  return pair<bool,Measurement1D>(false,Measurement1D(0.,0.));
211  }
212 
213  //get the Track line at origin
214  Line::PositionType posTrack(stateAtOrigin.globalPosition());
215  Line::DirectionType dirTrack((stateAtOrigin.globalMomentum()).unit());
216  Line trackLine(posTrack,dirTrack);
217  // get the Jet line
218  // Vertex vertex(vertex);
219  GlobalVector jetVector = direction.unit();
220  Line::PositionType posJet(GlobalPoint(vertex.x(),vertex.y(),vertex.z()));
221  Line::DirectionType dirJet(jetVector);
222  Line jetLine(posJet,dirJet);
223 
224  // now compute the distance between the two lines
225  // If the track has been used to refit the Primary vertex then sign it positively, otherwise negative
226  double theDistanceToJetAxis = (jetLine.distance(trackLine)).mag();
227  if (weight<1) theDistanceToJetAxis= -theDistanceToJetAxis;
228 
229  // ... and the flight distance along the Jet axis.
230  GlobalPoint V = jetLine.position();
231  GlobalVector Q = dirTrack - jetVector.dot(dirTrack) * jetVector;
232  GlobalVector P = jetVector - jetVector.dot(dirTrack) * dirTrack;
233  double theDistanceAlongJetAxis = P.dot(V-posTrack)/Q.dot(dirTrack);
234 
235  //
236  // get the covariance matrix of the vertex and compute the error on theDistanceToJetAxis
237  //
238 
240 
241  // build the vector of closest approach between lines
242 
243 
244  //FIXME: error not computed.
245  GlobalVector H((jetVector.cross(dirTrack).unit()));
246  CLHEP::HepVector Hh(3);
247  Hh[0] = H.x();
248  Hh[1] = H.y();
249  Hh[2] = H.z();
250 
251  // theLDist_err = sqrt(vertexError.similarity(Hh));
252 
253  // cout << "distance to jet axis : "<< theDistanceToJetAxis <<" and error : "<< theLDist_err<<endl;
254  // Now the impact parameter ...
255 
256  /* GlobalPoint T0 = track.position();
257  GlobalVector D = (T0-V)- (T0-V).dot(dir) * dir;
258  double IP = D.mag();
259  GlobalVector Dold = distance(aTSOS, aJet.vertex(), jetDirection);
260  double IPold = Dold.mag();
261  */
262 
263 
264 
265 
266  Measurement1D DTJA(theDistanceToJetAxis,theLDist_err);
267 
268  return pair<double,Measurement1D> (theDistanceAlongJetAxis,DTJA);
269  }
270 
271 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
GlobalVector linearImpactParameter(const TrajectoryStateOnSurface &aTSOS, const GlobalPoint &point)
Definition: IPTools.cc:191
Definition: Line.h:10
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
double y() const
y coordinate
Definition: Vertex.h:103
double sign(double x)
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:50
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::pair< bool, Measurement1D > signedDecayLength3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:94
const CartesianTrajectoryError cartesianError() const
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &fts, const Line &L) const
extrapolation from FreeTrajectoryState
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:71
GlobalPoint globalPosition() const
#define P
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:177
const MagneticField * field() const
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:116
const Point & position() const
position
Definition: Vertex.h:99
ROOT::Math::SVector< double, 6 > AlgebraicVector6
GlobalVector distance(const Line &aLine) const
Definition: Line.h:38
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:200
tuple result
Definition: mps_fire.py:83
PositionType position() const
Definition: Line.h:26
string unit
Definition: csvLumiCalc.py:46
T sqrt(T t)
Definition: SSEVec.h:18
ROOT::Math::SVector< double, 3 > AlgebraicVector3
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
double z() const
y coordinate
Definition: Vertex.h:105
const AlgebraicSymMatrix66 & matrix() const
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:27
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double x() const
x coordinate
Definition: Vertex.h:101
Measurement1D distance(const reco::Vertex &, const reco::Vertex &) const
const GlobalError position() const
Position error submatrix.
Error error() const
return SMatrix
Definition: Vertex.h:122
std::pair< bool, Measurement1D > absoluteTransverseImpactParameter(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:43
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:17
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
GlobalVector globalMomentum() const
std::pair< bool, Measurement1D > linearizedSignedImpactParameter3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:133
TrajectoryStateOnSurface impactPointState() const
dbl *** dir
Definition: mlp_gen.cc:35
std::pair< bool, Measurement1D > absoluteImpactParameter(const TrajectoryStateOnSurface &tsos, const reco::Vertex &vertex, VertexDistance &distanceComputer)
Impact parameter without direction (internally used)
Definition: IPTools.cc:26
T x() const
Definition: PV3DBase.h:62
*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