CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

TkRotation< T > Class Template Reference

#include <newTkRotation.h>

List of all members.

Public Types

typedef Basic3DVector< TBasicVector
typedef Vector3DBase< T,
GlobalTag
GlobalVector
typedef Vector3DBase< T,
GlobalTag
GlobalVector

Public Member Functions

Basic3DVector< TmultiplyInverse (const Basic3DVector< T > &v) const
template<class Scalar >
Basic3DVector< Scalar > multiplyInverse (const Basic3DVector< Scalar > &v) const
Basic3DVector< TmultiplyInverse (const Basic2DVector< T > &v) const
TkRotation multiplyInverse (const TkRotation &b) const
Basic3DVector< TmultiplyInverse (const Basic3DVector< T > &v) const
template<class Scalar >
Basic3DVector< Scalar > multiplyInverse (const Basic3DVector< Scalar > &v) const
Basic3DVector< TmultiplyInverse (const Basic2DVector< T > &v) const
TkRotation multiplyInverse (const TkRotation &b) const
Basic3DVector< Toperator* (const Basic2DVector< T > &v) const
TkRotation operator* (const TkRotation &b) const
Basic3DVector< Toperator* (const Basic3DVector< T > &v) const
Basic3DVector< Toperator* (const Basic2DVector< T > &v) const
Basic3DVector< Toperator* (const Basic3DVector< T > &v) const
TkRotation operator* (const TkRotation &b) const
TkRotationoperator*= (const TkRotation &b)
TkRotationoperator*= (const TkRotation &b)
Basic3DVector< Trotate (const Basic3DVector< T > &v) const
TkRotationrotateAxes (const Basic3DVector< T > &newX, const Basic3DVector< T > &newY, const Basic3DVector< T > &newZ)
TkRotationrotateAxes (const Basic3DVector< T > &newX, const Basic3DVector< T > &newY, const Basic3DVector< T > &newZ)
Basic3DVector< TrotateBack (const Basic3DVector< T > &v) const
 TkRotation (mathSSE::Rot3< T > const &irot)
 TkRotation (const BasicVector &uX, const BasicVector &uY, const BasicVector &uZ)
 TkRotation (const T *p)
 TkRotation (const Basic3DVector< T > &axis, T phi)
template<typename U >
 TkRotation (const TkRotation< U > &a)
 TkRotation ()
 TkRotation (T xx, T xy, T xz, T yx, T yy, T yz, T zx, T zy, T zz)
 TkRotation (const T *p)
 TkRotation (const GlobalVector &aX, const GlobalVector &aY)
template<typename U >
 TkRotation (const TkRotation< U > &a)
 TkRotation (const GlobalVector &aX, const GlobalVector &aY, const GlobalVector &aZ)
 TkRotation (const Basic3DVector< T > &axis, T phi)
 TkRotation (const GlobalVector &aX, const GlobalVector &aY)
 TkRotation ()
 TkRotation (const BasicVector &aX, const BasicVector &aY)
 TkRotation (T xx, T xy, T xz, T yx, T yy, T yz, T zx, T zy, T zz)
 TkRotation (const GlobalVector &uX, const GlobalVector &uY, const GlobalVector &uZ)
TkRotationtransform (const TkRotation &b)
TkRotationtransform (const TkRotation &b)
TkRotation transposed () const
TkRotation transposed () const
Basic3DVector< Tx () const
Basic3DVector< Tx () const
T xx () const
T const & xx () const
T xy () const
T const & xy () const
T xz () const
T const & xz () const
Basic3DVector< Ty () const
Basic3DVector< Ty () const
T yx () const
T const & yx () const
T yy () const
T const & yy () const
T const & yz () const
T yz () const
Basic3DVector< Tz () const
Basic3DVector< Tz () const
T const & zx () const
T zx () const
T const & zy () const
T zy () const
T const & zz () const
T zz () const

Private Attributes

T R11
T R12
T R13
T R21
T R22
T R23
T R31
T R32
T R33
mathSSE::Rot3< Trot

Detailed Description

template<class T>
class TkRotation< T >

Rotaion matrix used by Surface.

Definition at line 32 of file newTkRotation.h.


Member Typedef Documentation

template<class T>
typedef Basic3DVector<T> TkRotation< T >::BasicVector

Definition at line 36 of file newTkRotation.h.

template<class T>
typedef Vector3DBase< T, GlobalTag> TkRotation< T >::GlobalVector

Definition at line 31 of file oldTkRotation.h.

template<class T>
typedef Vector3DBase< T, GlobalTag> TkRotation< T >::GlobalVector

Definition at line 35 of file newTkRotation.h.


Constructor & Destructor Documentation

template<class T>
TkRotation< T >::TkRotation ( ) [inline]
template<class T>
TkRotation< T >::TkRotation ( mathSSE::Rot3< T > const &  irot) [inline]

Definition at line 39 of file newTkRotation.h.

: rot(irot){}
template<class T>
TkRotation< T >::TkRotation ( T  xx,
T  xy,
T  xz,
T  yx,
T  yy,
T  yz,
T  zx,
T  zy,
T  zz 
) [inline]

Definition at line 41 of file newTkRotation.h.

                                                                    :
    rot(xx,xy,xz, yx,yy,yz, zx, zy,zz){}
template<class T>
TkRotation< T >::TkRotation ( const T p) [inline]

Definition at line 44 of file newTkRotation.h.

                          : 
    rot(p[0],p[1],p[2], 
        p[3],p[4],p[5],
        p[6],p[7],p[8]) {}
template<class T>
TkRotation< T >::TkRotation ( const GlobalVector aX,
const GlobalVector aY 
) [inline]

Definition at line 49 of file newTkRotation.h.

                                                                 {
    
    GlobalVector uX = aX.unit();
    GlobalVector uY = aY.unit();
    GlobalVector uZ(uX.cross(uY));
    
    rot.axis[0]= uX.basicVector().v;
    rot.axis[1]= uY.basicVector().v;
    rot.axis[2]= uZ.basicVector().v;
    
  }
template<class T>
TkRotation< T >::TkRotation ( const BasicVector aX,
const BasicVector aY 
) [inline]

Definition at line 61 of file newTkRotation.h.

                                                               {
    
    BasicVector uX = aX.unit();
    BasicVector uY = aY.unit();
    BasicVector uZ(uX.cross(uY));
    
    rot.axis[0]= uX.v;
    rot.axis[1]= uY.v;
    rot.axis[2]= uZ.v;
    
  }
template<class T>
TkRotation< T >::TkRotation ( const GlobalVector uX,
const GlobalVector uY,
const GlobalVector uZ 
) [inline]

Construct from global vectors of the x, y and z axes. The axes are assumed to be unit vectors forming a right-handed orthonormal basis. No checks are performed!

Definition at line 78 of file newTkRotation.h.

                                       {
    rot.axis[0]= uX.basicVector().v;
    rot.axis[1]= uY.basicVector().v;
    rot.axis[2]= uZ.basicVector().v;
  }
template<class T>
TkRotation< T >::TkRotation ( const BasicVector uX,
const BasicVector uY,
const BasicVector uZ 
) [inline]

Definition at line 85 of file newTkRotation.h.

                                      {
    rot.axis[0]= uX.v;
    rot.axis[1]= uY.v;
    rot.axis[2]= uZ.v;
  }
template<class T>
TkRotation< T >::TkRotation ( const Basic3DVector< T > &  axis,
T  phi 
) [inline]

rotation around abritrary axis by the amount of phi: its constructed by O^-1(z<->axis) rot_z(phi) O(z<->axis) the frame is rotated such that the z-asis corresponds to the rotation axis desired. THen it's rotated round the "new" z-axis, and then the initial transformation is "taken back" again. unfortuately I'm too stupid to describe such thing directly by 3 Euler angles.. hence I have to construckt it this way...by brute force

Definition at line 101 of file newTkRotation.h.

                                                   :
    rot(  cos(phi), sin(phi), 0, 
          -sin(phi), cos(phi), 0,
          0,        0, 1) {
    
    //rotation around the z-axis by  phi
    TkRotation tmpRotz ( cos(axis.phi()), sin(axis.phi()), 0.,
                         -sin(axis.phi()), cos(axis.phi()), 0.,
                         0.,              0.,              1. );
    //rotation around y-axis by theta
    TkRotation tmpRoty ( cos(axis.theta()), 0.,-sin(axis.theta()),
                         0.,              1.,              0.,
                         sin(axis.theta()), 0., cos(axis.theta()) );
    (*this)*=tmpRoty;
    (*this)*=tmpRotz;      // =  this * tmpRoty * tmpRotz 
    
    // (tmpRoty * tmpRotz)^-1 * this * tmpRoty * tmpRotz
    
    *this = (tmpRoty*tmpRotz).multiplyInverse(*this);
    
  }
template<class T>
template<typename U >
TkRotation< T >::TkRotation ( const TkRotation< U > &  a) [inline]

Definition at line 156 of file newTkRotation.h.

                                      : 
    rot (a.xx(), a.xy(), a.xz(), 
         a.yx(), a.yy(), a.yz(),
         a.zx(), a.zy(), a.zz()) {}
template<class T>
TkRotation< T >::TkRotation ( ) [inline]

Definition at line 33 of file oldTkRotation.h.

               : 
    R11( 1), R12( 0), R13( 0), 
    R21( 0), R22( 1), R23( 0),
    R31( 0), R32( 0), R33( 1) {}
template<class T>
TkRotation< T >::TkRotation ( T  xx,
T  xy,
T  xz,
T  yx,
T  yy,
T  yz,
T  zx,
T  zy,
T  zz 
) [inline]

Definition at line 38 of file oldTkRotation.h.

                                                                    :
    R11(xx), R12(xy), R13(xz), 
    R21(yx), R22(yy), R23(yz),
    R31(zx), R32(zy), R33(zz) {}
template<class T>
TkRotation< T >::TkRotation ( const T p) [inline]

Definition at line 43 of file oldTkRotation.h.

                          :
    R11(p[0]), R12(p[1]), R13(p[2]), 
    R21(p[3]), R22(p[4]), R23(p[5]),
    R31(p[6]), R32(p[7]), R33(p[8]) {}
template<class T>
TkRotation< T >::TkRotation ( const GlobalVector aX,
const GlobalVector aY 
) [inline]

Definition at line 48 of file oldTkRotation.h.

                                                                 {

    GlobalVector uX = aX.unit();
    GlobalVector uY = aY.unit();
    GlobalVector uZ(uX.cross(uY));
 
    R11 = uX.x(); R12 = uX.y();  R13 = uX.z();
    R21 = uY.x(); R22 = uY.y();  R23 = uY.z();
    R31 = uZ.x(); R32 = uZ.y();  R33 = uZ.z();

  }
template<class T>
TkRotation< T >::TkRotation ( const GlobalVector aX,
const GlobalVector aY,
const GlobalVector aZ 
) [inline]

Construct from global vectors of the x, y and z axes. The axes are assumed to be unit vectors forming a right-handed orthonormal basis. No checks are performed!

Definition at line 64 of file oldTkRotation.h.

                                       :
    R11( aX.x()), R12( aX.y()), R13( aX.z()), 
    R21( aY.x()), R22( aY.y()), R23( aY.z()),
    R31( aZ.x()), R32( aZ.y()), R33( aZ.z()) {}
template<class T>
TkRotation< T >::TkRotation ( const Basic3DVector< T > &  axis,
T  phi 
) [inline]

rotation around abritrary axis by the amount of phi: its constructed by O^-1(z<->axis) rot_z(phi) O(z<->axis) the frame is rotated such that the z-asis corresponds to the rotation axis desired. THen it's rotated round the "new" z-axis, and then the initial transformation is "taken back" again. unfortuately I'm too stupid to describe such thing directly by 3 Euler angles.. hence I have to construckt it this way...by brute force

Definition at line 79 of file oldTkRotation.h.

                                                   :
    R11( cos(phi) ), R12( sin(phi)), R13( 0), 
    R21( -sin(phi)), R22( cos(phi)), R23( 0),
    R31( 0), R32( 0), R33( 1) {

    //rotation around the z-axis by  phi
    TkRotation tmpRotz ( cos(axis.phi()), sin(axis.phi()), 0.,
                        -sin(axis.phi()), cos(axis.phi()), 0.,
                         0.,              0.,              1. );
    //rotation around y-axis by theta
    TkRotation tmpRoty ( cos(axis.theta()), 0.,-sin(axis.theta()),
                         0.,              1.,              0.,
                         sin(axis.theta()), 0., cos(axis.theta()) );
    (*this)*=tmpRoty;
    (*this)*=tmpRotz;      // =  this * tmpRoty * tmpRotz 
   
    // (tmpRoty * tmpRotz)^-1 * this * tmpRoty * tmpRotz

    *this = (tmpRoty*tmpRotz).multiplyInverse(*this);

  }
template<class T>
template<typename U >
TkRotation< T >::TkRotation ( const TkRotation< U > &  a) [inline]

Definition at line 134 of file oldTkRotation.h.

                                      : 
    R11(a.xx()), R12(a.xy()), R13(a.xz()), 
    R21(a.yx()), R22(a.yy()), R23(a.yz()),
    R31(a.zx()), R32(a.zy()), R33(a.zz()) {}

Member Function Documentation

template<class T>
Basic3DVector<T> TkRotation< T >::multiplyInverse ( const Basic3DVector< T > &  v) const [inline]
template<class T>
Basic3DVector<T> TkRotation< T >::multiplyInverse ( const Basic2DVector< T > &  v) const [inline]

Definition at line 194 of file newTkRotation.h.

                                                                     {
    return Basic3DVector<T>( xx()*v.x() + yx()*v.y(),
                             xy()*v.x() + yy()*v.y(),
                             xz()*v.x() + yz()*v.y());
  }
template<class T>
TkRotation TkRotation< T >::multiplyInverse ( const TkRotation< T > &  b) const [inline]

Definition at line 191 of file oldTkRotation.h.

                                                         {
    return TkRotation(R11*b.R11 + R21*b.R21 + R31*b.R31,
                      R11*b.R12 + R21*b.R22 + R31*b.R32,
                      R11*b.R13 + R21*b.R23 + R31*b.R33,
                      R12*b.R11 + R22*b.R21 + R32*b.R31,
                      R12*b.R12 + R22*b.R22 + R32*b.R32,
                      R12*b.R13 + R22*b.R23 + R32*b.R33,
                      R13*b.R11 + R23*b.R21 + R33*b.R31,
                      R13*b.R12 + R23*b.R22 + R33*b.R32,
                      R13*b.R13 + R23*b.R23 + R33*b.R33);
  }
template<class T>
TkRotation TkRotation< T >::multiplyInverse ( const TkRotation< T > &  b) const [inline]

Definition at line 205 of file newTkRotation.h.

                                                         {
    return rot.transpose()*b.rot;
  }
template<class T>
template<class Scalar >
Basic3DVector<Scalar> TkRotation< T >::multiplyInverse ( const Basic3DVector< Scalar > &  v) const [inline]

Definition at line 183 of file newTkRotation.h.

                                                                               {
    return Basic3DVector<Scalar>( xx()*v.x() + yx()*v.y() + zx()*v.z(),
                                  xy()*v.x() + yy()*v.y() + zy()*v.z(),
                                  xz()*v.x() + yz()*v.y() + zz()*v.z());
  }
template<class T>
Basic3DVector<T> TkRotation< T >::multiplyInverse ( const Basic3DVector< T > &  v) const [inline]

Definition at line 151 of file oldTkRotation.h.

                                                                     {
    return Basic3DVector<T>( R11*v.x() + R21*v.y() + R31*v.z(),
                             R12*v.x() + R22*v.y() + R32*v.z(),
                             R13*v.x() + R23*v.y() + R33*v.z());
  }
template<class T>
template<class Scalar >
Basic3DVector<Scalar> TkRotation< T >::multiplyInverse ( const Basic3DVector< Scalar > &  v) const [inline]

Definition at line 159 of file oldTkRotation.h.

                                                                               {
    return Basic3DVector<Scalar>( R11*v.x() + R21*v.y() + R31*v.z(),
                                  R12*v.x() + R22*v.y() + R32*v.z(),
                                  R13*v.x() + R23*v.y() + R33*v.z());
  }
template<class T>
Basic3DVector<T> TkRotation< T >::multiplyInverse ( const Basic2DVector< T > &  v) const [inline]

Definition at line 171 of file oldTkRotation.h.

                                                                     {
    return Basic3DVector<T>( R11*v.x() + R21*v.y(),
                             R12*v.x() + R22*v.y(),
                             R13*v.x() + R23*v.y());
  }
template<class T>
Basic3DVector<T> TkRotation< T >::operator* ( const Basic3DVector< T > &  v) const [inline]

Definition at line 174 of file newTkRotation.h.

Referenced by TkRotation< align::Scalar >::operator*=().

                                                               {
    return rot.rotate(v.v);
  }
template<class T>
TkRotation TkRotation< T >::operator* ( const TkRotation< T > &  b) const [inline]

Definition at line 179 of file oldTkRotation.h.

                                                   {
    return TkRotation(R11*b.R11 + R12*b.R21 + R13*b.R31,
                      R11*b.R12 + R12*b.R22 + R13*b.R32,
                      R11*b.R13 + R12*b.R23 + R13*b.R33,
                      R21*b.R11 + R22*b.R21 + R23*b.R31,
                      R21*b.R12 + R22*b.R22 + R23*b.R32,
                      R21*b.R13 + R22*b.R23 + R23*b.R33,
                      R31*b.R11 + R32*b.R21 + R33*b.R31,
                      R31*b.R12 + R32*b.R22 + R33*b.R32,
                      R31*b.R13 + R32*b.R23 + R33*b.R33);
  }
template<class T>
TkRotation TkRotation< T >::operator* ( const TkRotation< T > &  b) const [inline]

Definition at line 202 of file newTkRotation.h.

                                                   {
    return rot*b.rot;
  }
template<class T>
Basic3DVector<T> TkRotation< T >::operator* ( const Basic3DVector< T > &  v) const [inline]

Definition at line 145 of file oldTkRotation.h.

                                                               {
    return Basic3DVector<T>( R11*v.x() + R12*v.y() + R13*v.z(),
                             R21*v.x() + R22*v.y() + R23*v.z(),
                             R31*v.x() + R32*v.y() + R33*v.z());
  }
template<class T>
Basic3DVector<T> TkRotation< T >::operator* ( const Basic2DVector< T > &  v) const [inline]

Definition at line 189 of file newTkRotation.h.

                                                               {
    return Basic3DVector<T>( xx()*v.x() + xy()*v.y(),
                             yx()*v.x() + yy()*v.y(),
                             zx()*v.x() + zy()*v.y());
  }
template<class T>
Basic3DVector<T> TkRotation< T >::operator* ( const Basic2DVector< T > &  v) const [inline]

Definition at line 166 of file oldTkRotation.h.

                                                               {
    return Basic3DVector<T>( R11*v.x() + R12*v.y(),
                             R21*v.x() + R22*v.y(),
                             R31*v.x() + R32*v.y());
  }
template<class T>
TkRotation& TkRotation< T >::operator*= ( const TkRotation< T > &  b) [inline]

Definition at line 203 of file oldTkRotation.h.

                                               {
    return *this = operator * (b);
  }
template<class T>
TkRotation& TkRotation< T >::operator*= ( const TkRotation< T > &  b) [inline]

Definition at line 209 of file newTkRotation.h.

                                               {
    return *this = operator * (b);
  }
template<class T>
Basic3DVector<T> TkRotation< T >::rotate ( const Basic3DVector< T > &  v) const [inline]

Definition at line 165 of file newTkRotation.h.

                                                            {
    return rot.rotate(v.v);
  }
template<class T>
TkRotation& TkRotation< T >::rotateAxes ( const Basic3DVector< T > &  newX,
const Basic3DVector< T > &  newY,
const Basic3DVector< T > &  newZ 
) [inline]

Definition at line 212 of file oldTkRotation.h.

                                                        {
    T del = 0.001;

    if (
        
        // the check for right-handedness is not needed since
        // we want to change the handedness when it's left in cmsim
        //
        //       fabs(newZ.x()-w.x()) > del ||
        //       fabs(newZ.y()-w.y()) > del ||
        //       fabs(newZ.z()-w.z()) > del ||
        fabs(newX.mag2()-1.) > del ||
        fabs(newY.mag2()-1.) > del || 
        fabs(newZ.mag2()-1.) > del ||
        fabs(newX.dot(newY)) > del ||
        fabs(newY.dot(newZ)) > del ||
        fabs(newZ.dot(newX)) > del) {
      geometryDetails::TkRotationErr2();
      return *this;
    } else {
      return transform(TkRotation(newX.x(), newY.x(), newZ.x(),
                                  newX.y(), newY.y(), newZ.y(),
                                  newX.z(), newY.z(), newZ.z()));
    }
  }
template<class T>
TkRotation& TkRotation< T >::rotateAxes ( const Basic3DVector< T > &  newX,
const Basic3DVector< T > &  newY,
const Basic3DVector< T > &  newZ 
) [inline]

Definition at line 218 of file newTkRotation.h.

Referenced by RPCGeometryBuilderFromCondDB::build(), CSCGeometryBuilder::buildChamber(), and RPCGeometryBuilderFromDDD::buildGeometry().

                                                        {
    T del = 0.001;
    
    if (
        
        // the check for right-handedness is not needed since
        // we want to change the handedness when it's left in cmsim
        //
        //       fabs(newZ.x()-w.x()) > del ||
        //       fabs(newZ.y()-w.y()) > del ||
        //       fabs(newZ.z()-w.z()) > del ||
        fabs(newX.mag2()-1.) > del ||
        fabs(newY.mag2()-1.) > del || 
        fabs(newZ.mag2()-1.) > del ||
        fabs(newX.dot(newY)) > del ||
        fabs(newY.dot(newZ)) > del ||
        fabs(newZ.dot(newX)) > del) {
      geometryDetails::TkRotationErr2();
      return *this;
    } else {
      return transform(TkRotation(newX.x(), newY.x(), newZ.x(),
                                  newX.y(), newY.y(), newZ.y(),
                                  newX.z(), newY.z(), newZ.z()));
    }
  }
template<class T>
Basic3DVector<T> TkRotation< T >::rotateBack ( const Basic3DVector< T > &  v) const [inline]

Definition at line 169 of file newTkRotation.h.

                                                                {
    return rot.rotateBack(v.v);
  }
template<class T>
TkRotation& TkRotation< T >::transform ( const TkRotation< T > &  b) [inline]

Definition at line 208 of file oldTkRotation.h.

                                             {
    return *this = b.operator * (*this);
  }
template<class T>
TkRotation& TkRotation< T >::transform ( const TkRotation< T > &  b) [inline]

Definition at line 214 of file newTkRotation.h.

Referenced by TkRotation< align::Scalar >::rotateAxes().

                                             {
    return *this = b.operator * (*this);
  }
template<class T>
TkRotation TkRotation< T >::transposed ( ) const [inline]
template<class T>
TkRotation TkRotation< T >::transposed ( ) const [inline]

Definition at line 139 of file oldTkRotation.h.

                                {
      return TkRotation( R11, R21, R31, 
                         R12, R22, R32,
                         R13, R23, R33);
  }
template<class T>
Basic3DVector<T> TkRotation< T >::x ( ) const [inline]

Definition at line 240 of file oldTkRotation.h.

{ return  Basic3DVector<T>(xx(),xy(),xz());}
template<class T>
Basic3DVector<T> TkRotation< T >::x ( ) const [inline]
template<class T>
T TkRotation< T >::xx ( ) const [inline]
template<class T>
T const& TkRotation< T >::xx ( ) const [inline]

Definition at line 245 of file oldTkRotation.h.

{ return R11;} 
template<class T>
T TkRotation< T >::xy ( ) const [inline]
template<class T>
T const& TkRotation< T >::xy ( ) const [inline]

Definition at line 246 of file oldTkRotation.h.

{ return R12;} 
template<class T>
T TkRotation< T >::xz ( ) const [inline]
template<class T>
T const& TkRotation< T >::xz ( ) const [inline]

Definition at line 247 of file oldTkRotation.h.

{ return R13;} 
template<class T>
Basic3DVector<T> TkRotation< T >::y ( ) const [inline]

Definition at line 248 of file newTkRotation.h.

Referenced by JacobianCurvilinearToLocal::compute(), and JacobianLocalToCurvilinear::compute().

{ return rot.axis[1];}
template<class T>
Basic3DVector<T> TkRotation< T >::y ( ) const [inline]

Definition at line 241 of file oldTkRotation.h.

{ return  Basic3DVector<T>(yx(),yy(),yz());}
template<class T>
T TkRotation< T >::yx ( ) const [inline]
template<class T>
T const& TkRotation< T >::yx ( ) const [inline]

Definition at line 248 of file oldTkRotation.h.

{ return R21;} 
template<class T>
T TkRotation< T >::yy ( ) const [inline]
template<class T>
T const& TkRotation< T >::yy ( ) const [inline]

Definition at line 249 of file oldTkRotation.h.

{ return R22;} 
template<class T>
T TkRotation< T >::yz ( ) const [inline]
template<class T>
T const& TkRotation< T >::yz ( ) const [inline]

Definition at line 250 of file oldTkRotation.h.

{ return R23;} 
template<class T>
Basic3DVector<T> TkRotation< T >::z ( ) const [inline]

Definition at line 249 of file newTkRotation.h.

Referenced by JacobianCurvilinearToLocal::compute().

{ return rot.axis[2];}
template<class T>
Basic3DVector<T> TkRotation< T >::z ( ) const [inline]

Definition at line 242 of file oldTkRotation.h.

{ return  Basic3DVector<T>(zx(),zy(),zz());}
template<class T>
T const& TkRotation< T >::zx ( ) const [inline]

Definition at line 251 of file oldTkRotation.h.

{ return R31;} 
template<class T>
T TkRotation< T >::zx ( ) const [inline]
template<class T>
T TkRotation< T >::zy ( ) const [inline]
template<class T>
T const& TkRotation< T >::zy ( ) const [inline]

Definition at line 252 of file oldTkRotation.h.

{ return R32;} 
template<class T>
T TkRotation< T >::zz ( ) const [inline]
template<class T>
T const& TkRotation< T >::zz ( ) const [inline]

Definition at line 253 of file oldTkRotation.h.

{ return R33;} 

Member Data Documentation

template<class T>
T TkRotation< T >::R11 [private]
template<class T>
T TkRotation< T >::R12 [private]
template<class T>
T TkRotation< T >::R13 [private]
template<class T>
T TkRotation< T >::R21 [private]
template<class T>
T TkRotation< T >::R22 [private]
template<class T>
T TkRotation< T >::R23 [private]
template<class T>
T TkRotation< T >::R31 [private]
template<class T>
T TkRotation< T >::R32 [private]
template<class T>
T TkRotation< T >::R33 [private]
template<class T>
mathSSE::Rot3<T> TkRotation< T >::rot [private]