CMS 3D CMS Logo

Typedefs | Functions

PerigeeConversions Namespace Reference

Typedefs

typedef FreeTrajectoryState FTS

Functions

CurvilinearTrajectoryError curvilinearError (const PerigeeTrajectoryError &perigeeError, const GlobalTrajectoryParameters &gtp)
PerigeeTrajectoryError ftsToPerigeeError (const FTS &originalFTS)
PerigeeTrajectoryParameters ftsToPerigeeParameters (const FTS &originalFTS, const GlobalPoint &referencePoint, double &pt)
AlgebraicMatrix55 jacobianCurvilinear2Perigee (const FreeTrajectoryState &fts)
AlgebraicMatrix66 jacobianParameters2Cartesian (const AlgebraicVector3 &momentum, const GlobalPoint &position, const TrackCharge &charge, const MagneticField *field)
AlgebraicMatrix55 jacobianPerigee2Curvilinear (const GlobalTrajectoryParameters &gtp)
GlobalVector momentumFromPerigee (const AlgebraicVector3 &momentum, const TrackCharge &charge, const GlobalPoint &referencePoint, const MagneticField *field)
GlobalVector momentumFromPerigee (const PerigeeTrajectoryParameters &parameters, double pt, const GlobalPoint &referencePoint)
GlobalPoint positionFromPerigee (const PerigeeTrajectoryParameters &parameters, const GlobalPoint &referencePoint)
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint (const AlgebraicVector3 &momentum, const GlobalPoint &referencePoint, const TrackCharge &charge, const AlgebraicSymMatrix66 &theCovarianceMatrix, const MagneticField *field)

Detailed Description

namespace provides several functions to transform perigee parameters to and from various other parametrisations.


Typedef Documentation

Definition at line 16 of file PerigeeConversions.h.


Function Documentation

CurvilinearTrajectoryError PerigeeConversions::curvilinearError ( const PerigeeTrajectoryError perigeeError,
const GlobalTrajectoryParameters gtp 
)

Definition at line 61 of file PerigeeConversions.cc.

References PerigeeTrajectoryError::covarianceMatrix(), and jacobianPerigee2Curvilinear().

Referenced by TrajectoryStateClosestToPoint::calculateFTS().

{
  AlgebraicMatrix55 perigee2curv = jacobianPerigee2Curvilinear(gtp);
  return CurvilinearTrajectoryError(ROOT::Math::Similarity(perigee2curv, perigeeError.covarianceMatrix()));
}
PerigeeTrajectoryError PerigeeConversions::ftsToPerigeeError ( const FTS originalFTS)
PerigeeTrajectoryParameters PerigeeConversions::ftsToPerigeeParameters ( const FTS originalFTS,
const GlobalPoint referencePoint,
double &  pt 
)

calculates the perigee parameters from a given FTS and a reference point.

Definition at line 8 of file PerigeeConversions.cc.

References abs, FreeTrajectoryState::charge(), epsilon, Exception, MagneticField::inInverseGeV(), M_PI, GlobalTrajectoryParameters::magneticField(), FreeTrajectoryState::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), phi, FreeTrajectoryState::position(), mathSSE::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 MatcherUsingTracksAlgorithm::getChi2(), and TrajectoryStateClosestToPoint::TrajectoryStateClosestToPoint().

