TVector3.cxx

Go to the documentation of this file.
00001 // @(#)root/physics:$Id: TVector3.cxx 24127 2008-06-04 07:59:20Z moneta $
00002 // Author: Pasha Murat, Peter Malzacher   12/02/99
00003 //    Aug 11 1999: added Pt == 0 guard to Eta()
00004 //    Oct  8 1999: changed Warning to Error and
00005 //                 return fX in Double_t & operator()
00006 //    Oct 20 1999: Bug fix: sign in PseudoRapidity
00007 //                 Warning-> Error in Double_t operator()
00008 
00009 //______________________________________________________________________________
00010 //*-*-*-*-*-*-*-*-*-*-*-*The Physics Vector package *-*-*-*-*-*-*-*-*-*-*-*
00011 //*-*                    ==========================                       *
00012 //*-* The Physics Vector package consists of five classes:                *
00013 //*-*   - TVector2                                                        *
00014 //*-*   - TVector3                                                        *
00015 //*-*   - TRotation                                                       *
00016 //*-*   - TLorentzVector                                                  *
00017 //*-*   - TLorentzRotation                                                *
00018 //*-* It is a combination of CLHEPs Vector package written by             *
00019 //*-* Leif Lonnblad, Andreas Nilsson and Evgueni Tcherniaev               *
00020 //*-* and a ROOT package written by Pasha Murat.                          *
00021 //*-* for CLHEP see:  http://wwwinfo.cern.ch/asd/lhc++/clhep/             *
00022 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00023 //BEGIN_HTML <!--
00024 /* -->
00025 <H2>
00026 TVector3</H2>
00027 <TT>TVector3</TT> is a general three vector class, which can be used for
00028 the description of different vectors in 3D.
00029 <H3>
00030 Declaration / Access to the components</H3>
00031 <TT>TVector3</TT> has been implemented as a vector of three <TT>Double_t</TT>
00032 variables, representing the cartesian coordinates. By default all components
00033 are initialized to zero:
00034 
00035 <P><TT>&nbsp; TVector3 v1;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; //
00036 v1 = (0,0,0)</TT>
00037 <BR><TT>&nbsp; TVector3 v2(1);&nbsp;&nbsp;&nbsp;&nbsp; // v2 = (1,0,0)</TT>
00038 <BR><TT>&nbsp; TVector3 v3(1,2,3); // v3 = (1,2,3)</TT>
00039 <BR><TT>&nbsp; TVector3 v4(v2);&nbsp;&nbsp;&nbsp; // v4 = v2</TT>
00040 
00041 <P>It is also possible (but not recommended) to initialize a <TT>TVector3</TT>
00042 with a <TT>Double_t</TT> or <TT>Float_t</TT> C array.
00043 
00044 <P>You can get the basic components either by name or by index using <TT>operator()</TT>:
00045 
00046 <P><TT>&nbsp; xx = v1.X();&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp; xx =
00047 v1(0);</TT>
00048 <BR><TT>&nbsp; yy = v1.Y();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
00049 yy = v1(1);</TT>
00050 <BR><TT>&nbsp; zz = v1.Z();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
00051 zz = v1(2);</TT>
00052 
00053 <P>The memberfunctions <TT>SetX()</TT>, <TT>SetY()</TT>, <TT>SetZ()</TT>
00054 and<TT> SetXYZ()</TT> allow to set the components:
00055 
00056 <P><TT>&nbsp; v1.SetX(1.); v1.SetY(2.); v1.SetZ(3.);</TT>
00057 <BR><TT>&nbsp; v1.SetXYZ(1.,2.,3.);</TT>
00058 <BR>&nbsp;
00059 <H3>
00060 Noncartesian coordinates</H3>
00061 To get information on the <TT>TVector3</TT> in spherical (rho,phi,theta)
00062 or cylindrical (z,r,theta) coordinates, the
00063 <BR>the member functions <TT>Mag()</TT> (=magnitude=rho in spherical coordinates),
00064 <TT>Mag2()</TT>, <TT>Theta()</TT>, <TT>CosTheta()</TT>, <TT>Phi()</TT>,
00065 <TT>Perp()</TT> (the transverse component=r in cylindrical coordinates),
00066 <TT>Perp2()</TT> can be used:
00067 
00068 <P><TT>&nbsp; Double_t m&nbsp; = v.Mag();&nbsp;&nbsp;&nbsp; // get magnitude
00069 (=rho=Sqrt(x*x+y*y+z*z)))</TT>
00070 <BR><TT>&nbsp; Double_t m2 = v.Mag2();&nbsp;&nbsp; // get magnitude squared</TT>
00071 <BR><TT>&nbsp; Double_t t&nbsp; = v.Theta();&nbsp; // get polar angle</TT>
00072 <BR><TT>&nbsp; Double_t ct = v.CosTheta();// get cos of theta</TT>
00073 <BR><TT>&nbsp; Double_t p&nbsp; = v.Phi();&nbsp;&nbsp;&nbsp; // get azimuth
00074 angle</TT>
00075 <BR><TT>&nbsp; Double_t pp = v.Perp();&nbsp;&nbsp; // get transverse component</TT>
00076 <BR><TT>&nbsp; Double_t pp2= v.Perp2();&nbsp; // get transvers component
00077 squared</TT>
00078 
00079 <P>It is also possible to get the transverse component with respect to
00080 another vector:
00081 
00082 <P><TT>&nbsp; Double_t ppv1 = v.Perp(v1);</TT>
00083 <BR><TT>&nbsp; Double_t pp2v1 = v.Perp2(v1);</TT>
00084 
00085 <P>The pseudo-rapidity ( eta=-ln (tan (theta/2)) ) can be obtained by <TT>Eta()</TT>
00086 or <TT>PseudoRapidity()</TT>:
00087 <BR>&nbsp;
00088 <BR><TT>&nbsp; Double_t eta = v.PseudoRapidity();</TT>
00089 
00090 <P>There are set functions to change one of the noncartesian coordinates:
00091 
00092 <P><TT>&nbsp; v.SetTheta(.5); // keeping rho and phi</TT>
00093 <BR><TT>&nbsp; v.SetPhi(.8);&nbsp;&nbsp; // keeping rho and theta</TT>
00094 <BR><TT>&nbsp; v.SetMag(10.);&nbsp; // keeping theta and phi</TT>
00095 <BR><TT>&nbsp; v.SetPerp(3.);&nbsp; // keeping z and phi</TT>
00096 <BR>&nbsp;
00097 <H3>
00098 Arithmetic / Comparison</H3>
00099 The <TT>TVector3</TT> class provides the operators to add, subtract, scale and compare
00100 vectors:
00101 
00102 <P><TT>&nbsp; v3&nbsp; = -v1;</TT>
00103 <BR><TT>&nbsp; v1&nbsp; = v2+v3;</TT>
00104 <BR><TT>&nbsp; v1 += v3;</TT>
00105 <BR><TT>&nbsp; v1&nbsp; = v1 - v3</TT>
00106 <BR><TT>&nbsp; v1 -= v3;</TT>
00107 <BR><TT>&nbsp; v1 *= 10;</TT>
00108 <BR><TT>&nbsp; v1&nbsp; = 5*v2;</TT>
00109 
00110 <P><TT>&nbsp; if(v1==v2) {...}</TT>
00111 <BR><TT>&nbsp; if(v1!=v2) {...}</TT>
00112 <BR>&nbsp;
00113 <H3>
00114 Related Vectors</H3>
00115 <TT>&nbsp; v2 = v1.Unit();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // get unit
00116 vector parallel to v1</TT>
00117 <BR><TT>&nbsp; v2 = v1.Orthogonal(); // get vector orthogonal to v1</TT>
00118 <H3>
00119 Scalar and vector products</H3>
00120 <TT>&nbsp; s = v1.Dot(v2);&nbsp;&nbsp; // scalar product</TT>
00121 <BR><TT>&nbsp; s = v1 * v2;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // scalar product</TT>
00122 <BR><TT>&nbsp; v = v1.Cross(v2); // vector product</TT>
00123 <H3>
00124 &nbsp;Angle between two vectors</H3>
00125 <TT>&nbsp; Double_t a = v1.Angle(v2);</TT>
00126 <H3>
00127 Rotations</H3>
00128 
00129 <H5>
00130 Rotation around axes</H5>
00131 <TT>&nbsp; v.RotateX(.5);</TT>
00132 <BR><TT>&nbsp; v.RotateY(TMath::Pi());</TT>
00133 <BR><TT>&nbsp; v.RotateZ(angle);</TT>
00134 <H5>
00135 Rotation around a vector</H5>
00136 <TT>&nbsp; v1.Rotate(TMath::Pi()/4, v2); // rotation around v2</TT>
00137 <H5>
00138 Rotation by TRotation</H5>
00139 <TT>TVector3</TT> objects can be rotated by objects of the <TT>TRotation</TT>
00140 class using the <TT>Transform()</TT> member functions,
00141 <BR>the <TT>operator *=</TT> or the <TT>operator *</TT> of the TRotation
00142 class:
00143 
00144 <P><TT>&nbsp; TRotation m;</TT>
00145 <BR><TT>&nbsp; ...</TT>
00146 <BR><TT>&nbsp; v1.transform(m);</TT>
00147 <BR><TT>&nbsp; v1 = m*v1;</TT>
00148 <BR><TT>&nbsp; v1 *= m; // Attention v1 = m*v1</TT>
00149 <H5>
00150 Transformation from rotated frame</H5>
00151 <TT>&nbsp; TVector3 direction = v.Unit()</TT>
00152 <BR><TT>&nbsp; v1.RotateUz(direction); // direction must be TVector3 of
00153 unit length</TT>
00154 
00155 <P>transforms v1 from the rotated frame (z' parallel to direction, x' in
00156 the theta plane and y' in the xy plane as well as perpendicular to the
00157 theta plane) to the (x,y,z) frame.
00158 
00159 <!--*/
00160 // -->END_HTML
00161 //
00162 
00163 #include "TVector3.h"
00164 #include "TRotation.h"
00165 #include "TMath.h"
00166 #include "TClass.h"
00167 
00168 ClassImp(TVector3)
00169 
00170 //______________________________________________________________________________
00171 TVector3::TVector3(const TVector3 & p) : TObject(p),
00172   fX(p.fX), fY(p.fY), fZ(p.fZ) {}
00173 
00174 TVector3::TVector3(Double_t xx, Double_t yy, Double_t zz)
00175 : fX(xx), fY(yy), fZ(zz) {}
00176 
00177 TVector3::TVector3(const Double_t * x0)
00178 : fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
00179 
00180 TVector3::TVector3(const Float_t * x0)
00181 : fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
00182 
00183 TVector3::~TVector3() {}
00184 
00185 //______________________________________________________________________________
00186 Double_t TVector3::operator () (int i) const {
00187    //dereferencing operator const
00188    switch(i) {
00189       case 0:
00190          return fX;
00191       case 1:
00192          return fY;
00193       case 2:
00194          return fZ;
00195       default:
00196          Error("operator()(i)", "bad index (%d) returning 0",i);
00197    }
00198    return 0.;
00199 }
00200 
00201 //______________________________________________________________________________
00202 Double_t & TVector3::operator () (int i) {
00203    //dereferencing operator
00204    switch(i) {
00205       case 0:
00206          return fX;
00207       case 1:
00208          return fY;
00209       case 2:
00210          return fZ;
00211       default:
00212          Error("operator()(i)", "bad index (%d) returning &fX",i);
00213    }
00214    return fX;
00215 }
00216 
00217 //______________________________________________________________________________
00218 TVector3 & TVector3::operator *= (const TRotation & m){
00219    //multiplication operator
00220    return *this = m * (*this);
00221 }
00222 
00223 //______________________________________________________________________________
00224 TVector3 & TVector3::Transform(const TRotation & m) {
00225    //transform this vector with a TRotation
00226    return *this = m * (*this);
00227 }
00228 
00229 //______________________________________________________________________________
00230 Double_t TVector3::Angle(const TVector3 & q) const 
00231 {
00232    // return the angle w.r.t. another 3-vector
00233    Double_t ptot2 = Mag2()*q.Mag2();
00234    if(ptot2 <= 0) {
00235       return 0.0;
00236    } else {
00237       Double_t arg = Dot(q)/TMath::Sqrt(ptot2);
00238       if(arg >  1.0) arg =  1.0;
00239       if(arg < -1.0) arg = -1.0;
00240       return TMath::ACos(arg);
00241    }
00242 }
00243 
00244 //______________________________________________________________________________
00245 Double_t TVector3::Mag() const 
00246 { 
00247    // return the magnitude (rho in spherical coordinate system)
00248    
00249    return TMath::Sqrt(Mag2()); 
00250 }
00251 
00252 //______________________________________________________________________________
00253 Double_t TVector3::Perp() const 
00254 { 
00255    //return the transverse component  (R in cylindrical coordinate system)
00256 
00257    return TMath::Sqrt(Perp2()); 
00258 }
00259 
00260 
00261 //______________________________________________________________________________
00262 Double_t TVector3::Perp(const TVector3 & p) const
00263 { 
00264    //return the transverse component (R in cylindrical coordinate system)
00265 
00266    return TMath::Sqrt(Perp2(p)); 
00267 }
00268 
00269 //______________________________________________________________________________
00270 Double_t TVector3::Phi() const 
00271 {
00272    //return the  azimuth angle. returns phi from -pi to pi
00273    return fX == 0.0 && fY == 0.0 ? 0.0 : TMath::ATan2(fY,fX);
00274 }
00275 
00276 //______________________________________________________________________________
00277 Double_t TVector3::Theta() const 
00278 {
00279    //return the polar angle
00280    return fX == 0.0 && fY == 0.0 && fZ == 0.0 ? 0.0 : TMath::ATan2(Perp(),fZ);
00281 }
00282 
00283 //______________________________________________________________________________
00284 TVector3 TVector3::Unit() const 
00285 {
00286    // return unit vector parallel to this.
00287    Double_t  tot = Mag2();
00288    TVector3 p(fX,fY,fZ);
00289    return tot > 0.0 ? p *= (1.0/TMath::Sqrt(tot)) : p;
00290 }
00291 
00292 //______________________________________________________________________________
00293 void TVector3::RotateX(Double_t angle) {
00294    //rotate vector around X
00295    Double_t s = TMath::Sin(angle);
00296    Double_t c = TMath::Cos(angle);
00297    Double_t yy = fY;
00298    fY = c*yy - s*fZ;
00299    fZ = s*yy + c*fZ;
00300 }
00301 
00302 //______________________________________________________________________________
00303 void TVector3::RotateY(Double_t angle) {
00304    //rotate vector around Y
00305    Double_t s = TMath::Sin(angle);
00306    Double_t c = TMath::Cos(angle);
00307    Double_t zz = fZ;
00308    fZ = c*zz - s*fX;
00309    fX = s*zz + c*fX;
00310 }
00311 
00312 //______________________________________________________________________________
00313 void TVector3::RotateZ(Double_t angle) {
00314    //rotate vector around Z
00315    Double_t s = TMath::Sin(angle);
00316    Double_t c = TMath::Cos(angle);
00317    Double_t xx = fX;
00318    fX = c*xx - s*fY;
00319    fY = s*xx + c*fY;
00320 }
00321 
00322 //______________________________________________________________________________
00323 void TVector3::Rotate(Double_t angle, const TVector3 & axis){
00324    //rotate vector
00325    TRotation trans;
00326    trans.Rotate(angle, axis);
00327    operator*=(trans);
00328 }
00329 
00330 //______________________________________________________________________________
00331 void TVector3::RotateUz(const TVector3& NewUzVector) {
00332    // NewUzVector must be normalized !
00333 
00334    Double_t u1 = NewUzVector.fX;
00335    Double_t u2 = NewUzVector.fY;
00336    Double_t u3 = NewUzVector.fZ;
00337    Double_t up = u1*u1 + u2*u2;
00338 
00339    if (up) {
00340       up = TMath::Sqrt(up);
00341       Double_t px = fX,  py = fY,  pz = fZ;
00342       fX = (u1*u3*px - u2*py + u1*up*pz)/up;
00343       fY = (u2*u3*px + u1*py + u2*up*pz)/up;
00344       fZ = (u3*u3*px -    px + u3*up*pz)/up;
00345    } else if (u3 < 0.) { fX = -fX; fZ = -fZ; }      // phi=0  teta=pi
00346    else {};
00347 }
00348 
00349 //______________________________________________________________________________
00350 Double_t TVector3::PseudoRapidity() const {
00351    //Double_t m = Mag();
00352    //return 0.5*log( (m+fZ)/(m-fZ) );
00353    // guard against Pt=0
00354    double cosTheta = CosTheta();
00355    if (cosTheta*cosTheta < 1) return -0.5* TMath::Log( (1.0-cosTheta)/(1.0+cosTheta) );
00356    Warning("PseudoRapidity","transvers momentum = 0! return +/- 10e10");
00357    if (fZ > 0) return 10e10;
00358    else        return -10e10;
00359 }
00360 
00361 //______________________________________________________________________________
00362 void TVector3::SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi) {
00363    //set Pt, Eta and Phi
00364    Double_t apt = TMath::Abs(pt);
00365    SetXYZ(apt*TMath::Cos(phi), apt*TMath::Sin(phi), apt/TMath::Tan(2.0*TMath::ATan(TMath::Exp(-eta))) );
00366 }
00367 
00368 //______________________________________________________________________________
00369 void TVector3::SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi) {
00370    //set Pt, Theta and Phi
00371    fX = pt * TMath::Cos(phi);
00372    fY = pt * TMath::Sin(phi); 
00373    Double_t tanTheta = TMath::Tan(theta);
00374    fZ = tanTheta ? pt / tanTheta : 0;
00375 }
00376 
00377 //______________________________________________________________________________
00378 void TVector3::SetTheta(Double_t th) 
00379 {
00380    // Set theta keeping mag and phi constant (BaBar).
00381    Double_t ma   = Mag();
00382    Double_t ph   = Phi();
00383    SetX(ma*TMath::Sin(th)*TMath::Cos(ph));
00384    SetY(ma*TMath::Sin(th)*TMath::Sin(ph));
00385    SetZ(ma*TMath::Cos(th));
00386 }
00387 
00388 //______________________________________________________________________________
00389 void TVector3::SetPhi(Double_t ph) 
00390 {
00391    // Set phi keeping mag and theta constant (BaBar).
00392    Double_t xy   = Perp();
00393    SetX(xy*TMath::Cos(ph));
00394    SetY(xy*TMath::Sin(ph));
00395 }
00396 
00397 //______________________________________________________________________________
00398 Double_t TVector3::DeltaR(const TVector3 & v) const 
00399 {
00400    //return deltaR with respect to v
00401    Double_t deta = Eta()-v.Eta();
00402    Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
00403    return TMath::Sqrt( deta*deta+dphi*dphi );
00404 }
00405 
00406 //______________________________________________________________________________
00407 void TVector3::SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi) 
00408 {
00409    //setter with mag, theta, phi
00410    Double_t amag = TMath::Abs(mag);
00411    fX = amag * TMath::Sin(theta) * TMath::Cos(phi);
00412    fY = amag * TMath::Sin(theta) * TMath::Sin(phi);
00413    fZ = amag * TMath::Cos(theta);
00414 }
00415 
00416 //______________________________________________________________________________
00417 void TVector3::Streamer(TBuffer &R__b)
00418 {
00419    // Stream an object of class TVector3.
00420 
00421    if (R__b.IsReading()) {
00422       UInt_t R__s, R__c;
00423       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
00424       if (R__v > 2) {
00425          R__b.ReadClassBuffer(TVector3::Class(), this, R__v, R__s, R__c);
00426          return;
00427       }
00428       //====process old versions before automatic schema evolution
00429       if (R__v < 2) TObject::Streamer(R__b);
00430       R__b >> fX;
00431       R__b >> fY;
00432       R__b >> fZ;
00433       R__b.CheckByteCount(R__s, R__c, TVector3::IsA());
00434       //====end of old versions
00435       
00436    } else {
00437       R__b.WriteClassBuffer(TVector3::Class(),this);
00438    }
00439 }
00440 
00441 TVector3 operator + (const TVector3 & a, const TVector3 & b) {
00442    return TVector3(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
00443 }
00444 
00445 TVector3 operator - (const TVector3 & a, const TVector3 & b) {
00446    return TVector3(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
00447 }
00448 
00449 TVector3 operator * (const TVector3 & p, Double_t a) {
00450    return TVector3(a*p.X(), a*p.Y(), a*p.Z());
00451 }
00452 
00453 TVector3 operator * (Double_t a, const TVector3 & p) {
00454    return TVector3(a*p.X(), a*p.Y(), a*p.Z());
00455 }
00456 
00457 Double_t operator * (const TVector3 & a, const TVector3 & b) {
00458    return a.Dot(b);
00459 }
00460 
00461 TVector3 operator * (const TMatrix & m, const TVector3 & v ) {
00462    return TVector3( m(0,0)*v.X()+m(0,1)*v.Y()+m(0,2)*v.Z(),
00463                     m(1,0)*v.X()+m(1,1)*v.Y()+m(1,2)*v.Z(),
00464                     m(2,0)*v.X()+m(2,1)*v.Y()+m(2,2)*v.Z());
00465 }
00466 
00467 
00468 //const TVector3 kXHat(1.0, 0.0, 0.0);
00469 //const TVector3 kYHat(0.0, 1.0, 0.0);
00470 //const TVector3 kZHat(0.0, 0.0, 1.0);
00471 
00472 void TVector3::Print(Option_t*)const
00473 {
00474    //print vector parameters
00475    Printf("%s %s (x,y,z)=(%f,%f,%f) (rho,theta,phi)=(%f,%f,%f)",GetName(),GetTitle(),X(),Y(),Z(),
00476                                           Mag(),Theta()*TMath::RadToDeg(),Phi()*TMath::RadToDeg());
00477 }

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