CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
JacobianLocalToCurvilinear Class Reference

#include <JacobianLocalToCurvilinear.h>

Public Member Functions

const AlgebraicMatrix55jacobian () 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)
 

Private Attributes

AlgebraicMatrix55 theJacobian
 

Detailed Description

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 21 of file JacobianLocalToCurvilinear.h.

Constructor & Destructor Documentation

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 8 of file JacobianLocalToCurvilinear.cc.

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

11  : theJacobian(ROOT::Math::SMatrixNoInit()) {
12  GlobalPoint x = surface.toGlobal(localParameters.position());
13  GlobalVector h = magField.inInverseGeV(x);
14  GlobalVector hq = h * localParameters.signedInverseMomentum(); // changed sign
15 
16  LocalVector tnl = localParameters.direction();
17  GlobalVector tn = surface.toGlobal(tnl);
18 
19  Surface::RotationType const& rot = surface.rotation();
20 
21  compute(rot, tnl, tn, hq);
22 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:81
LocalVector direction() const
Momentum vector unit in the local frame.
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
LocalPoint position() const
Local x and y position coordinates.
float signedInverseMomentum() const
Signed inverse momentum q/p (zero for neutrals).
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:36
void compute(Surface::RotationType const &rot, LocalVector const &tnl, GlobalVector const &tn, GlobalVector const &hq)
const RotationType & rotation() const
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 24 of file JacobianLocalToCurvilinear.cc.

References compute(), LocalTrajectoryParameters::direction(), h, GlobalTrajectoryParameters::magneticFieldInInverseGeV(), makeMuonMisalignmentScenario::rot, GloballyPositioned< T >::rotation(), LocalTrajectoryParameters::signedInverseMomentum(), and Surface::toGlobal().

28  : theJacobian(ROOT::Math::SMatrixNoInit()) {
29  GlobalVector h = globalParameters.magneticFieldInInverseGeV();
30  GlobalVector hq = h * localParameters.signedInverseMomentum(); // changed sign
31 
32  LocalVector tnl = localParameters.direction();
33  GlobalVector tn = surface.toGlobal(tnl); // globalParameters.momentum().unit();
34 
35  Surface::RotationType const& rot = surface.rotation();
36 
37  compute(rot, tnl, tn, hq);
38 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:81
LocalVector direction() const
Momentum vector unit in the local frame.
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
GlobalVector magneticFieldInInverseGeV(const GlobalPoint &x) const
float signedInverseMomentum() const
Signed inverse momentum q/p (zero for neutrals).
void compute(Surface::RotationType const &rot, LocalVector const &tnl, GlobalVector const &tn, GlobalVector const &hq)
const RotationType & rotation() const

Member Function Documentation

void JacobianLocalToCurvilinear::compute ( Surface::RotationType const &  rot,
LocalVector const &  tnl,
GlobalVector const &  tn,
GlobalVector const &  hq 
)
private

Definition at line 40 of file JacobianLocalToCurvilinear.cc.

References Vector3DBase< T, FrameTag >::dot(), MillePedeFileConverter_cfg::e, mps_fire::i, 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 jacobian(), and JacobianLocalToCurvilinear().

43  {
44  // Origin: TRSDSC
45 
46  GlobalVector dj(rot.x());
47  GlobalVector dk(rot.y());
48 
49  // rotate coordinates because of wrong coordinate system in orca
50  double tvwX = tnl.z(), tvwY = tnl.x(), tvwZ = tnl.y();
51  double cosl = tn.perp();
52  if (cosl < 1.e-30)
53  cosl = 1.e-30;
54  double cosl1 = 1. / cosl;
55 
56  GlobalVector un(-tn.y() * cosl1, tn.x() * cosl1, 0.);
57  double uj = un.dot(dj);
58  double uk = un.dot(dk);
59  double sinz = -un.dot(hq);
60 
61  GlobalVector vn(-tn.z() * un.y(), tn.z() * un.x(), cosl);
62  double vj = vn.dot(dj);
63  double vk = vn.dot(dk);
64  double cosz = vn.dot(hq);
65 
66  theJacobian(0, 0) = 1.;
67  for (auto i = 1; i < 5; ++i)
68  theJacobian(0, i) = 0.;
69 
70  theJacobian(1, 0) = 0.;
71  theJacobian(2, 0) = 0.;
72 
73  theJacobian(1, 1) = tvwX * vj;
74  theJacobian(1, 2) = tvwX * vk;
75  theJacobian(2, 1) = tvwX * uj * cosl1;
76  theJacobian(2, 2) = tvwX * uk * cosl1;
77 
78  for (auto i = 0; i < 3; ++i) {
79  theJacobian(3, i) = 0.;
80  theJacobian(4, i) = 0.;
81  }
82 
83  theJacobian(3, 3) = uj;
84  theJacobian(3, 4) = uk;
85  theJacobian(4, 3) = vj;
86  theJacobian(4, 4) = vk;
87 
88  theJacobian(1, 3) = tvwY * sinz;
89  theJacobian(1, 4) = tvwZ * sinz;
90  theJacobian(2, 3) = tvwY * (cosz * cosl1);
91  theJacobian(2, 4) = tvwZ * (cosz * cosl1);
92  // end of TRSDSC
93 
94  //dbg::dbg_trace(1,"Loc2Cu", localParameters.vector(),x,dj,dk,theJacobian);
95 }
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
T z() const
Definition: PV3DBase.h:61
const AlgebraicMatrix55& JacobianLocalToCurvilinear::jacobian ( ) const
inline

Access to Jacobian.

Definition at line 42 of file JacobianLocalToCurvilinear.h.

References compute(), dso_internal, makeMuonMisalignmentScenario::rot, and theJacobian.

Referenced by BasicTrajectoryState::checkCurvilinError().

42 { return theJacobian; }

Member Data Documentation

AlgebraicMatrix55 JacobianLocalToCurvilinear::theJacobian
private

Definition at line 48 of file JacobianLocalToCurvilinear.h.

Referenced by compute(), and jacobian().