{
  GlobalVector impactDistance = originalFTS.position() - referencePoint;

  pt = originalFTS.momentum().perp();
  if (pt==0.) throw cms::Exception("PerigeeConversions", "Track with pt=0");
  
  double theta = originalFTS.momentum().theta();
  double phi = originalFTS.momentum().phi();
  double field  = originalFTS.parameters().magneticField().inInverseGeV(originalFTS.position()).z();
//   if (field==0.) throw cms::Exception("PerigeeConversions", "Field is 0") << " at " << originalFTS.position() << "\n" ;

  double positiveMomentumPhi = ( (phi>0) ? phi : (2*M_PI+phi) );
  double positionPhi = impactDistance.phi();
  double positivePositionPhi =
   ( (positionPhi>0) ? positionPhi : (2*M_PI+positionPhi) );
  double phiDiff = positiveMomentumPhi - positivePositionPhi;
  if (phiDiff<0.0) phiDiff+= (2*M_PI);
  double signEpsilon = ( (phiDiff > M_PI) ? -1.0 : 1.0);

  double epsilon = signEpsilon *
                     sqrt ( impactDistance.x()*impactDistance.x() +
                            impactDistance.y()*impactDistance.y() );

  // The track parameters:
  AlgebraicVector5 theTrackParameters;

  double signTC = - originalFTS.charge();
  bool isCharged = (signTC!=0) && (std::abs(field)>1.e-10);
  if (isCharged) {
    theTrackParameters[0] = field / pt*signTC;
  } else {
    theTrackParameters[0] = 1 / pt;
  }
  theTrackParameters[1] = theta;
  theTrackParameters[2] = phi;
  theTrackParameters[3] = epsilon;
  theTrackParameters[4] = impactDistance.z();
  return PerigeeTrajectoryParameters(theTrackParameters, isCharged);
}
AlgebraicMatrix55 PerigeeConversions::jacobianCurvilinear2Perigee ( const FreeTrajectoryState fts)

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 156 of file PerigeeConversions.cc.

References abs, alpha, Vector3DBase< T, FrameTag >::cross(), Vector3DBase< T, FrameTag >::dot(), alignCSCRings::e, Exhume::I, MagneticField::inInverseGeV(), M_PI, PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), FreeTrajectoryState::momentum(), N, AlCaHLTBitMon_ParallelJobs::p, FreeTrajectoryState::parameters(), FreeTrajectoryState::position(), FreeTrajectoryState::signedInverseMomentum(), PV3DBase< T, PVType, FrameType >::theta(), FreeTrajectoryState::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and Gflash::Z.

Referenced by ftsToPerigeeError(), and PerigeeKinematicState::PerigeeKinematicState().

                                                                             {

  GlobalVector p = fts.momentum();

  GlobalVector Z = GlobalVector(0.,0.,1.);
  GlobalVector T = p.unit();
  GlobalVector U = Z.cross(T).unit();; 
  GlobalVector V = T.cross(U);

  GlobalVector I = GlobalVector(-p.x(), -p.y(), 0.); //opposite to track dir.
  I = I.unit();
  GlobalVector J(-I.y(), I.x(),0.); //counterclockwise rotation
  GlobalVector K(Z);
  GlobalPoint x = fts.position();
  GlobalVector B  = fts.parameters().magneticField().inInverseGeV(x);
  GlobalVector H = B.unit();
  GlobalVector HxT = H.cross(T);
  GlobalVector N = HxT.unit();
  double alpha = HxT.mag();
  double qbp = fts.signedInverseMomentum();
  double Q = -B.mag() * qbp;
  double alphaQ = alpha * Q;

  double lambda = 0.5 * M_PI - p.theta();
  double sinlambda, coslambda;  vdt::fast_sincos(lambda, sinlambda, coslambda);
  double seclambda = 1./coslambda;

  double ITI = 1./T.dot(I);
  double NU = N.dot(U);
  double NV = N.dot(V);
  double UI = U.dot(I);
  double VI = V.dot(I);
  double UJ = U.dot(J);
  double VJ = V.dot(J);
  double UK = U.dot(K);
  double VK = V.dot(K);

  AlgebraicMatrix55 jac;

  if( fabs(fts.transverseCurvature())<1.e-10 ) {
    jac(0,0) = seclambda;
    jac(0,1) = sinlambda*seclambda*seclambda*std::abs(qbp);
  }else{
    double Bz = B.z();
    jac(0,0) = -Bz * seclambda;
    jac(0,1) = -Bz * sinlambda*seclambda*seclambda*qbp;
    jac(1,3) = alphaQ * NV * UI*ITI;
    jac(1,4) = alphaQ * NV * VI*ITI;
    jac(0,3) = -jac(0,1) * jac(1,3);
    jac(0,4) = -jac(0,1) * jac(1,4);
    jac(2,3) = -alphaQ*seclambda * NU * UI*ITI;
    jac(2,4) = -alphaQ*seclambda * NU * VI*ITI;
  }
  jac(1,1) = -1.;
  jac(2,2) = 1.;
  jac(3,3) = VK*ITI;
  jac(3,4) = -UK*ITI;
  jac(4,3) = -VJ*ITI;
  jac(4,4) = UJ*ITI;
  
  return jac;
  
}
AlgebraicMatrix66 PerigeeConversions::jacobianParameters2Cartesian ( const AlgebraicVector3 momentum,
const GlobalPoint position,
const TrackCharge charge,
const MagneticField field 
)

