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