#include <JacobianLocalToCurvilinear.h>
Public Member Functions | |
const AlgebraicMatrix55 & | jacobian () const |
JacobianLocalToCurvilinear (const Surface &surface, const LocalTrajectoryParameters &localParameters, const MagneticField &magField) | |
JacobianLocalToCurvilinear (const Surface &surface, const LocalTrajectoryParameters &localParameters, const GlobalTrajectoryParameters &globalParameters, const MagneticField &magField) | |
Private Member Functions | |
void | compute (Surface::RotationType const &rot, LocalVector const &tnl, GlobalVector const &tn, GlobalVector const &hq) dso_internal |
Private Attributes | |
AlgebraicMatrix55 | theJacobian |
Class which calculates the Jacobian matrix of the transformation from the local to the curvilinear frame. The Jacobian is calculated during construction and thereafter cached, enabling reuse of the same Jacobian without calculating it again.
Definition at line 22 of file JacobianLocalToCurvilinear.h.
JacobianLocalToCurvilinear::JacobianLocalToCurvilinear | ( | 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 JacobianLocalToCurvilinear.cc.
References compute(), LocalTrajectoryParameters::direction(), h, MagneticField::inInverseGeV(), LocalTrajectoryParameters::position(), makeMuonMisalignmentScenario::rot, GloballyPositioned< T >::rotation(), LocalTrajectoryParameters::signedInverseMomentum(), Surface::toGlobal(), and x.
: theJacobian() { GlobalPoint x = surface.toGlobal(localParameters.position()); GlobalVector h = magField.inInverseGeV(x); GlobalVector hq = h*localParameters.signedInverseMomentum(); // changed sign LocalVector tnl = localParameters.direction(); GlobalVector tn = surface.toGlobal(tnl); // GlobalVector dj = surface.toGlobal(LocalVector(1., 0., 0.)); // GlobalVector dk = surface.toGlobal(LocalVector(0., 1., 0.)); Surface::RotationType const & rot = surface.rotation(); compute(rot, tnl, tn, hq); }
JacobianLocalToCurvilinear::JacobianLocalToCurvilinear | ( | const Surface & | surface, |
const LocalTrajectoryParameters & | localParameters, | ||
const GlobalTrajectoryParameters & | globalParameters, | ||
const MagneticField & | magField | ||
) |
Constructor from local and global trajectory parameters and surface defining the local frame.
Definition at line 30 of file JacobianLocalToCurvilinear.cc.
References compute(), LocalTrajectoryParameters::direction(), h, MagneticField::inInverseGeV(), GlobalTrajectoryParameters::position(), makeMuonMisalignmentScenario::rot, GloballyPositioned< T >::rotation(), LocalTrajectoryParameters::signedInverseMomentum(), Surface::toGlobal(), and x.
: theJacobian() { GlobalPoint x = globalParameters.position(); GlobalVector h = magField.inInverseGeV(x); GlobalVector hq = h*localParameters.signedInverseMomentum(); // changed sign LocalVector tnl = localParameters.direction(); GlobalVector tn = surface.toGlobal(tnl); // globalParameters.momentum().unit(); // GlobalVector dj = surface.toGlobal(LocalVector(1., 0., 0.)); // GlobalVector dk = surface.toGlobal(LocalVector(0., 1., 0.)); Surface::RotationType const & rot = surface.rotation(); compute(rot, tnl, tn, hq); }
void JacobianLocalToCurvilinear::compute | ( | Surface::RotationType const & | rot, |
LocalVector const & | tnl, | ||
GlobalVector const & | tn, | ||
GlobalVector const & | hq | ||
) | [private] |
Definition at line 51 of file JacobianLocalToCurvilinear.cc.
References Vector3DBase< T, FrameTag >::dot(), alignCSCRings::e, PV3DBase< T, PVType, FrameType >::perp(), theJacobian, PV3DBase< T, PVType, FrameType >::x(), TkRotation< T >::x(), PV3DBase< T, PVType, FrameType >::y(), TkRotation< T >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by JacobianLocalToCurvilinear().
{ // Origin: TRSDSC GlobalVector dj(rot.x()); GlobalVector dk(rot.y()); // GlobalVector p = surface.toGlobal(localParameters.momentum()); // GlobalVector pt(p.x(), p.y(), 0.); // pt = pt.unit(); // GlobalVector tn = p.unit(); // GlobalVector di = tsos.surface().toGlobal(LocalVector(0., 0., 1.)); // 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.); double uj = un.dot(dj); double uk = un.dot(dk); double sinz =-un.dot(hq); GlobalVector vn(-tn.z()*un.y(), tn.z()*un.x(), cosl); double vj = vn.dot(dj); double vk = vn.dot(dk); double cosz = vn.dot(hq); theJacobian(0,0) = 1.; theJacobian(1,1) = tvw.x()*vj; theJacobian(1,2) = tvw.x()*vk; theJacobian(2,1) = tvw.x()*uj*cosl1; theJacobian(2,2) = tvw.x()*uk*cosl1; theJacobian(3,3) = uj; theJacobian(3,4) = uk; theJacobian(4,3) = vj; theJacobian(4,4) = vk; theJacobian(1,3) = tvw.y()*sinz; theJacobian(1,4) = tvw.z()*sinz; theJacobian(2,3) = tvw.y()*(cosz*cosl1); theJacobian(2,4) = tvw.z()*(cosz*cosl1); // end of TRSDSC //dbg::dbg_trace(1,"Loc2Cu", localParameters.vector(),x,dj,dk,theJacobian); }
const AlgebraicMatrix55& JacobianLocalToCurvilinear::jacobian | ( | ) | const [inline] |
Access to Jacobian.
Definition at line 45 of file JacobianLocalToCurvilinear.h.
References theJacobian.
Referenced by BasicTrajectoryState::checkCurvilinError().
{return theJacobian;}
Definition at line 51 of file JacobianLocalToCurvilinear.h.
Referenced by compute(), and jacobian().