Jacobians of tranformations between the parametrixation (x, y, z, transverse curvature, theta, phi) to Cartesian

Definition at line 124 of file PerigeeConversions.cc.

References abs, ecalTB2006H4_GenSimDigiReco_cfg::bField, alignmentValidation::c1, DeDxDiscriminatorTools::charge(), alignCSCRings::e, f, python::connectstrParser::f1, python::connectstrParser::f2, MagneticField::inInverseGeV(), indexGen::s2, and PV3DBase< T, PVType, FrameType >::z().

Referenced by KinematicPerigeeConversions::jacobianParameters2Kinematic(), and trajectoryStateClosestToPoint().

{
  float factor = -1.;
  float bField = field->inInverseGeV(position).z();
  if (charge!=0 && std::abs(bField)>1.e-10f)
    factor = bField*charge;
 

  float s1,c1; vdt::fast_sincosf(momentum[1],s1,c1);
  float s2,c2; vdt::fast_sincosf(momentum[2],s2,c2);
  float f1 = factor/(momentum[0]*momentum[0]);
  float f2 = factor/momentum[0];

  AlgebraicMatrix66 frameTransJ;
  frameTransJ(0,0) = 1;
  frameTransJ(1,1) = 1;
  frameTransJ(2,2) = 1;
  frameTransJ(3,3) =  (f1 * c2); 
  frameTransJ(4,3) =  (f1 * s2);
  frameTransJ(5,3) =  (f1*c1/s1);  

  frameTransJ(3,5) = (f2 * s2);
  frameTransJ(4,5) = -(f2 * c2);
  frameTransJ(5,4) = (f2/(s1*s1));

  return frameTransJ;
}
AlgebraicMatrix55 PerigeeConversions::jacobianPerigee2Curvilinear ( const GlobalTrajectoryParameters gtp)

Definition at line 223 of file PerigeeConversions.cc.

References alpha, Vector3DBase< T, FrameTag >::cross(), Vector3DBase< T, FrameTag >::dot(), alignCSCRings::e, f, Exhume::I, MagneticField::inInverseGeV(), M_PI, PV3DBase< T, PVType, FrameType >::mag(), GlobalTrajectoryParameters::magneticField(), GlobalTrajectoryParameters::momentum(), N, AlCaHLTBitMon_ParallelJobs::p, GlobalTrajectoryParameters::position(), GlobalTrajectoryParameters::signedInverseMomentum(), PV3DBase< T, PVType, FrameType >::theta(), GlobalTrajectoryParameters::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and Gflash::Z.

Referenced by curvilinearError().

                                                                                     {

  GlobalVector p = gtp.momentum();

  GlobalVector Z = GlobalVector(0.f,0.f,1.f);
  GlobalVector T = p.unit();
  GlobalVector U = Z.cross(T).unit();; 
  GlobalVector V = T.cross(U);

  GlobalVector I = GlobalVector(-p.x(), -p.y(), 0.f); //opposite to track dir.
  I = I.unit();
  GlobalVector J(-I.y(), I.x(),0.f); //counterclockwise rotation
  GlobalVector K(Z);
  GlobalPoint x = gtp.position();
  GlobalVector B  = gtp.magneticField().inInverseGeV(x);
  GlobalVector H = B.unit();
  GlobalVector HxT = H.cross(T);
  GlobalVector N = HxT.unit();
  double alpha = HxT.mag();
  double qbp = gtp.signedInverseMomentum();
  double Q = -B.mag() * qbp;
  double alphaQ = alpha * Q;


  double lambda = 0.5 * M_PI - p.theta();
  double sinlambda, coslambda;  vdt::fast_sincos(lambda, sinlambda, coslambda);
  double seclambda = 1./coslambda;

  double mqbpt = -1./coslambda * qbp;

  double TJ = T.dot(J);
  double TK = T.dot(K);
  double NU = N.dot(U);
  double NV = N.dot(V);
  double UJ = U.dot(J);
  double VJ = V.dot(J);
  double UK = U.dot(K);
  double VK = V.dot(K);

  AlgebraicMatrix55 jac;

  if( fabs(gtp.transverseCurvature())<1.e-10f ) {
    jac(0,0) = coslambda;
    jac(0,1) = sinlambda/coslambda/gtp.momentum().mag();
  }else{
    jac(0,0) = -coslambda/B.z();
    jac(0,1) = -sinlambda * mqbpt;
    jac(1,3) = -alphaQ * NV * TJ;
    jac(1,4) = -alphaQ * NV * TK;
    jac(2,3) = -alphaQ*seclambda * NU * TJ;
    jac(2,4) = -alphaQ*seclambda * NU * TK;
  }
  jac(1,1) = -1.;
  jac(2,2) = 1.;
  jac(3,3) = UJ;
  jac(3,4) = UK;
  jac(4,3) = VJ;
  jac(4,4) = VK;
  
  return jac;
}
GlobalVector PerigeeConversions::momentumFromPerigee ( const AlgebraicVector3 momentum,
const TrackCharge charge,
const GlobalPoint referencePoint,
const MagneticField field 
)

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 87 of file PerigeeConversions.cc.

