#include <TrackingTools/TrajectoryState/interface/PerigeeConversions.h>
Public Member Functions | |
TrackCharge | chargeFromPerigee (const PerigeeTrajectoryParameters &perigee) const |
This method returns the charge. | |
CurvilinearTrajectoryError | curvilinearError (const PerigeeTrajectoryError &perigeeError, const GlobalTrajectoryParameters >p) const |
PerigeeTrajectoryError | ftsToPerigeeError (const FTS &originalFTS) const |
PerigeeTrajectoryParameters | ftsToPerigeeParameters (const FTS &originalFTS, const GlobalPoint &referencePoint, double &pt) const |
This method calculates the perigee parameters from a given FTS and a reference point. | |
AlgebraicMatrix55 | jacobianCurvilinear2Perigee (const FreeTrajectoryState &fts) const |
AlgebraicMatrix | jacobianCurvilinear2Perigee_old (const FreeTrajectoryState &fts) const |
Jacobians of tranformations between curvilinear frame at point of closest approach in transverse plane and perigee frame. | |
AlgebraicMatrix66 | jacobianParameters2Cartesian (const AlgebraicVector3 &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field) const |
Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian. | |
AlgebraicMatrix | jacobianParameters2Cartesian_old (const AlgebraicVector &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field) const |
Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian. | |
AlgebraicMatrix55 | jacobianPerigee2Curvilinear (const GlobalTrajectoryParameters >p) const |
AlgebraicMatrix | jacobianPerigee2Curvilinear_old (const GlobalTrajectoryParameters >p) const |
GlobalVector | momentumFromPerigee (const PerigeeTrajectoryParameters ¶meters, double pt, const GlobalPoint &referencePoint) const |
This method returns the (Cartesian) momentum from the PerigeeTrajectoryParameters. | |
GlobalVector | momentumFromPerigee (const AlgebraicVector3 &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field) const |
GlobalVector | momentumFromPerigee (const AlgebraicVector &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field) const |
This method returns the (Cartesian) momentum. | |
GlobalPoint | positionFromPerigee (const PerigeeTrajectoryParameters ¶meters, const GlobalPoint &referencePoint) const |
This method returns the position (on the helix) at which the parameters are defined. | |
TrajectoryStateClosestToPoint | trajectoryStateClosestToPoint (const AlgebraicVector3 &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicSymMatrix66 &theCovarianceMatrix, const MagneticField *field) const |
TrajectoryStateClosestToPoint | trajectoryStateClosestToPoint (const AlgebraicVector &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicMatrix &theCovarianceMatrix, const MagneticField *field) const |
Public constructor. | |
Private Types | |
typedef FreeTrajectoryState | FTS |
Definition at line 16 of file PerigeeConversions.h.
typedef FreeTrajectoryState PerigeeConversions::FTS [private] |
Definition at line 18 of file PerigeeConversions.h.
TrackCharge PerigeeConversions::chargeFromPerigee | ( | const PerigeeTrajectoryParameters & | perigee | ) | const |
This method returns the charge.
Definition at line 149 of file PerigeeConversions.cc.
References PerigeeTrajectoryParameters::charge().
Referenced by MuonStandaloneAlgorithm::produceFreeTrajectoryState().
00150 { 00151 return parameters.charge(); 00152 }
CurvilinearTrajectoryError PerigeeConversions::curvilinearError | ( | const PerigeeTrajectoryError & | perigeeError, | |
const GlobalTrajectoryParameters & | gtp | |||
) | const |
Definition at line 101 of file PerigeeConversions.cc.
References PerigeeTrajectoryError::covarianceMatrix().
Referenced by TrajectoryStateClosestToPoint::calculateFTS().
00102 { 00103 AlgebraicMatrix55 perigee2curv = jacobianPerigee2Curvilinear(gtp); 00104 return CurvilinearTrajectoryError(ROOT::Math::Similarity(perigee2curv, perigeeError.covarianceMatrix())); 00105 }
PerigeeTrajectoryError PerigeeConversions::ftsToPerigeeError | ( | const FTS & | originalFTS | ) | const |
Definition at line 62 of file PerigeeConversions.cc.
References FreeTrajectoryState::curvilinearError(), and CurvilinearTrajectoryError::matrix().
Referenced by TrackerSeedValidator::analyze(), and TrajectoryStateClosestToPoint::TrajectoryStateClosestToPoint().
00063 { 00064 AlgebraicSymMatrix55 errorMatrix = originalFTS.curvilinearError().matrix(); 00065 AlgebraicMatrix55 curv2perigee = jacobianCurvilinear2Perigee(originalFTS); 00066 return PerigeeTrajectoryError(ROOT::Math::Similarity(curv2perigee,errorMatrix)); 00067 }
PerigeeTrajectoryParameters PerigeeConversions::ftsToPerigeeParameters | ( | const FTS & | originalFTS, | |
const GlobalPoint & | referencePoint, | |||
double & | pt | |||
) | const |
This method calculates the perigee parameters from a given FTS and a reference point.
Definition at line 7 of file PerigeeConversions.cc.
References FreeTrajectoryState::charge(), geometryDiff::epsilon, Exception, MagneticField::inInverseGeV(), GlobalTrajectoryParameters::magneticField(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), phi, FreeTrajectoryState::position(), funct::sqrt(), theta, PV3DBase< T, PVType, FrameType >::theta(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and z.
Referenced by TrajectoryStateClosestToPoint::TrajectoryStateClosestToPoint().
00009 { 00010 GlobalVector impactDistance = originalFTS.position() - referencePoint; 00011 00012 pt = originalFTS.momentum().perp(); 00013 if (pt==0.) throw cms::Exception("PerigeeConversions", "Track with pt=0"); 00014 00015 double theta = originalFTS.momentum().theta(); 00016 double phi = originalFTS.momentum().phi(); 00017 double field = originalFTS.parameters().magneticField().inInverseGeV(originalFTS.position()).z(); 00018 00019 double positiveMomentumPhi = ( (phi>0) ? phi : (2*M_PI+phi) ); 00020 double positionPhi = impactDistance.phi(); 00021 double positivePositionPhi = 00022 ( (positionPhi>0) ? positionPhi : (2*M_PI+positionPhi) ); 00023 double phiDiff = positiveMomentumPhi - positivePositionPhi; 00024 if (phiDiff<0.0) phiDiff+= (2*M_PI); 00025 double signEpsilon = ( (phiDiff > M_PI) ? -1.0 : 1.0); 00026 00027 double epsilon = signEpsilon * 00028 sqrt ( impactDistance.x()*impactDistance.x() + 00029 impactDistance.y()*impactDistance.y() ); 00030 00031 // The track parameters: 00032 AlgebraicVector5 theTrackParameters; 00033 00034 double signTC = - originalFTS.charge(); 00035 bool isCharged = (signTC!=0); 00036 if (isCharged) { 00037 theTrackParameters[0] = field / pt*signTC; 00038 } else { 00039 theTrackParameters[0] = 1 / pt; 00040 } 00041 theTrackParameters[1] = theta; 00042 theTrackParameters[2] = phi; 00043 theTrackParameters[3] = epsilon; 00044 theTrackParameters[4] = impactDistance.z(); 00045 return PerigeeTrajectoryParameters(theTrackParameters, isCharged); 00046 }
AlgebraicMatrix55 PerigeeConversions::jacobianCurvilinear2Perigee | ( | const FreeTrajectoryState & | fts | ) | const |
Definition at line 221 of file PerigeeConversions.cc.
References funct::cos(), Vector3DBase< T, FrameTag >::cross(), Vector3DBase< T, FrameTag >::dot(), e, I, MagneticField::inInverseGeV(), K, PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), FreeTrajectoryState::momentum(), N, p, FreeTrajectoryState::parameters(), FreeTrajectoryState::position(), FreeTrajectoryState::signedInverseMomentum(), funct::tan(), PV3DBase< T, PVType, FrameType >::theta(), FreeTrajectoryState::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), V, PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by jacobianCurvilinear2Perigee_old(), and PerigeeKinematicState::PerigeeKinematicState().
00221 { 00222 00223 GlobalVector p = fts.momentum(); 00224 00225 GlobalVector Z = GlobalVector(0.,0.,1.); 00226 GlobalVector T = p.unit(); 00227 GlobalVector U = Z.cross(T).unit();; 00228 GlobalVector V = T.cross(U); 00229 00230 GlobalVector I = GlobalVector(-p.x(), -p.y(), 0.); //opposite to track dir. 00231 I = I.unit(); 00232 GlobalVector J(-I.y(), I.x(),0.); //counterclockwise rotation 00233 GlobalVector K(Z); 00234 GlobalPoint x = fts.position(); 00235 GlobalVector B = fts.parameters().magneticField().inInverseGeV(x); 00236 GlobalVector H = B.unit(); 00237 GlobalVector HxT = H.cross(T); 00238 GlobalVector N = HxT.unit(); 00239 double alpha = HxT.mag(); 00240 double qbp = fts.signedInverseMomentum(); 00241 double Q = -B.mag() * qbp; 00242 double alphaQ = alpha * Q; 00243 00244 double lambda = 0.5 * M_PI - p.theta(); 00245 double coslambda = cos(lambda), tanlambda = tan(lambda); 00246 00247 double TI = T.dot(I); 00248 double NU = N.dot(U); 00249 double NV = N.dot(V); 00250 double UI = U.dot(I); 00251 double VI = V.dot(I); 00252 double UJ = U.dot(J); 00253 double VJ = V.dot(J); 00254 double UK = U.dot(K); 00255 double VK = V.dot(K); 00256 00257 AlgebraicMatrix55 jac; 00258 00259 if( fabs(fts.transverseCurvature())<1.e-10 ) { 00260 jac(0,0) = 1/coslambda; 00261 jac(0,1) = tanlambda/coslambda/fts.parameters().momentum().mag(); 00262 }else{ 00263 double Bz = B.z(); 00264 jac(0,0) = -Bz/coslambda; 00265 jac(0,1) = -Bz * tanlambda/coslambda*qbp; 00266 jac(1,3) = alphaQ * NV * UI/TI; 00267 jac(1,4) = alphaQ * NV * VI/TI; 00268 jac(0,3) = -jac(0,1) * jac(1,3); 00269 jac(0,4) = -jac(0,1) * jac(1,4); 00270 jac(2,3) = -alphaQ/coslambda * NU * UI/TI; 00271 jac(2,4) = -alphaQ/coslambda * NU * VI/TI; 00272 } 00273 jac(1,1) = -1.; 00274 jac(2,2) = 1.; 00275 jac(3,3) = VK/TI; 00276 jac(3,4) = -UK/TI; 00277 jac(4,3) = -VJ/TI; 00278 jac(4,4) = UJ/TI; 00279 00280 return jac; 00281 00282 }
AlgebraicMatrix PerigeeConversions::jacobianCurvilinear2Perigee_old | ( | const FreeTrajectoryState & | fts | ) | const |
Jacobians of tranformations between curvilinear frame at point of closest approach in transverse plane and perigee frame.
The fts must therefore be given at exactly this point in order to yield the correct Jacobians.
Definition at line 215 of file PerigeeConversions.cc.
References asHepMatrix(), and jacobianCurvilinear2Perigee().
00215 { 00216 return asHepMatrix(jacobianCurvilinear2Perigee(fts)); 00217 }
AlgebraicMatrix66 PerigeeConversions::jacobianParameters2Cartesian | ( | const AlgebraicVector3 & | momentum, | |
const GlobalPoint & | position, | |||
const TrackCharge & | charge, | |||
const MagneticField * | field | |||
) | const |
Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian.
Definition at line 189 of file PerigeeConversions.cc.
References funct::cos(), Exception, MagneticField::inInverseGeV(), funct::sin(), funct::tan(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by jacobianParameters2Cartesian_old(), and KinematicPerigeeConversions::jacobianParameters2Kinematic().
00192 { 00193 if (momentum[0]==0.) throw cms::Exception("PerigeeConversions", "Track with rho=0"); 00194 double factor = 1.; 00195 if (charge!=0) { 00196 double bField = field->inInverseGeV(position).z(); 00197 factor = -bField*charge; 00198 } 00199 AlgebraicMatrix66 frameTransJ; 00200 frameTransJ(0,0) = 1; 00201 frameTransJ(1,1) = 1; 00202 frameTransJ(2,2) = 1; 00203 frameTransJ(3,3) = - factor * cos(momentum[2]) / (momentum[0]*momentum[0]); 00204 frameTransJ(4,3) = - factor * sin(momentum[2]) / (momentum[0]*momentum[0]); 00205 frameTransJ(5,3) = - factor / tan(momentum[1]) / (momentum[0]*momentum[0]); 00206 00207 frameTransJ(3,5) = - factor * sin(momentum[2]) / (momentum[0]); 00208 frameTransJ(4,5) = factor * cos(momentum[2]) / (momentum[0]); 00209 frameTransJ(5,4) = - factor / (momentum[0]*sin(momentum[1])*sin(momentum[1])); 00210 00211 return frameTransJ; 00212 }
AlgebraicMatrix PerigeeConversions::jacobianParameters2Cartesian_old | ( | const AlgebraicVector & | momentum, | |
const GlobalPoint & | position, | |||
const TrackCharge & | charge, | |||
const MagneticField * | field | |||
) | const |
Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian.
Definition at line 182 of file PerigeeConversions.cc.
References asHepMatrix(), and jacobianParameters2Cartesian().
00184 { 00185 return asHepMatrix(jacobianParameters2Cartesian(asSVector<3>(momentum), position, charge, field)); 00186 }
AlgebraicMatrix55 PerigeeConversions::jacobianPerigee2Curvilinear | ( | const GlobalTrajectoryParameters & | gtp | ) | const |
Definition at line 291 of file PerigeeConversions.cc.
References funct::cos(), Vector3DBase< T, FrameTag >::cross(), Vector3DBase< T, FrameTag >::dot(), e, I, MagneticField::inInverseGeV(), K, PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), N, p, GlobalTrajectoryParameters::position(), GlobalTrajectoryParameters::signedInverseMomentum(), funct::sin(), PV3DBase< T, PVType, FrameType >::theta(), GlobalTrajectoryParameters::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), V, PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by jacobianPerigee2Curvilinear_old().
00291 { 00292 00293 GlobalVector p = gtp.momentum(); 00294 00295 GlobalVector Z = GlobalVector(0.,0.,1.); 00296 GlobalVector T = p.unit(); 00297 GlobalVector U = Z.cross(T).unit();; 00298 GlobalVector V = T.cross(U); 00299 00300 GlobalVector I = GlobalVector(-p.x(), -p.y(), 0.); //opposite to track dir. 00301 I = I.unit(); 00302 GlobalVector J(-I.y(), I.x(),0.); //counterclockwise rotation 00303 GlobalVector K(Z); 00304 GlobalPoint x = gtp.position(); 00305 GlobalVector B = gtp.magneticField().inInverseGeV(x); 00306 GlobalVector H = B.unit(); 00307 GlobalVector HxT = H.cross(T); 00308 GlobalVector N = HxT.unit(); 00309 double alpha = HxT.mag(); 00310 double qbp = gtp.signedInverseMomentum(); 00311 double Q = -B.mag() * qbp; 00312 double alphaQ = alpha * Q; 00313 00314 double lambda = 0.5 * M_PI - p.theta(); 00315 double coslambda = cos(lambda), sinlambda = sin(lambda); 00316 double mqbpt = -1./coslambda * qbp; 00317 00318 double TJ = T.dot(J); 00319 double TK = T.dot(K); 00320 double NU = N.dot(U); 00321 double NV = N.dot(V); 00322 double UJ = U.dot(J); 00323 double VJ = V.dot(J); 00324 double UK = U.dot(K); 00325 double VK = V.dot(K); 00326 00327 AlgebraicMatrix55 jac; 00328 00329 if( fabs(gtp.transverseCurvature())<1.e-10 ) { 00330 jac(0,0) = coslambda; 00331 jac(0,1) = sinlambda/coslambda/gtp.momentum().mag(); 00332 }else{ 00333 jac(0,0) = -coslambda/B.z(); 00334 jac(0,1) = -sinlambda * mqbpt; 00335 jac(1,3) = -alphaQ * NV * TJ; 00336 jac(1,4) = -alphaQ * NV * TK; 00337 jac(2,3) = -alphaQ/coslambda * NU * TJ; 00338 jac(2,4) = -alphaQ/coslambda * NU * TK; 00339 } 00340 jac(1,1) = -1.; 00341 jac(2,2) = 1.; 00342 jac(3,3) = UJ; 00343 jac(3,4) = UK; 00344 jac(4,3) = VJ; 00345 jac(4,4) = VK; 00346 00347 return jac; 00348 }
AlgebraicMatrix PerigeeConversions::jacobianPerigee2Curvilinear_old | ( | const GlobalTrajectoryParameters & | gtp | ) | const |
Definition at line 286 of file PerigeeConversions.cc.
References asHepMatrix(), and jacobianPerigee2Curvilinear().
00286 { 00287 return asHepMatrix(jacobianPerigee2Curvilinear(gtp)); 00288 }
GlobalVector PerigeeConversions::momentumFromPerigee | ( | const PerigeeTrajectoryParameters & | parameters, | |
double | pt, | |||
const GlobalPoint & | referencePoint | |||
) | const |
This method returns the (Cartesian) momentum from the PerigeeTrajectoryParameters.
Definition at line 118 of file PerigeeConversions.cc.
References funct::cos(), PerigeeTrajectoryParameters::phi(), funct::sin(), funct::tan(), and PerigeeTrajectoryParameters::theta().
00119 { 00120 return GlobalVector(cos(parameters.phi()) * pt, 00121 sin(parameters.phi()) * pt, 00122 pt / tan(parameters.theta())); 00123 }
GlobalVector PerigeeConversions::momentumFromPerigee | ( | const AlgebraicVector3 & | momentum, | |
const TrackCharge & | charge, | |||
const GlobalPoint & | referencePoint, | |||
const MagneticField * | field | |||
) | const |
Definition at line 132 of file PerigeeConversions.cc.
References funct::abs(), funct::cos(), Exception, MagneticField::inInverseGeV(), funct::sin(), funct::tan(), and PV3DBase< T, PVType, FrameType >::z().
00134 { 00135 double pt; 00136 if (momentum[0]==0.) throw cms::Exception("PerigeeConversions", "Track with rho=0"); 00137 00138 if (charge!=0) { 00139 pt = std::abs(field->inInverseGeV(referencePoint).z() / momentum[0]); 00140 } else { 00141 pt = 1 / momentum[0]; 00142 } 00143 return GlobalVector(cos(momentum[2]) * pt, 00144 sin(momentum[2]) * pt, 00145 pt/tan(momentum[1])); 00146 }
GlobalVector PerigeeConversions::momentumFromPerigee | ( | const AlgebraicVector & | momentum, | |
const TrackCharge & | charge, | |||
const GlobalPoint & | referencePoint, | |||
const MagneticField * | field | |||
) | const |
This method returns the (Cartesian) momentum.
The parameters need not be the full perigee parameters, as long as the first 3 parameters are the transverse curvature, theta and phi.
Definition at line 126 of file PerigeeConversions.cc.
Referenced by TrajectoryStateClosestToPoint::calculateFTS(), TrajectoryStateClosestToPoint::momentum(), and MuonStandaloneAlgorithm::produceFreeTrajectoryState().
00127 { 00128 return momentumFromPerigee(asSVector<3>(momentum), charge, referencePoint, field); 00129 }
GlobalPoint PerigeeConversions::positionFromPerigee | ( | const PerigeeTrajectoryParameters & | parameters, | |
const GlobalPoint & | referencePoint | |||
) | const |
This method returns the position (on the helix) at which the parameters are defined.
Definition at line 108 of file PerigeeConversions.cc.
References funct::cos(), funct::sin(), PerigeeTrajectoryParameters::vector(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by TrajectoryStateClosestToPoint::calculateFTS(), TrajectoryStateClosestToPoint::position(), and MuonStandaloneAlgorithm::produceFreeTrajectoryState().
00109 { 00110 AlgebraicVector5 theVector = parameters.vector(); 00111 return GlobalPoint(theVector[3]*sin(theVector[2])+referencePoint.x(), 00112 -theVector[3]*cos(theVector[2])+referencePoint.y(), 00113 theVector[4]+referencePoint.z()); 00114 }
TrajectoryStateClosestToPoint PerigeeConversions::trajectoryStateClosestToPoint | ( | const AlgebraicVector3 & | momentum, | |
const GlobalPoint & | referencePoint, | |||
const TrackCharge & | charge, | |||
const AlgebraicSymMatrix66 & | theCovarianceMatrix, | |||
const MagneticField * | field | |||
) | const |
Definition at line 166 of file PerigeeConversions.cc.
00169 { 00170 AlgebraicMatrix66 param2cart = jacobianParameters2Cartesian 00171 (momentum, referencePoint, charge, field); 00172 CartesianTrajectoryError cartesianTrajErr(ROOT::Math::Similarity(param2cart, theCovarianceMatrix)); 00173 00174 FTS theFTS(GlobalTrajectoryParameters(referencePoint, 00175 momentumFromPerigee(momentum, charge, referencePoint, field), charge, 00176 field), cartesianTrajErr); 00177 00178 return TrajectoryStateClosestToPoint(theFTS, referencePoint); 00179 }
TrajectoryStateClosestToPoint PerigeeConversions::trajectoryStateClosestToPoint | ( | const AlgebraicVector & | momentum, | |
const GlobalPoint & | referencePoint, | |||
const TrackCharge & | charge, | |||
const AlgebraicMatrix & | theCovarianceMatrix, | |||
const MagneticField * | field | |||
) | const |
Public constructor.
This constructor takes a momentum, with parameters (transverse curvature, theta, phi) and a position, which is both the reference position and the position at which the momentum is defined. The covariance matrix is defined for these 6 parameters, in the order (x, y, z, transverse curvature, theta, phi).
field | FIXME !!! why not Sym !!?? |
Definition at line 155 of file PerigeeConversions.cc.
Referenced by PerigeeLinearizedTrackState::createRefittedTrackState(), and PerigeeMultiLTS::createRefittedTrackState().
00157 { 00158 AlgebraicSymMatrix sym; sym.assign(theCovarianceMatrix); // below, this was used for Matrix => SymMatrix 00159 return trajectoryStateClosestToPoint(asSVector<3>(momentum), referencePoint, 00160 charge, asSMatrix<6>(sym), field); 00161 00162 }