#include <TrackingTools/AnalyticalJacobians/interface/JacobianCurvilinearToLocal.h>
Public Member Functions | |
const AlgebraicMatrix55 & | jacobian () const |
Access to Jacobian. | |
const AlgebraicMatrix | jacobian_old () const |
JacobianCurvilinearToLocal (const Surface &surface, const LocalTrajectoryParameters &localParameters, const MagneticField &magField) | |
Constructor from local trajectory parameters and surface defining the local frame. | |
Private Attributes | |
AlgebraicMatrix55 | theJacobian |
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(), e, h, MagneticField::inTesla(), PV3DBase< T, PVType, FrameType >::mag(), LocalTrajectoryParameters::momentum(), PV3DBase< T, PVType, FrameType >::perp(), LocalTrajectoryParameters::position(), LocalTrajectoryParameters::signedInverseMomentum(), theJacobian, Surface::toGlobal(), GloballyPositioned< T >::toLocal(), interactiveExample::ui, Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
00009 : theJacobian() { 00010 00011 // Origin: TRSCSD 00012 GlobalVector tn = surface.toGlobal(localParameters.momentum()).unit(); 00013 GlobalVector dj = surface.toGlobal(LocalVector(1., 0., 0.)); 00014 GlobalVector dk = surface.toGlobal(LocalVector(0., 1., 0.)); 00015 GlobalVector di = surface.toGlobal(LocalVector(0., 0., 1.)); 00016 LocalVector tnl = surface.toLocal(tn); 00017 // rotate coordinates because of wrong coordinate system in orca 00018 LocalVector tvw(tnl.z(), tnl.x(), tnl.y()); 00019 double cosl = tn.perp(); if (cosl < 1.e-30) cosl = 1.e-30; 00020 double cosl1 = 1./cosl; 00021 GlobalVector un(-tn.y()*cosl1, tn.x()*cosl1, 0.); 00022 GlobalVector vn(-tn.z()*un.y(), tn.z()*un.x(), cosl); 00023 double uj = un.dot(dj); 00024 double uk = un.dot(dk); 00025 double vj = vn.dot(dj); 00026 double vk = vn.dot(dk); 00027 double t1r = 1./tvw.x(); 00028 double t2r = t1r*t1r; 00029 double t3r = t1r*t2r; 00030 theJacobian(0,0) = 1.; 00031 theJacobian(1,1) = -uk*t2r; 00032 theJacobian(1,2) = vk*cosl*t2r; 00033 theJacobian(2,1) = uj*t2r; 00034 theJacobian(2,2) = -vj*cosl*t2r; 00035 theJacobian(3,3) = vk*t1r; 00036 theJacobian(3,4) = -uk*t1r; 00037 theJacobian(4,3) = -vj*t1r; 00038 theJacobian(4,4) = uj*t1r; 00039 GlobalPoint x = surface.toGlobal(localParameters.position()); 00040 // GlobalVector h = MagneticField::inInverseGeV(x); 00041 GlobalVector h = magField.inTesla(x) * 2.99792458e-3; 00042 double q = -h.mag() * localParameters.signedInverseMomentum(); 00043 double sinz =-un.dot(h.unit()); 00044 double cosz = vn.dot(h.unit()); 00045 double ui = un.dot(di); 00046 double vi = vn.dot(di); 00047 theJacobian(1,3) =-q*ui*(vk*cosz-uk*sinz)*t3r; 00048 theJacobian(1,4) =-q*vi*(vk*cosz-uk*sinz)*t3r; 00049 theJacobian(2,3) = q*ui*(vj*cosz-uj*sinz)*t3r; 00050 theJacobian(2,4) = q*vi*(vj*cosz-uj*sinz)*t3r; 00051 // end of TRSCSD 00052 //dbg::dbg_trace(1,"Cu2L", localParameters.vector(),di,dj,dk,theJacobian); 00053 }
const AlgebraicMatrix55 & JacobianCurvilinearToLocal::jacobian | ( | ) | const |
Access to Jacobian.
Definition at line 54 of file JacobianCurvilinearToLocal.cc.
References theJacobian.
Referenced by BasicSingleTrajectoryState::createLocalErrorFromCurvilinearError().
00054 { 00055 return theJacobian; 00056 }
const AlgebraicMatrix JacobianCurvilinearToLocal::jacobian_old | ( | ) | const |
Definition at line 57 of file JacobianCurvilinearToLocal.cc.
References asHepMatrix(), and theJacobian.
00057 { 00058 return asHepMatrix(theJacobian); 00059 }
Definition at line 37 of file JacobianCurvilinearToLocal.h.
Referenced by jacobian(), jacobian_old(), and JacobianCurvilinearToLocal().