References abs, MagneticField::inInverseGeV(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by TrajectoryStateClosestToPoint::calculateFTS(), TrajectoryStateClosestToPoint::momentum(), and trajectoryStateClosestToPoint().

{
  double pt;

  double bz = fabs(field->inInverseGeV(referencePoint).z());
  if ( charge!=0 && std::abs(bz)>1.e-10 ) {
    pt = std::abs(bz/momentum[0]);
    // if (pt<1.e-10) throw cms::Exception("PerigeeConversions", "pt is 0");
  } else {
    pt = 1 / momentum[0];
  }

  return GlobalVector(vdt::fast_cos(momentum[2]) * pt,
                      vdt::fast_sin(momentum[2]) * pt,
                      pt/vdt::fast_tan(momentum[1]));
}
GlobalVector PerigeeConversions::momentumFromPerigee ( const PerigeeTrajectoryParameters parameters,
double  pt,
const GlobalPoint referencePoint 
)

returns the (Cartesian) momentum from the PerigeeTrajectoryParameters

Definition at line 78 of file PerigeeConversions.cc.

References PerigeeTrajectoryParameters::phi(), and PerigeeTrajectoryParameters::theta().

{
  return GlobalVector(vdt::fast_cos(parameters.phi()) * pt,
                      vdt::fast_sin(parameters.phi()) * pt,
                      pt /vdt::fast_tan(parameters.theta()));
}
GlobalPoint PerigeeConversions::positionFromPerigee ( const PerigeeTrajectoryParameters parameters,
const GlobalPoint referencePoint 
)

returns the position (on the helix) at which the parameters are defined

Definition at line 68 of file PerigeeConversions.cc.

References PerigeeTrajectoryParameters::vector(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by TrajectoryStateClosestToPoint::calculateFTS(), and TrajectoryStateClosestToPoint::position().

{
  AlgebraicVector5 theVector = parameters.vector();
  return GlobalPoint(theVector[3]*vdt::fast_sin(theVector[2])+referencePoint.x(),
                     -theVector[3]*vdt::fast_cos(theVector[2])+referencePoint.y(),
                     theVector[4]+referencePoint.z());
}
TrajectoryStateClosestToPoint PerigeeConversions::trajectoryStateClosestToPoint ( const AlgebraicVector3 momentum,
const GlobalPoint referencePoint,
const TrackCharge charge,
const AlgebraicSymMatrix66 theCovarianceMatrix,
const MagneticField field 
)

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).

Definition at line 107 of file PerigeeConversions.cc.

References jacobianParameters2Cartesian(), and momentumFromPerigee().

Referenced by PerigeeMultiLTS::createRefittedTrackState().

{
  AlgebraicMatrix66 param2cart = jacobianParameters2Cartesian
        (momentum, referencePoint, charge, field);
  CartesianTrajectoryError cartesianTrajErr(ROOT::Math::Similarity(param2cart, theCovarianceMatrix));

  FTS theFTS(GlobalTrajectoryParameters(referencePoint,
            momentumFromPerigee(momentum, charge, referencePoint, field), charge,
            field), cartesianTrajErr);

  return TrajectoryStateClosestToPoint(theFTS, referencePoint);
}