CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

JacobianCurvilinearToLocal Class Reference

#include <JacobianCurvilinearToLocal.h>

List of all members.

Public Member Functions

const AlgebraicMatrix55jacobian () const
 JacobianCurvilinearToLocal (const Surface &surface, const LocalTrajectoryParameters &localParameters, const MagneticField &magField)
 JacobianCurvilinearToLocal (const Surface &surface, const LocalTrajectoryParameters &localParameters, const GlobalTrajectoryParameters &globalParameters, const MagneticField &magField)

Private Member Functions

void compute (Surface::RotationType const &rot, GlobalVector const &tn, GlobalVector const &qh, double lz) dso_internal

Private Attributes

AlgebraicMatrix55 theJacobian

Detailed Description

Class which calculates the Jacobian matrix of the transformation from the curvilinear to the local frame. The Jacobian is calculated during construction and thereafter cached, enabling reuse of the same Jacobian without calculating it again.

Definition at line 19 of file JacobianCurvilinearToLocal.h.


Constructor & Destructor Documentation

JacobianCurvilinearToLocal::JacobianCurvilinearToLocal ( const Surface surface,
const LocalTrajectoryParameters localParameters,
const MagneticField magField 
)

Constructor from local trajectory parameters and surface defining the local frame. NB!! No default constructor exists!

Definition at line 9 of file JacobianCurvilinearToLocal.cc.

References compute(), LocalTrajectoryParameters::direction(), h, MagneticField::inInverseGeV(), LocalTrajectoryParameters::position(), makeMuonMisalignmentScenario::rot, GloballyPositioned< T >::rotation(), LocalTrajectoryParameters::signedInverseMomentum(), Surface::toGlobal(), x, and PV3DBase< T, PVType, FrameType >::z().

                                                          : theJacobian() {
 
  GlobalPoint  x = surface.toGlobal(localParameters.position());
  GlobalVector h  = magField.inInverseGeV(x);
  GlobalVector qh = h*localParameters.signedInverseMomentum();  // changed sign


  //  GlobalVector  hdir =  h.unit();
  //double q = -h.mag() * localParameters.signedInverseMomentum();
  
  LocalVector tnl = localParameters.direction();
  GlobalVector tn = surface.toGlobal(tnl);
  double t1r = 1./tnl.z();

  // GlobalVector dj = surface.toGlobal(LocalVector(1., 0., 0.));
  // GlobalVector dk = surface.toGlobal(LocalVector(0., 1., 0.));
  //  GlobalVector di = surface.toGlobal(LocalVector(0., 0., 1.));
  Surface::RotationType const & rot = surface.rotation();

  compute(rot, tn, qh, t1r);
}
JacobianCurvilinearToLocal::JacobianCurvilinearToLocal ( const Surface surface,
const LocalTrajectoryParameters localParameters,
const GlobalTrajectoryParameters globalParameters,
const MagneticField magField 
)

Definition at line 34 of file JacobianCurvilinearToLocal.cc.

References compute(), LocalTrajectoryParameters::direction(), h, MagneticField::inInverseGeV(), GlobalTrajectoryParameters::position(), makeMuonMisalignmentScenario::rot, GloballyPositioned< T >::rotation(), GlobalTrajectoryParameters::signedInverseMomentum(), Surface::toGlobal(), x, and PV3DBase< T, PVType, FrameType >::z().

                                                          : theJacobian() {
 
  GlobalPoint  x =  globalParameters.position();
  GlobalVector h  = magField.inInverseGeV(x);
  GlobalVector qh = h*globalParameters.signedInverseMomentum();  // changed sign

  //GlobalVector  hdir =  h.unit();
  //double q = -h.mag() * localParameters.signedInverseMomentum();

 
  //  GlobalVector tn = globalParameters.momentum().unit();
  //  LocalVector tnl = localParameters.momentum().unit();

  LocalVector tnl = localParameters.direction();
  GlobalVector tn = surface.toGlobal(tnl);
  double t1r = 1./tnl.z();
 
 
  Surface::RotationType const & rot = surface.rotation();

  compute(rot, tn, qh, t1r);
}

Member Function Documentation

void JacobianCurvilinearToLocal::compute ( Surface::RotationType const &  rot,
GlobalVector const &  tn,
GlobalVector const &  qh,
double  lz 
) [private]

Definition at line 61 of file JacobianCurvilinearToLocal.cc.

References Vector3DBase< T, FrameTag >::dot(), alignCSCRings::e, PV3DBase< T, PVType, FrameType >::perp(), theJacobian, interactiveExample::ui, PV3DBase< T, PVType, FrameType >::x(), TkRotation< T >::x(), PV3DBase< T, PVType, FrameType >::y(), TkRotation< T >::y(), PV3DBase< T, PVType, FrameType >::z(), and TkRotation< T >::z().

Referenced by JacobianCurvilinearToLocal().

                                                                                                                                       {
  // Origin: TRSCSD

  GlobalVector dj(rot.x());
  GlobalVector dk(rot.y());
  GlobalVector di(rot.z());

  double cosl = tn.perp(); if (cosl < 1.e-30) cosl = 1.e-30;
  double cosl1 = 1./cosl;
  GlobalVector un(-tn.y()*cosl1, tn.x()*cosl1, 0.);
  GlobalVector vn(-tn.z()*un.y(), tn.z()*un.x(), cosl);

  double uj = un.dot(dj);
  double uk = un.dot(dk);
  double vj = vn.dot(dj);
  double vk = vn.dot(dk);

  //  double t1r = 1./tvw.x();
  double t2r = t1r*t1r;
  double t3r = t1r*t2r;

  theJacobian(0,0) = 1.;
  theJacobian(1,1) = -uk*t2r;
  theJacobian(1,2) = vk*(cosl*t2r);
  theJacobian(2,1) = uj*t2r;
  theJacobian(2,2) = -vj*(cosl*t2r);
  theJacobian(3,3) = vk*t1r;
  theJacobian(3,4) = -uk*t1r;
  theJacobian(4,3) = -vj*t1r;
  theJacobian(4,4) = uj*t1r;


  double sinz = un.dot(qh);
  double cosz =-vn.dot(qh);
  double ui = un.dot(di)*(t3r);
  double vi = vn.dot(di)*(t3r);
  theJacobian(1,3) =-ui*(vk*cosz-uk*sinz);
  theJacobian(1,4) =-vi*(vk*cosz-uk*sinz);
  theJacobian(2,3) = ui*(vj*cosz-uj*sinz);
  theJacobian(2,4) = vi*(vj*cosz-uj*sinz);
  // end of TRSCSD
  //dbg::dbg_trace(1,"Cu2L", localParameters.vector(),di,dj,dk,theJacobian);
}
const AlgebraicMatrix55& JacobianCurvilinearToLocal::jacobian ( ) const [inline]

Member Data Documentation

Definition at line 46 of file JacobianCurvilinearToLocal.h.

Referenced by compute(), and jacobian().