TRotation.cxx

Go to the documentation of this file.
00001 // @(#)root/physics:$Id: TRotation.cxx 20882 2007-11-19 11:31:26Z rdm $
00002 // Author: Peter Malzacher   19/06/99
00003 
00004 //______________________________________________________________________________
00005 //*-*-*-*-*-*-*-*-*-*-*-*The Physics Vector package *-*-*-*-*-*-*-*-*-*-*-*
00006 //*-*                    ==========================                       *
00007 //*-* The Physics Vector package consists of five classes:                *
00008 //*-*   - TVector2                                                        *
00009 //*-*   - TVector3                                                        *
00010 //*-*   - TRotation                                                       *
00011 //*-*   - TLorentzVector                                                  *
00012 //*-*   - TLorentzRotation                                                *
00013 //*-* It is a combination of CLHEPs Vector package written by             *
00014 //*-* Leif Lonnblad, Andreas Nilsson and Evgueni Tcherniaev               *
00015 //*-* and a ROOT package written by Pasha Murat.                          *
00016 //*-* for CLHEP see:  http://wwwinfo.cern.ch/asd/lhc++/clhep/             *
00017 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00018 //BEGIN_HTML <!--
00019 /* -->
00020 <H2>
00021 TRotation</H2>
00022 The TRotation class describes a rotation of objects of the TVector3 class.
00023 It is a 3*3 matrix of Double_t:
00024 
00025 <P><TT>| xx&nbsp; xy&nbsp; xz |</TT>
00026 <BR><TT>| yx&nbsp; yy&nbsp; yz |</TT>
00027 <BR><TT>| zx&nbsp; zy&nbsp; zz |</TT>
00028 
00029 <P>It describes a so called active rotation, i.e. rotation of objects inside
00030 a static system of coordinates. In case you want to rotate the frame and
00031 want to know the coordinates of objects in the rotated system, you should
00032 apply the inverse rotation to the objects. If you want to transform coordinates
00033 from the rotated frame to the original frame you have to apply the direct
00034 transformation.
00035 
00036 <P>A rotation around a specified axis means counterclockwise rotation around
00037 the positive direction of the axis.
00038 <BR>&nbsp;
00039 <H3>
00040 Declaration, Access, Comparisons</H3>
00041 <TT>&nbsp; TRotation r;&nbsp;&nbsp;&nbsp; // r initialized as identity</TT>
00042 <BR><TT>&nbsp; TRotation m(r); // m = r</TT>
00043 
00044 <P>There is no direct way to to set the matrix elements - to ensure that
00045 a <TT>TRotation</TT> object always describes a real rotation. But you can get the
00046 values by the member functions <TT>XX()..ZZ()</TT> or the<TT> (,)</TT>
00047 operator:
00048 
00049 <P><TT>&nbsp; Double_t xx = r.XX();&nbsp;&nbsp;&nbsp;&nbsp; //&nbsp; the
00050 same as xx=r(0,0)</TT>
00051 <BR><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; xx
00052 = r(0,0);</TT>
00053 
00054 <P><TT>&nbsp; if (r==m) {...}&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;// test for equality</TT>
00055 <BR><TT>&nbsp; if (r!=m) {..}&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;// test for inequality</TT>
00056 <BR><TT>&nbsp; if (r.IsIdentity()) {...} // test for identity</TT>
00057 <BR>&nbsp;
00058 <H3>
00059 Rotation around axes</H3>
00060 The following matrices desrcibe counterclockwise rotations around coordinate
00061 axes
00062 
00063 <P><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; | 1&nbsp;&nbsp; 0&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
00064 0&nbsp;&nbsp;&nbsp; |</TT>
00065 <BR><TT>Rx(a) = | 0 cos(a) -sin(a) |</TT>
00066 <BR><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; | 0 sin(a) cos(a)&nbsp;
00067 |</TT>
00068 
00069 <P><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; | cos(a)&nbsp; 0 sin(a)
00070 |</TT>
00071 <BR><TT>Ry(a) = |&nbsp;&nbsp; 0&nbsp;&nbsp;&nbsp;&nbsp; 1&nbsp;&nbsp;&nbsp;
00072 0&nbsp;&nbsp; |</TT>
00073 <BR><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; | -sin(a) 0 cos(a) |</TT>
00074 
00075 <P><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; | cos(a) -sin(a) 0 |</TT>
00076 <BR><TT>Rz(a) = | sin(a) cos(a) 0 |</TT>
00077 <BR><TT>&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; |&nbsp;&nbsp; 0&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
00078 0&nbsp;&nbsp;&nbsp;&nbsp; 1 |</TT>
00079 <BR>and are implemented as member functions <TT>RotateX()</TT>, <TT>RotateY()</TT>
00080 and <TT>RotateZ()</TT>:
00081 
00082 <P><TT>&nbsp; r.RotateX(TMath::Pi()); // rotation around the x-axis</TT>
00083 <H3>
00084 Rotation around arbitary axis</H3>
00085 The member function <TT>Rotate()</TT> allows to rotate around an arbitary vector
00086 (not neccessary a unit one) and returns the result.
00087 
00088 <P><TT>&nbsp; r.Rotate(TMath::Pi()/3,TVector3(3,4,5));</TT>
00089 
00090 <P>It is possible to find a unit vector and an angle, which describe the
00091 same rotation as the current one:
00092 
00093 <P><TT>&nbsp; Double_t angle;</TT>
00094 <BR><TT>&nbsp; TVector3 axis;</TT>
00095 <BR><TT>&nbsp; r.GetAngleAxis(angle,axis);</TT>
00096 <H3>
00097 Rotation of local axes</H3>
00098 Member function <TT>RotateAxes()</TT> adds a rotation of local axes to
00099 the current rotation and returns the result:
00100 
00101 <P><TT>&nbsp; TVector3 newX(0,1,0);</TT>
00102 <BR><TT>&nbsp; TVector3 newY(0,0,1);</TT>
00103 <BR><TT>&nbsp; TVector3 newZ(1,0,0);</TT>
00104 <BR><TT>&nbsp; a.RotateAxes(newX,newY,newZ);</TT>
00105 
00106 <P>Member functions <TT>ThetaX()</TT>, <TT>ThetaY()</TT>, <TT>ThetaZ()</TT>,
00107 <TT>PhiX()</TT>, <TT>PhiY()</TT>,<TT>PhiZ()</TT> return azimuth and polar
00108 angles of the rotated axes:
00109 
00110 <P><TT>&nbsp; Double_t tx,ty,tz,px,py,pz;</TT>
00111 <BR><TT>&nbsp; tx= a.ThetaX();</TT>
00112 <BR><TT>&nbsp; ...</TT>
00113 <BR><TT>&nbsp; pz= a.PhiZ();</TT>
00114 
00115 <H3>
00116 Setting The Rotations</H3>
00117 The member function <TT>SetToIdentity()</TT> will set the rotation object 
00118 to the identity (no rotation).
00119 
00120 With a minor caveat, the Euler angles of the rotation may be set using 
00121 <TT>SetXEulerAngles()</TT> or individually set with <TT>SetXPhi()</TT>, 
00122 <TT>SetXTheta()</TT>, and <TT>SetXPsi()</TT>.  These routines set the Euler 
00123 angles using the X-convention which is defined by a rotation about the Z-axis,
00124 about the new X-axis, and about the new Z-axis.  This is the convention used
00125 in Landau and Lifshitz, Goldstein and other common physics texts.  The 
00126 Y-convention euler angles can be set with <TT>SetYEulerAngles()</TT>,
00127 <TT>SetYPhi()</TT>, <TT>SetYTheta()</TT>, and <TT>SetYPsi()</TT>.  The caveat 
00128 is that Euler angles usually define the rotation of the new coordinate system 
00129 with respect to the original system, however, the TRotation class specifies 
00130 the rotation of the object in the original system (an active rotation).  To 
00131 recover the usual Euler rotations (ie. rotate the system not the object), you 
00132 must take the inverse of the rotation.
00133 
00134 The member functions <TT>SetXAxis()</TT>, <TT>SetYAxis()</TT>, and 
00135 <TT>SetZAxis()</TT> will create a rotation which rotates the requested axis
00136 of the object to be parallel to a vector.  If used with one argument, the 
00137 rotation about that axis is arbitrary.  If used with two arguments, the
00138 second variable defines the <TT>XY</TT>, <TT>YZ</TT>, or <TT>ZX</TT> 
00139 respectively.
00140 
00141 <H3>
00142 Inverse rotation</H3>
00143 <TT>&nbsp; TRotation a,b;</TT>
00144 <BR><TT>&nbsp; ...</TT>
00145 <BR><TT>&nbsp; b = a.Inverse();&nbsp; // b is inverse of a, a is unchanged</TT>
00146 <BR><TT>&nbsp; b = a.Invert();&nbsp;&nbsp; // invert a and set b = a</TT>
00147 <H3>
00148 Compound Rotations</H3>
00149 The <TT>operator *</TT> has been implemented in a way that follows the 
00150 mathematical notation of a product of the two matrices which describe the 
00151 two consecutive rotations. Therefore the second rotation should be placed 
00152 first:
00153 
00154 <P><TT>&nbsp; r = r2 * r1;</TT>
00155 <H3>
00156 Rotation of TVector3</H3>
00157 The TRotation class provides an <TT>operator *</TT> which allows to express
00158 a rotation of a <TT>TVector3</TT> analog to the mathematical notation
00159 
00160 <P><TT>&nbsp; | x' |&nbsp;&nbsp; | xx xy xz | | x |</TT>
00161 <BR><TT>&nbsp; | y' | = | yx yy yz | | y |</TT>
00162 <BR><TT>&nbsp; | z' |&nbsp;&nbsp; | zx zy zz | | z |</TT><TT></TT>
00163 
00164 <P>e.g.:
00165 
00166 <P><TT>&nbsp; TVector3 v(1,1,1);</TT>
00167 <BR><TT>&nbsp; v = r * v;</TT><TT></TT>
00168 
00169 <P>You can also use the <TT>Transform()</TT> member function or the o<TT>perator
00170 *=</TT> of the
00171 <BR>TVector3 class:<TT></TT>
00172 
00173 <P><TT>&nbsp; TVector3 v;</TT>
00174 <BR><TT>&nbsp; TRotation r;</TT>
00175 <BR><TT>&nbsp; v.Transform(r);</TT>
00176 <BR><TT>&nbsp; v *= r;&nbsp; //Attention v = r * v</TT>
00177 <!--*/
00178 // -->END_HTML
00179 //
00180 
00181 #include "TRotation.h"
00182 #include "TMath.h"
00183 #include "TQuaternion.h"
00184 #include "TError.h"
00185 
00186 ClassImp(TRotation)
00187 
00188 #define TOLERANCE (1.0E-6)
00189 
00190 TRotation::TRotation()
00191 : fxx(1.0), fxy(0.0), fxz(0.0), fyx(0.0), fyy(1.0), fyz(0.0),
00192   fzx(0.0), fzy(0.0), fzz(1.0) {}
00193 
00194 TRotation::TRotation(const TRotation & m) : TObject(m),
00195   fxx(m.fxx), fxy(m.fxy), fxz(m.fxz), fyx(m.fyx), fyy(m.fyy), fyz(m.fyz),
00196   fzx(m.fzx), fzy(m.fzy), fzz(m.fzz) {}
00197 
00198 TRotation::TRotation(Double_t mxx, Double_t mxy, Double_t mxz,
00199                          Double_t myx, Double_t myy, Double_t myz,
00200                          Double_t mzx, Double_t mzy, Double_t mzz)
00201 : fxx(mxx), fxy(mxy), fxz(mxz), fyx(myx), fyy(myy), fyz(myz),
00202   fzx(mzx), fzy(mzy), fzz(mzz) {}
00203 
00204 
00205 Double_t TRotation::operator() (int i, int j) const {
00206    //dereferencing operator const
00207    if (i == 0) {
00208       if (j == 0) { return fxx; }
00209       if (j == 1) { return fxy; }
00210       if (j == 2) { return fxz; }
00211    } else if (i == 1) {
00212       if (j == 0) { return fyx; }
00213       if (j == 1) { return fyy; }
00214       if (j == 2) { return fyz; }
00215    } else if (i == 2) {
00216       if (j == 0) { return fzx; }
00217       if (j == 1) { return fzy; }
00218       if (j == 2) { return fzz; }
00219    }
00220   
00221    Warning("operator()(i,j)", "bad indices (%d , %d)",i,j);
00222 
00223    return 0.0;
00224 }
00225 
00226 TRotation TRotation::operator* (const TRotation & b) const {
00227    //multiplication operator
00228    return TRotation(fxx*b.fxx + fxy*b.fyx + fxz*b.fzx,
00229                     fxx*b.fxy + fxy*b.fyy + fxz*b.fzy,
00230                     fxx*b.fxz + fxy*b.fyz + fxz*b.fzz,
00231                     fyx*b.fxx + fyy*b.fyx + fyz*b.fzx,
00232                     fyx*b.fxy + fyy*b.fyy + fyz*b.fzy,
00233                     fyx*b.fxz + fyy*b.fyz + fyz*b.fzz,
00234                     fzx*b.fxx + fzy*b.fyx + fzz*b.fzx,
00235                     fzx*b.fxy + fzy*b.fyy + fzz*b.fzy,
00236                     fzx*b.fxz + fzy*b.fyz + fzz*b.fzz);
00237 }
00238 
00239 //_____________________________________
00240 TRotation::TRotation(const TQuaternion & Q) {
00241    // Constructor for a rotation based on a Quaternion
00242    // if magnitude of quaternion is null, creates identity rotation
00243    // if quaternion is non-unit, creates rotation corresponding to the normalized (unit) quaternion
00244 
00245 
00246    double two_r2 = 2 * Q.fRealPart * Q.fRealPart;
00247    double two_x2 = 2 * Q.fVectorPart.X() * Q.fVectorPart.X();
00248    double two_y2 = 2 * Q.fVectorPart.Y() * Q.fVectorPart.Y();
00249    double two_z2 = 2 * Q.fVectorPart.Z() * Q.fVectorPart.Z();
00250    double two_xy = 2 * Q.fVectorPart.X() * Q.fVectorPart.Y();
00251    double two_xz = 2 * Q.fVectorPart.X() * Q.fVectorPart.Z();
00252    double two_xr = 2 * Q.fVectorPart.X() * Q.fRealPart;
00253    double two_yz = 2 * Q.fVectorPart.Y() * Q.fVectorPart.Z();
00254    double two_yr = 2 * Q.fVectorPart.Y() * Q.fRealPart;
00255    double two_zr = 2 * Q.fVectorPart.Z() * Q.fRealPart;
00256 
00257    // protect agains zero quaternion
00258    double mag2 = Q.QMag2();
00259    if (mag2 > 0) {
00260 
00261       // diago + identity
00262       fxx = two_r2 + two_x2;
00263       fyy = two_r2 + two_y2;
00264       fzz = two_r2 + two_z2;
00265 
00266       //        line 0 column 1 and conjugate
00267       fxy = two_xy - two_zr;
00268       fyx = two_xy + two_zr;
00269 
00270       //        line 0 column 2 and conjugate
00271       fxz = two_xz + two_yr;
00272       fzx = two_xz - two_yr;
00273 
00274       //        line 1 column 2 and conjugate
00275       fyz = two_yz - two_xr;
00276       fzy = two_yz + two_xr;
00277 
00278       // protect agains non-unit quaternion 
00279       if (TMath::Abs(mag2-1) > 1e-10) {
00280          fxx /= mag2;
00281          fyy /= mag2;
00282          fzz /= mag2;
00283          fxy /= mag2;
00284          fyx /= mag2;
00285          fxz /= mag2;
00286          fzx /= mag2;
00287          fyz /= mag2;
00288          fzy /= mag2;
00289       }
00290 
00291       // diago : remove identity
00292       fxx -= 1;
00293       fyy -= 1;
00294       fzz -= 1;
00295 
00296 
00297    } else {
00298       // Identity
00299 
00300       fxx = fyy = fzz = 1;
00301       fxy = fyx = fxz = fzx = fyz = fzy = 0;
00302 
00303    }
00304 
00305 }
00306 
00307 TRotation & TRotation::Rotate(Double_t a, const TVector3& axis) {
00308    //rotate along an axis
00309    if (a != 0.0) {
00310       Double_t ll = axis.Mag();
00311       if (ll == 0.0) {
00312          Warning("Rotate(angle,axis)"," zero axis");
00313       } else {
00314          Double_t sa = TMath::Sin(a), ca = TMath::Cos(a);
00315          Double_t dx = axis.X()/ll, dy = axis.Y()/ll, dz = axis.Z()/ll;
00316          TRotation m(
00317              ca+(1-ca)*dx*dx,          (1-ca)*dx*dy-sa*dz,    (1-ca)*dx*dz+sa*dy,
00318              (1-ca)*dy*dx+sa*dz, ca+(1-ca)*dy*dy,          (1-ca)*dy*dz-sa*dx,
00319              (1-ca)*dz*dx-sa*dy,    (1-ca)*dz*dy+sa*dx, ca+(1-ca)*dz*dz );
00320          Transform(m);
00321       }
00322    }
00323    return *this;
00324 }
00325 
00326 TRotation & TRotation::RotateX(Double_t a) {
00327    //rotate around x
00328    Double_t c = TMath::Cos(a);
00329    Double_t s = TMath::Sin(a);
00330    Double_t x = fyx, y = fyy, z = fyz;
00331    fyx = c*x - s*fzx;
00332    fyy = c*y - s*fzy;
00333    fyz = c*z - s*fzz;
00334    fzx = s*x + c*fzx;
00335    fzy = s*y + c*fzy;
00336    fzz = s*z + c*fzz;
00337    return *this;
00338 }
00339 
00340 TRotation & TRotation::RotateY(Double_t a){
00341    //rotate around y
00342    Double_t c = TMath::Cos(a);
00343    Double_t s = TMath::Sin(a);
00344    Double_t x = fzx, y = fzy, z = fzz;
00345    fzx = c*x - s*fxx;
00346    fzy = c*y - s*fxy;
00347    fzz = c*z - s*fxz;
00348    fxx = s*x + c*fxx;
00349    fxy = s*y + c*fxy;
00350    fxz = s*z + c*fxz;
00351    return *this;
00352 }
00353 
00354 TRotation & TRotation::RotateZ(Double_t a) {
00355    //rotate around z
00356    Double_t c = TMath::Cos(a);
00357    Double_t s = TMath::Sin(a);
00358    Double_t x = fxx, y = fxy, z = fxz;
00359    fxx = c*x - s*fyx;
00360    fxy = c*y - s*fyy;
00361    fxz = c*z - s*fyz;
00362    fyx = s*x + c*fyx;
00363    fyy = s*y + c*fyy;
00364    fyz = s*z + c*fyz;
00365    return *this;
00366 }
00367 
00368 TRotation & TRotation::RotateAxes(const TVector3 &newX,
00369                                   const TVector3 &newY,
00370                                   const TVector3 &newZ) {
00371    //rotate axes
00372    Double_t del = 0.001;
00373    TVector3 w = newX.Cross(newY);
00374 
00375    if (TMath::Abs(newZ.X()-w.X()) > del ||
00376        TMath::Abs(newZ.Y()-w.Y()) > del ||
00377        TMath::Abs(newZ.Z()-w.Z()) > del ||
00378        TMath::Abs(newX.Mag2()-1.) > del ||
00379        TMath::Abs(newY.Mag2()-1.) > del ||
00380        TMath::Abs(newZ.Mag2()-1.) > del ||
00381        TMath::Abs(newX.Dot(newY)) > del ||
00382        TMath::Abs(newY.Dot(newZ)) > del ||
00383        TMath::Abs(newZ.Dot(newX)) > del) {
00384       Warning("RotateAxes","bad axis vectors");
00385       return *this;
00386    } else {
00387       return Transform(TRotation(newX.X(), newY.X(), newZ.X(),
00388                                  newX.Y(), newY.Y(), newZ.Y(),
00389                                  newX.Z(), newY.Z(), newZ.Z()));
00390    }
00391 }
00392 
00393 Double_t TRotation::PhiX() const {
00394    //return Phi
00395    return (fyx == 0.0 && fxx == 0.0) ? 0.0 : TMath::ATan2(fyx,fxx);
00396 }
00397 
00398 Double_t TRotation::PhiY() const {
00399    //return Phi
00400    return (fyy == 0.0 && fxy == 0.0) ? 0.0 : TMath::ATan2(fyy,fxy);
00401 }
00402 
00403 Double_t TRotation::PhiZ() const {
00404    //return Phi
00405    return (fyz == 0.0 && fxz == 0.0) ? 0.0 : TMath::ATan2(fyz,fxz);
00406 }
00407 
00408 Double_t TRotation::ThetaX() const {
00409    //return Phi
00410    return TMath::ACos(fzx);
00411 }
00412 
00413 Double_t TRotation::ThetaY() const {
00414    //return Theta
00415    return TMath::ACos(fzy);
00416 }
00417 
00418 Double_t TRotation::ThetaZ() const {
00419    //return Theta
00420    return TMath::ACos(fzz);
00421 }
00422 
00423 void TRotation::AngleAxis(Double_t &angle, TVector3 &axis) const {
00424    //rotation defined by an angle and a vector
00425    Double_t cosa  = 0.5*(fxx+fyy+fzz-1);
00426    Double_t cosa1 = 1-cosa;
00427    if (cosa1 <= 0) {
00428       angle = 0;
00429       axis  = TVector3(0,0,1);
00430    } else {
00431       Double_t x=0, y=0, z=0;
00432       if (fxx > cosa) x = TMath::Sqrt((fxx-cosa)/cosa1);
00433       if (fyy > cosa) y = TMath::Sqrt((fyy-cosa)/cosa1);
00434       if (fzz > cosa) z = TMath::Sqrt((fzz-cosa)/cosa1);
00435       if (fzy < fyz)  x = -x;
00436       if (fxz < fzx)  y = -y;
00437       if (fyx < fxy)  z = -z;
00438       angle = TMath::ACos(cosa);
00439       axis  = TVector3(x,y,z);
00440    }
00441 }
00442 
00443 TRotation & TRotation::SetXEulerAngles(Double_t phi,
00444                                       Double_t theta,
00445                                       Double_t psi) {
00446    // Rotate using the x-convention (Landau and Lifshitz, Goldstein, &c) by 
00447    // doing the explicit rotations.  This is slightly less efficient than 
00448    // directly applying the rotation, but makes the code much clearer.  My
00449    // presumption is that this code is not going to be a speed bottle neck.
00450 
00451    SetToIdentity();
00452    RotateZ(phi);
00453    RotateX(theta);
00454    RotateZ(psi);
00455   
00456    return *this;
00457 }
00458 
00459 TRotation & TRotation::SetYEulerAngles(Double_t phi,
00460                                        Double_t theta,
00461                                        Double_t psi) {
00462    // Rotate using the y-convention.
00463     
00464    SetToIdentity();
00465    RotateZ(phi);
00466    RotateY(theta);
00467    RotateZ(psi);
00468    return *this;
00469 }
00470 
00471 TRotation & TRotation::RotateXEulerAngles(Double_t phi,
00472                                          Double_t theta,
00473                                          Double_t psi) {
00474    // Rotate using the x-convention.
00475    TRotation euler;
00476    euler.SetXEulerAngles(phi,theta,psi);
00477    return Transform(euler);
00478 }
00479 
00480 TRotation & TRotation::RotateYEulerAngles(Double_t phi,
00481                                           Double_t theta,
00482                                           Double_t psi) {
00483    // Rotate using the y-convention.
00484    TRotation euler;
00485    euler.SetYEulerAngles(phi,theta,psi);
00486    return Transform(euler);
00487 }
00488 
00489 void TRotation::SetXPhi(Double_t phi) {
00490    //set XPhi
00491    SetXEulerAngles(phi,GetXTheta(),GetXPsi());
00492 }
00493 
00494 void TRotation::SetXTheta(Double_t theta) {
00495    //set XTheta
00496    SetXEulerAngles(GetXPhi(),theta,GetXPsi());
00497 }
00498 
00499 void TRotation::SetXPsi(Double_t psi) {
00500    //set XPsi
00501    SetXEulerAngles(GetXPhi(),GetXTheta(),psi);
00502 }
00503 
00504 void TRotation::SetYPhi(Double_t phi) {
00505    //set YPhi
00506    SetYEulerAngles(phi,GetYTheta(),GetYPsi());
00507 }
00508 
00509 void TRotation::SetYTheta(Double_t theta) {
00510    //set YTheta
00511    SetYEulerAngles(GetYPhi(),theta,GetYPsi());
00512 }
00513 
00514 void TRotation::SetYPsi(Double_t psi) {
00515    //set YPsi
00516    SetYEulerAngles(GetYPhi(),GetYTheta(),psi);
00517 }
00518 
00519 Double_t TRotation::GetXPhi(void) const {
00520    //return phi angle
00521    Double_t finalPhi;
00522 
00523    Double_t s2 =  1.0 - fzz*fzz;
00524    if (s2 < 0) {
00525       Warning("GetPhi()"," |fzz| > 1 ");
00526       s2 = 0;
00527    }
00528    const Double_t sinTheta = TMath::Sqrt(s2);
00529 
00530    if (sinTheta != 0) {
00531       const Double_t cscTheta = 1/sinTheta;
00532       Double_t cosAbsPhi =  fzy * cscTheta;
00533       if ( TMath::Abs(cosAbsPhi) > 1 ) {        // NaN-proofing
00534          Warning("GetPhi()","finds | cos phi | > 1");
00535          cosAbsPhi = 1;
00536       }
00537       const Double_t absPhi = TMath::ACos(cosAbsPhi);
00538       if (fzx > 0) {
00539          finalPhi = absPhi;
00540       } else if (fzx < 0) {
00541          finalPhi = -absPhi;
00542       } else if (fzy > 0) {
00543          finalPhi = 0.0;
00544       } else {
00545          finalPhi = TMath::Pi();
00546       }
00547    } else {              // sinTheta == 0 so |Fzz| = 1
00548       const Double_t absPhi = .5 * TMath::ACos (fxx);
00549       if (fxy > 0) {
00550          finalPhi =  -absPhi;
00551       } else if (fxy < 0) {
00552          finalPhi =   absPhi;
00553       } else if (fxx>0) {
00554          finalPhi = 0.0;
00555       } else {
00556          finalPhi = fzz * TMath::PiOver2();
00557       }
00558    }
00559    return finalPhi;
00560 }
00561 
00562 Double_t TRotation::GetYPhi(void) const {
00563    //return YPhi
00564    return GetXPhi() + TMath::Pi()/2.0;
00565 }
00566 
00567 Double_t TRotation::GetXTheta(void) const {
00568    //return XTheta
00569    return  ThetaZ();
00570 }
00571 
00572 Double_t TRotation::GetYTheta(void) const {
00573    //return YTheta
00574    return  ThetaZ();
00575 }
00576 
00577 Double_t TRotation::GetXPsi(void) const {
00578    //Get psi angle
00579    double finalPsi = 0.0;
00580 
00581    Double_t s2 =  1.0 - fzz*fzz;
00582    if (s2 < 0) {
00583       Warning("GetPsi()"," |fzz| > 1 ");
00584       s2 = 0;
00585    }
00586    const Double_t sinTheta = TMath::Sqrt(s2);
00587 
00588    if (sinTheta != 0) {
00589       const Double_t cscTheta = 1/sinTheta;
00590       Double_t cosAbsPsi =  - fyz * cscTheta;
00591       if ( TMath::Abs(cosAbsPsi) > 1 ) {        // NaN-proofing
00592          Warning("GetPsi()","| cos psi | > 1 ");
00593          cosAbsPsi = 1;
00594       }
00595       const Double_t absPsi = TMath::ACos(cosAbsPsi);
00596       if (fxz > 0) {
00597          finalPsi = absPsi;
00598       } else if (fxz < 0) {
00599          finalPsi = -absPsi;
00600       } else {
00601          finalPsi = (fyz < 0) ? 0 : TMath::Pi();
00602       }
00603    } else {              // sinTheta == 0 so |Fzz| = 1
00604       Double_t absPsi = fxx;
00605       if ( TMath::Abs(fxx) > 1 ) {        // NaN-proofing
00606          Warning("GetPsi()","| fxx | > 1 ");
00607          absPsi = 1;
00608       }
00609       absPsi = .5 * TMath::ACos (absPsi);
00610       if (fyx > 0) {
00611          finalPsi = absPsi;
00612       } else if (fyx < 0) {
00613          finalPsi = -absPsi;
00614       } else {
00615          finalPsi = (fxx > 0) ? 0 : TMath::PiOver2();
00616       }
00617    }
00618    return finalPsi;
00619 }
00620 
00621 Double_t TRotation::GetYPsi(void) const {
00622    //return YPsi
00623    return GetXPsi() - TMath::Pi()/2;
00624 }
00625 
00626 TRotation & TRotation::SetXAxis(const TVector3& axis, 
00627                                 const TVector3& xyPlane) {
00628    //set X axis
00629    TVector3 xAxis(xyPlane);
00630    TVector3 yAxis;
00631    TVector3 zAxis(axis);
00632    MakeBasis(xAxis,yAxis,zAxis);
00633    fxx = zAxis.X();  fyx = zAxis.Y();  fzx = zAxis.Z();
00634    fxy = xAxis.X();  fyy = xAxis.Y();  fzy = xAxis.Z();
00635    fxz = yAxis.X();  fyz = yAxis.Y();  fzz = yAxis.Z();
00636    return *this;
00637 }
00638 
00639 TRotation & TRotation::SetXAxis(const TVector3& axis) {
00640    //set X axis
00641    TVector3 xyPlane(0.0,1.0,0.0);
00642    return SetXAxis(axis,xyPlane);
00643 }
00644 
00645 TRotation & TRotation::SetYAxis(const TVector3& axis, 
00646                                 const TVector3& yzPlane) {
00647    //set Y axis
00648    TVector3 xAxis(yzPlane);
00649    TVector3 yAxis;
00650    TVector3 zAxis(axis);
00651    MakeBasis(xAxis,yAxis,zAxis);
00652    fxx = yAxis.X();  fyx = yAxis.Y();  fzx = yAxis.Z();
00653    fxy = zAxis.X();  fyy = zAxis.Y();  fzy = zAxis.Z();
00654    fxz = xAxis.X();  fyz = xAxis.Y();  fzz = xAxis.Z();
00655    return *this;
00656 }
00657 
00658 TRotation & TRotation::SetYAxis(const TVector3& axis) {
00659    //set Y axis
00660    TVector3 yzPlane(0.0,0.0,1.0);
00661    return SetYAxis(axis,yzPlane);
00662 }
00663 
00664 TRotation & TRotation::SetZAxis(const TVector3& axis, 
00665                                 const TVector3& zxPlane) {
00666    //set Z axis
00667    TVector3 xAxis(zxPlane);
00668    TVector3 yAxis;
00669    TVector3 zAxis(axis);
00670    MakeBasis(xAxis,yAxis,zAxis);
00671    fxx = xAxis.X();  fyx = xAxis.Y();  fzx = xAxis.Z();
00672    fxy = yAxis.X();  fyy = yAxis.Y();  fzy = yAxis.Z();
00673    fxz = zAxis.X();  fyz = zAxis.Y();  fzz = zAxis.Z();
00674    return *this;
00675 }
00676 
00677 TRotation & TRotation::SetZAxis(const TVector3& axis) {
00678    //set Z axis
00679    TVector3 zxPlane(1.0,0.0,0.0);
00680    return SetZAxis(axis,zxPlane);
00681 }
00682 
00683 void TRotation::MakeBasis(TVector3& xAxis,
00684                           TVector3& yAxis,
00685                           TVector3& zAxis) const {
00686    // Make the Z axis into a unit variable. 
00687    Double_t zmag = zAxis.Mag();
00688    if (zmag<TOLERANCE) {
00689       Warning("MakeBasis(X,Y,Z)","non-zero Z Axis is required");
00690    }
00691    zAxis *= (1.0/zmag);
00692 
00693    Double_t xmag = xAxis.Mag();
00694    if (xmag<TOLERANCE*zmag) {
00695       xAxis = zAxis.Orthogonal();
00696       xmag = 1.0;
00697    }
00698 
00699    // Find the Y axis
00700    yAxis = zAxis.Cross(xAxis)*(1.0/xmag);
00701    Double_t ymag = yAxis.Mag();
00702    if (ymag<TOLERANCE*zmag) {
00703       yAxis = zAxis.Orthogonal();
00704    } else {
00705       yAxis *= (1.0/ymag);
00706    }
00707 
00708    xAxis = yAxis.Cross(zAxis);
00709 }

Generated on Tue Jul 5 14:37:48 2011 for ROOT_528-00b_version by  doxygen 1.5.1