#include <JacobianCurvilinearToLocal.h>
Public Member Functions | |
const AlgebraicMatrix55 & | jacobian () const |
const AlgebraicMatrix | jacobian_old () const |
JacobianCurvilinearToLocal (const Surface &surface, const LocalTrajectoryParameters &localParameters, const MagneticField &magField) | |
Private Attributes | |
AlgebraicMatrix55 | theJacobian |
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 16 of file JacobianCurvilinearToLocal.h.
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 7 of file JacobianCurvilinearToLocal.cc.
References Vector3DBase< T, FrameTag >::dot(), h, MagneticField::inTesla(), PV3DBase< T, PVType, FrameType >::mag(), LocalTrajectoryParameters::momentum(), PV3DBase< T, PVType, FrameType >::perp(), LocalTrajectoryParameters::position(), lumiQueryAPI::q, GloballyPositioned< T >::rotation(), LocalTrajectoryParameters::signedInverseMomentum(), theJacobian, Surface::toGlobal(), interactiveExample::ui, Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), TkRotation< T >::x(), x, PV3DBase< T, PVType, FrameType >::y(), TkRotation< T >::y(), PV3DBase< T, PVType, FrameType >::z(), and TkRotation< T >::z().
: theJacobian() { // Origin: TRSCSD GlobalPoint x = surface.toGlobal(localParameters.position()); // GlobalVector h = MagneticField::inInverseGeV(x); GlobalVector h = magField.inTesla(x) * 2.99792458e-3; GlobalVector hdir = h.unit(); LocalVector tnl = localParameters.momentum().unit(); GlobalVector tn = surface.toGlobal(tnl); // 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(); GlobalVector dj(rot.x()); GlobalVector dk(rot.y()); GlobalVector di(rot.z()); // rotate coordinates because of wrong coordinate system in orca // LocalVector tvw(tnl.z(), tnl.x(), tnl.y()); 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 t1r = 1./tnl.z(); 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 q = -h.mag() * localParameters.signedInverseMomentum(); double sinz =-un.dot(hdir); double cosz = vn.dot(hdir); double ui = un.dot(di)*(q*t3r); double vi = vn.dot(di)*(q*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 |
Access to Jacobian.
Definition at line 68 of file JacobianCurvilinearToLocal.cc.
References theJacobian.
Referenced by BasicSingleTrajectoryState::createLocalErrorFromCurvilinearError(), and TwoBodyDecayTrajectoryState::propagateSingleState().
{ return theJacobian; }
const AlgebraicMatrix JacobianCurvilinearToLocal::jacobian_old | ( | ) | const |
Definition at line 71 of file JacobianCurvilinearToLocal.cc.
References asHepMatrix(), and theJacobian.
{ return asHepMatrix(theJacobian); }
Definition at line 37 of file JacobianCurvilinearToLocal.h.
Referenced by jacobian(), jacobian_old(), and JacobianCurvilinearToLocal().