ATLAS Offline Software
Loading...
Searching...
No Matches
Trk::JacobianLocalAnglesPhiTheta Class Reference

This jacobian transforms the local direction representation as given in the Trk::LocalDirection class into the corresponding polar and azimutal angle \( \phi, \theta \), therefor it is of dimension (2x2). More...

#include <JacobianLocalAnglesPhiTheta.h>

Inheritance diagram for Trk::JacobianLocalAnglesPhiTheta:
Collaboration diagram for Trk::JacobianLocalAnglesPhiTheta:

Public Member Functions

 JacobianLocalAnglesPhiTheta (const double angleXZ, const double angleYZ, const Amg::RotationMatrix3D &rot)
 ~JacobianLocalAnglesPhiTheta ()

Detailed Description

This jacobian transforms the local direction representation as given in the Trk::LocalDirection class into the corresponding polar and azimutal angle \( \phi, \theta \), therefor it is of dimension (2x2).

The Rotation given for construction is the one declared by the Surface frame:

Surface::transform()::getRotation()

Author
Andre.nosp@m.as.S.nosp@m.alzbu.nosp@m.rger.nosp@m.@cern.nosp@m..ch

Definition at line 31 of file JacobianLocalAnglesPhiTheta.h.

Constructor & Destructor Documentation

◆ JacobianLocalAnglesPhiTheta()

Trk::JacobianLocalAnglesPhiTheta::JacobianLocalAnglesPhiTheta ( const double angleXZ,
const double angleYZ,
const Amg::RotationMatrix3D & rot )

Definition at line 12 of file JacobianLocalAnglesPhiTheta.cxx.

12 :
13 AmgMatrix(2,2)()
14{
15
16 this->setIdentity();
17
18// Rot is a 3x3 rotation matrix from local 3-vector to global 3-vector
19
20 double cx = std::cos(angleXZ);
21 double sx = std::sin(angleXZ);
22 double cy = std::cos(angleYZ);
23 double sy = std::sin(angleYZ);
24 double norm = std::sqrt(cx*sy*cx*sy + cy*sx*cy*sx + sx*sy*sx*sy);
25 double zinv = norm/(sx*sy) ;
26//
27// local direction x = cx * sy / norm
28// y = cy * sx / norm
29// z = sx * sy / norm
30//
31
32// Determine Global theta and phi
33
34 Amg::Vector3D locdir(cx*sy/norm,cy*sx/norm,sx*sy/norm);
35 Amg::Vector3D gdir = rot*locdir;
36 double phi = gdir.phi();
37 double theta = gdir.theta();
38// std::cout << " phi " << phi << " theta " << theta << std::endl;
39 double cp = std::cos(phi);
40 double sp = std::sin(phi);
41 double ct = std::cos(theta);
42 double st = std::sin(theta);
43 AmgMatrix(3,2) J2 = AmgMatrix(3,2)::Zero(); // row, column
44//
45// derivatives dx/dphi dx/dtheta
46// dy/dphi dy/dtheta
47// dz/dphi dz/dtheta
48//
49// x = std::sin(theta)*std::cos(phi)
50// y = std::sin(theta)*std::sin(phi)
51// z = std::cos(theta)
52//
53
54 J2(0,0) = -st*sp;
55 J2(1,0) = st*cp;
56 J2(0,1) = ct*cp;
57 J2(1,1) = ct*sp;
58 J2(2,1) = -st;
59
60 AmgMatrix(2,3) J1 = AmgMatrix(2,3)::Zero();
61//
62// derivatives d alphaXZ/dx d alphaXZ/dy d alphaXZ/dz
63// d alphaYZ/dx d alphaYZ/dy d alphaYZ/dz
64//
65// tan(alphaXZ) = z / x
66// tan(alhpaYZ) = z / y
67//
68
69 J1(0,0) = -sx*sx*zinv;
70 J1(1,1) = -sy*sy*zinv;
71 J1(0,2) = cx*sx*zinv;
72 J1(1,2) = cy*sy*zinv;
73
74 const AmgMatrix(3,3)& Rot = rot;
75 AmgMatrix(2,2) Jinv = AmgMatrix(2,2)::Zero();
76//
77// Jacobian Inverse: lefthand local righthand global
78//
79 AmgMatrix(3,3) RI = Rot.inverse();
80 Jinv = J1*(RI*J2);
81//
82// Default Matrix
83//
84 (*this)(0,0) = 1.;
85 (*this)(0,1) = 0.;
86 (*this)(1,0) = 0.;
87 (*this)(1,1) = 1.;
88//
89// Store Jacobian J: right hand local
90//
91// if (succes==0) {
92// AmgMatrix(2,2) J = Jinv.inverse(succes);
93// if (succes==0) {
95// (*this)(0,0) = J(0,0);
96// (*this)(0,1) = J(0,1);
97// (*this)(1,0) = J(1,0);
98// (*this)(1,1) = J(1,1);
99// }
100// }
101}
#define AmgMatrix(rows, cols)
static Double_t sp
Eigen::Matrix< double, 3, 1 > Vector3D
@ theta
Definition ParamDefs.h:66
@ phi
Definition ParamDefs.h:75

◆ ~JacobianLocalAnglesPhiTheta()

Trk::JacobianLocalAnglesPhiTheta::~JacobianLocalAnglesPhiTheta ( )
inline

Definition at line 34 of file JacobianLocalAnglesPhiTheta.h.

34{}

The documentation for this class was generated from the following files: