test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
JacobianLocalToCartesian Class Reference

#include <JacobianLocalToCartesian.h>

Public Member Functions

const AlgebraicMatrix65jacobian () const
 
 JacobianLocalToCartesian (const Surface &surface, const LocalTrajectoryParameters &localParameters)
 

Private Attributes

AlgebraicMatrix65 theJacobian
 

Detailed Description

Class which calculates the Jacobian matrix of the transformation from the local to the Caresian frame. 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 JacobianLocalToCartesian.h.

Constructor & Destructor Documentation

JacobianLocalToCartesian::JacobianLocalToCartesian ( const Surface surface,
const LocalTrajectoryParameters localParameters 
)

Constructor from local trajectory parameters and surface defining the local frame. NB!! No default constructor exists!

Definition at line 5 of file JacobianLocalToCartesian.cc.

References LocalTrajectoryParameters::charge(), LocalTrajectoryParameters::dxdz(), LocalTrajectoryParameters::dydz(), LocalTrajectoryParameters::pzSign(), lumiQueryAPI::q, LocalTrajectoryParameters::qbp(), dttmaxenums::R, makeMuonMisalignmentScenario::rot, GloballyPositioned< T >::rotation(), funct::sqr(), mathSSE::sqrt(), theJacobian, TkRotation< T >::xx(), TkRotation< T >::xy(), TkRotation< T >::xz(), TkRotation< T >::yx(), TkRotation< T >::yy(), TkRotation< T >::yz(), TkRotation< T >::zx(), TkRotation< T >::zy(), and TkRotation< T >::zz().

6  : theJacobian() {
7 
8  //LocalCoordinates 1 = (x, y, z, px, py, pz)
9  //LocalCoordinates 2 = (q/|p|, dx/dz, dy/dz, x, y)
10  //Transformation: x = x
11  // y = y
12  // px = (q/(q/|p|)) * (dx/dz) * sqrt(1./(1.+(dx/dz)^2+(dy/dz)^2))
13  // py = (q/(q/|p|)) * (dy/dz) * sqrt(1./(1.+(dx/dz)^2+(dy/dz)^2))
14  // pz = (q/(q/|p|)) * sqrt(1./(1.+(dx/dz)^2+(dy/dz)^2))
15  //Jacobian J((x, y, px, py, pz) = f(q/|p|, dx/dz, dy/dz, x, y))
16 
17  // AlgebraicVector5 localTrackParams = localParameters.mixedFormatVector();
18  double qbp = localParameters.qbp();
19  double dxdz = localParameters.dxdz();
20  double dydz = localParameters.dydz();
21  TrackCharge iq = localParameters.charge();
22  // for neutrals: qbp is 1/p instead of q/p -
23  // equivalent to charge 1
24  if ( iq==0 ) iq = 1;
25  double pzSign = localParameters.pzSign();
26  double q = iq*pzSign;
27  double sqr = sqrt(dxdz*dxdz + dydz*dydz + 1);
28  double den = -q/(sqr*sqr*sqr*qbp);
29 
30  // no difference between local and data member
31  AlgebraicMatrix65 & lJacobian = theJacobian;
32  lJacobian(0,3) = 1.;
33  lJacobian(1,4) = 1.;
34  lJacobian(3,0) = ( dxdz*(-q/(sqr*qbp*qbp)) );
35  lJacobian(3,1) = ( q/(sqr*qbp) + (den*dxdz*dxdz) );
36  lJacobian(3,2) = ( (den*dxdz*dydz) );
37  lJacobian(4,0) = ( dydz*(-q/(sqr*qbp*qbp)) );
38  lJacobian(4,1) = ( (den*dxdz*dydz) );
39  lJacobian(4,2) = ( q/(sqr*qbp) + (den*dydz*dydz) );
40  lJacobian(5,0) = ( -q/(sqr*qbp*qbp) );
41  lJacobian(5,1) = ( (den*dxdz) );
42  lJacobian(5,2) = ( (den*dydz) );
43 
44  /*
45  GlobalVector g1 = surface.toGlobal(LocalVector(1., 0., 0.));
46  GlobalVector g2 = surface.toGlobal(LocalVector(0., 1., 0.));
47  GlobalVector g3 = surface.toGlobal(LocalVector(0., 0., 1.));
48  */
49  AlgebraicMatrix33 Rsub;
50  /*
51  Rsub(0,0) = g1.x(); Rsub(0,1) = g2.x(); Rsub(0,2) = g3.x();
52  Rsub(1,0) = g1.y(); Rsub(1,1) = g2.y(); Rsub(1,2) = g3.y();
53  Rsub(2,0) = g1.z(); Rsub(2,1) = g2.z(); Rsub(2,2) = g3.z();
54  */
55  // need to be copied anhhow to go from float to double...
56  Surface::RotationType const & rot = surface.rotation();
57  Rsub(0,0) = rot.xx(); Rsub(0,1) = rot.yx(); Rsub(0,2) = rot.zx();
58  Rsub(1,0) = rot.xy(); Rsub(1,1) = rot.yy(); Rsub(1,2) = rot.zy();
59  Rsub(2,0) = rot.xz(); Rsub(2,1) = rot.yz(); Rsub(2,2) = rot.zz();
60 
62  R.Place_at(Rsub, 0,0);
63  R.Place_at(Rsub, 3,3);
64  theJacobian = R * lJacobian;
65  //dbg::dbg_trace(1,"Loc2Ca", localParameters.vector(),g1,g2,g3,theJacobian);
66 }
T xx() const
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepStd< double, 6, 6 > > AlgebraicMatrix66
T yx() const
T zx() const
T xy() const
int TrackCharge
Definition: TrackCharge.h:4
T zz() const
T sqrt(T t)
Definition: SSEVec.h:18
T zy() const
T yy() const
TrackCharge charge() const
Charge (-1, 0 or 1)
ROOT::Math::SMatrix< double, 6, 5, ROOT::Math::MatRepStd< double, 6, 5 > > AlgebraicMatrix65
Square< F >::type sqr(const F &f)
Definition: Square.h:13
T xz() const
const RotationType & rotation() const
float pzSign() const
Sign of the z-component of the momentum in the local frame.
T yz() const
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepStd< double, 3, 3 > > AlgebraicMatrix33

Member Function Documentation

const AlgebraicMatrix65& JacobianLocalToCartesian::jacobian ( ) const
inline

Access to Jacobian.

Definition at line 29 of file JacobianLocalToCartesian.h.

References theJacobian.

29 { return theJacobian;}

Member Data Documentation

AlgebraicMatrix65 JacobianLocalToCartesian::theJacobian
private

Definition at line 33 of file JacobianLocalToCartesian.h.

Referenced by jacobian(), and JacobianLocalToCartesian().