Main Page | Modules | Namespace List | Class Hierarchy | Class List | File List | Namespace Members | Class Members | File Members

TMBVector3.hpp

Go to the documentation of this file.
00001 #ifndef INCLUDE_TMBVECTOR3
00002 #define INCLUDE_TMBVECTOR3
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 TVector3</H2>
00022 <TT>TVector3</TT> is a general three vector class, which can be used for
00023 the description of different vectors in 3D.
00024 <H3>
00025 Declaration / Access to the components</H3>
00026 <TT>TVector3</TT> has been implemented as a vector of three <TT>Double_t</TT>
00027 variables, representing the cartesian coordinates. By default all components
00028 are initialized to zero:
00029 
00030 <P><TT>&nbsp; TVector3 v1;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; //
00031 v1 = (0,0,0)</TT>
00032 <BR><TT>&nbsp; TVector3 v2(1);&nbsp;&nbsp;&nbsp;&nbsp; // v2 = (1,0,0)</TT>
00033 <BR><TT>&nbsp; TVector3 v3(1,2,3); // v3 = (1,2,3)</TT>
00034 <BR><TT>&nbsp; TVector3 v4(v2);&nbsp;&nbsp;&nbsp; // v4 = v2</TT>
00035 
00036 <P>It is also possible (but not recommended) to initialize a <TT>TVector3</TT>
00037 with a <TT>Double_t</TT> or <TT>Float_t</TT> C array.
00038 
00039 <P>You can get the basic components either by name or by index using <TT>operator()</TT>:
00040 
00041 <P><TT>&nbsp; xx = v1.X();&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp; xx =
00042 v1(0);</TT>
00043 <BR><TT>&nbsp; yy = v1.Y();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
00044 yy = v1(1);</TT>
00045 <BR><TT>&nbsp; zz = v1.Z();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
00046 zz = v1(2);</TT>
00047 
00048 <P>The memberfunctions <TT>SetX()</TT>, <TT>SetY()</TT>, <TT>SetZ()</TT>
00049 and<TT> SetXYZ()</TT> allow to set the components:
00050 
00051 <P><TT>&nbsp; v1.SetX(1.); v1.SetY(2.); v1.SetZ(3.);</TT>
00052 <BR><TT>&nbsp; v1.SetXYZ(1.,2.,3.);</TT>
00053 <BR>&nbsp;
00054 <H3>
00055 Noncartesian coordinates</H3>
00056 To get information on the <TT>TVector3</TT> in spherical (rho,phi,theta)
00057 or cylindrical (z,r,theta) coordinates, the
00058 <BR>the member functions <TT>Mag3()</TT> (=magnitude=rho in spherical coordinates),
00059 <TT>Mag32()</TT>, <TT>Theta()</TT>, <TT>CosTheta()</TT>, <TT>Phi()</TT>,
00060 <TT>Perp()</TT> (the transverse component=r in cylindrical coordinates),
00061 <TT>Perp2()</TT> can be used:
00062 
00063 <P><TT>&nbsp; Double_t m&nbsp; = v.Mag3();&nbsp;&nbsp;&nbsp; // get magnitude
00064 (=rho=Sqrt(x*x+y*y+z*z)))</TT>
00065 <BR><TT>&nbsp; Double_t m2 = v.Mag32();&nbsp;&nbsp; // get magnitude squared</TT>
00066 <BR><TT>&nbsp; Double_t t&nbsp; = v.Theta();&nbsp; // get polar angle</TT>
00067 <BR><TT>&nbsp; Double_t ct = v.CosTheta();// get cos of theta</TT>
00068 <BR><TT>&nbsp; Double_t p&nbsp; = v.Phi();&nbsp;&nbsp;&nbsp; // get azimuth angle</TT>
00069 <BR><TT>&nbsp; Double_t pp = v.Perp();&nbsp;&nbsp; // get transverse component</TT>
00070 <BR><TT>&nbsp; Double_t pp2= v.Perp2();&nbsp; // get transvers component
00071 squared</TT>
00072 
00073 <P>It is also possible to get the transverse component with respect to
00074 another vector:
00075 
00076 <P><TT>&nbsp; Double_t ppv1 = v.Perp(v1);</TT>
00077 <BR><TT>&nbsp; Double_t pp2v1 = v.Perp2(v1);</TT>
00078 
00079 <P>The pseudorapiditiy ( eta=-ln (tan (phi/2)) ) can be get by <TT>Eta()</TT>
00080 or <TT>PseudoRapidity()</TT>:
00081 <BR>&nbsp;
00082 <BR><TT>&nbsp; Double_t eta = v.PseudoRapidity();</TT>
00083 
00084 <P>There are set functions to change one of the noncartesian coordinates:
00085 
00086 <P><TT>&nbsp; v.SetTheta(.5); // keeping rho and phi</TT>
00087 <BR><TT>&nbsp; v.SetPhi(.8);&nbsp;&nbsp; // keeping rho and theta</TT>
00088 <BR><TT>&nbsp; v.SetMag3(10.);&nbsp; // keeping theta and phi</TT>
00089 <BR><TT>&nbsp; v.SetPerp(3.);&nbsp; // keeping z and phi</TT>
00090 <BR>&nbsp;
00091 <H3>
00092 Arithmetic / Comparison</H3>
00093 The <TT>TVector3</TT> class provides the operators to add, subtract, scale and compare
00094 vectors:
00095 
00096 <P><TT>&nbsp; v3&nbsp; = -v1;</TT>
00097 <BR><TT>&nbsp; v1&nbsp; = v2+v3;</TT>
00098 <BR><TT>&nbsp; v1 += v3;</TT>
00099 <BR><TT>&nbsp; v1&nbsp; = v1 - v3</TT>
00100 <BR><TT>&nbsp; v1 -= v3;</TT>
00101 <BR><TT>&nbsp; v1 *= 10;</TT>
00102 <BR><TT>&nbsp; v1&nbsp; = 5*v2;</TT>
00103 
00104 <P><TT>&nbsp; if(v1==v2) {...}</TT>
00105 <BR><TT>&nbsp; if(v1!=v2) {...}</TT>
00106 <BR>&nbsp;
00107 <H3>
00108 Related Vectors</H3>
00109 <TT>&nbsp; v2 = v1.Unit();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // get unit
00110 vector parallel to v1</TT>
00111 <BR><TT>&nbsp; v2 = v1.Orthogonal(); // get vector orthogonal to v1</TT>
00112 <H3>
00113 Scalar and vector products</H3>
00114 <TT>&nbsp; s = v1.Dot(v2);&nbsp;&nbsp; // scalar product</TT>
00115 <BR><TT>&nbsp; s = v1 * v2;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // scalar product</TT>
00116 <BR><TT>&nbsp; v = v1.Cross(v2); // vector product</TT>
00117 <H3>
00118 &nbsp;Angle between two vectors</H3>
00119 <TT>&nbsp; Double_t a = v1.Angle(v2);</TT>
00120 <H3>
00121 Rotations</H3>
00122 
00123 <H5>
00124 Rotation around axes</H5>
00125 <TT>&nbsp; v.RotateX(.5);</TT>
00126 <BR><TT>&nbsp; v.RotateY(TMath::Pi());</TT>
00127 <BR><TT>&nbsp; v.RotateZ(angle);</TT>
00128 <H5>
00129 Rotation around a vector</H5>
00130 <TT>&nbsp; v1.Rotate(TMath::Pi()/4, v2); // rotation around v2</TT>
00131 <H5>
00132 Rotation by TRotation</H5>
00133 <TT>TVector3</TT> objects can be rotated by objects of the <TT>TRotation</TT>
00134 class using the <TT>Transform()</TT> member functions,
00135 <BR>the <TT>operator *=</TT> or the <TT>operator *</TT> of the TRotation
00136 class:
00137 
00138 <P><TT>&nbsp; TRotation m;</TT>
00139 <BR><TT>&nbsp; ...</TT>
00140 <BR><TT>&nbsp; v1.transform(m);</TT>
00141 <BR><TT>&nbsp; v1 = m*v1;</TT>
00142 <BR><TT>&nbsp; v1 *= m; // Attention v1 = m*v1</TT>
00143 <H5>
00144 Transformation from rotated frame</H5>
00145 <TT>&nbsp; TVector3 direction = v.Unit()</TT>
00146 <BR><TT>&nbsp; v1.RotateUz(direction); // direction must be TVector3 of
00147 unit length</TT>
00148 
00149 <P>transforms v1 from the rotated frame (z' parallel to direction, x' in
00150 the theta plane and y' in the xy plane as well as perpendicular to the
00151 theta plane) to the (x,y,z) frame.'
00152 
00153 <!--*/
00154 // -->END_HTML
00155 //
00156 
00157 // See copyright statements at end of file.
00158 
00159 #include "TMath.h"
00160 #include "TError.h"
00161 #include "TVector2.h"
00162 #include "TMatrix.h"
00163 #include "TRotation.h"
00164 
00165 class TMBVector3 : public TObject {
00166 public:
00167 
00168    TMBVector3(Double_t x = 0.0, Double_t y = 0.0, Double_t z = 0.0):
00169        fX(x), fY(y), fZ(z) {  
00170    }
00171 
00172    TMBVector3(const Double_t *a): fX(a[0]), fY(a[1]), fZ(a[2]) 
00173    { 
00174    }
00175 
00176    TMBVector3(const Float_t *a): fX(a[0]), fY(a[1]), fZ(a[2]) 
00177    { 
00178    }
00179 
00180    TMBVector3(const TMBVector3 &p): TObject(p),
00181       fX(p.fX), fY(p.fY), fZ(p.fZ)
00182    {  
00183    }
00184 
00185    TMBVector3(const TVector3& v): fX(v.X()), fY(v.Y()), fZ(v.Z())
00186    {
00188    }
00189 
00190   virtual ~TMBVector3() {}
00191   // Destructor
00192 
00193   Double_t operator () (int i) const {
00195      switch(i) {
00196      case 0: return fX;
00197      case 1: return fY;
00198      case 2: return fZ;
00199      default: Error("operator()(i)", "bad index (%d) returning &fX",i);
00200      }
00201      return 0.;
00202   }
00203      
00204   inline Double_t operator [] (int i) const { 
00205      return operator()(i); }
00206 
00207 
00208   Double_t & operator () (int i) {
00210 
00211      switch(i) {
00212      case 0: return fX;
00213      case 1: return fY;
00214      case 2: return fZ;
00215      default: Error("operator()(i)", "bad index (%d) returning &fX",i);
00216      }
00217      return fX;
00218   }
00219   inline Double_t & operator [] (int i) { 
00220      return operator()(i); }
00221 
00222    Double_t x() const { 
00223       return fX; }
00224    Double_t X() const { 
00225       return x(); }
00226    Double_t Px() const { 
00227       return x(); }
00228 
00229    Double_t y() const { 
00230       return fY; }
00231    Double_t Y() const { 
00232       return y(); }
00233    Double_t Py() const { 
00234       return y(); }
00235 
00236    Double_t z() const { 
00237       return fZ; }
00238    Double_t Z() const { 
00239       return z(); }
00240    Double_t Pz() const { 
00241       return z(); }
00242 
00243    void SetX(Double_t a) { 
00244       fX=a; }
00245    void SetPx(Double_t a) { 
00246       SetX(a); }
00247    void SetY(Double_t a) { 
00248       fY=a; }
00249    void SetPy(Double_t a) { 
00250       SetY(a); }
00251    void SetZ(Double_t a) { 
00252       fZ=a; }
00253    void SetPz(Double_t a) { 
00254       SetZ(a); }
00255    inline void SetXYZ(Double_t x, Double_t y, Double_t z) {
00257        fX=x; fY=y; fZ=z; }
00258     
00259    void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi) {
00261       Double_t apt = TMath::Abs(pt);
00262       Double_t d=TMath::Tan(2.0*TMath::ATan(TMath::Exp(-eta)));
00263       if (d==0.) 
00264          Warning("SetPtEtaPhi", "tan(2*atan(exp(-eta)))==0.! Setting z=0.");
00265       SetXYZ(apt*TMath::Cos(phi), apt*TMath::Sin(phi), d==0.?0.:apt/d );
00266    }
00267    void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi) {
00269       Double_t apt = TMath::Abs(pt);
00270       fX = apt * TMath::Cos(phi);
00271       fY = apt * TMath::Sin(phi); 
00272       Double_t tanTheta = TMath::Tan(theta);
00273       if (tanTheta==0.) 
00274          Warning("SetPtThetaPhi", "tan(theta)==0.! Setting z=0.");
00275       fZ = tanTheta ? apt / tanTheta : 0;
00276    }
00277    
00278    inline void GetXYZ(Double_t *carray) const {
00280       carray[0]=fX; carray[1]=fY; carray[2]=fZ; }
00281    inline void GetXYZ(Float_t *carray) const {
00283       carray[0]=fX; carray[1]=fY; carray[2]=fZ; }
00284 
00285    Double_t Eta() const { 
00286       return PseudoRapidity(); }
00287    Double_t Theta() const { 
00288       return fX==.0 && fY==.0 && fZ==.0 ? .0 : TMath::ATan2(Perp(),fZ); }
00289    Double_t CosTheta() const { 
00290       Double_t ptot = Mag3();
00291       return ptot == 0.0 ? 1.0 : fZ/ptot;
00292    }
00293    Double_t Phi() const { 
00294       return (fX==.0 && fY==.0 ? .0 : TVector2::Phi_0_2pi(TMath::ATan2(fY,fX))); }
00295 
00296    Double_t Rho() const { 
00297       return Mag3(); }
00298 
00299    virtual Double_t Mag32() const {
00301       return fX*fX + fY*fY + fZ*fZ; }
00302    virtual Double_t Mag3() const {
00304       return TMath::Sqrt(Mag32()); }
00305    Double_t P()  const { 
00306       return Mag3(); }
00307 
00308    Double_t DeltaPhi(const TMBVector3 &v) const { 
00310       return TVector2::Phi_mpi_pi(Phi()-v.Phi()); }
00311 
00312    Double_t DeltaR(const TMBVector3 &v) const {
00314       Double_t deta = Eta()-v.Eta();
00315       Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
00316       return TMath::Sqrt( deta*deta+dphi*dphi );
00317    }
00318    Double_t DrEtaPhi(const TMBVector3 &lv) const {
00320       return DeltaR(lv); }
00321 
00322     //   TVector2 EtaPhiVector() { /// return TVector2(eta,phi)
00323     //      return TVector2 (Eta(),Phi()); }
00324 
00325    Double_t Angle(const TVector3 & q) const { 
00326       Double_t ptot2 = Mag32()*q.Mag2();
00327       if(ptot2 <= 0.) return .0;
00328       Double_t arg = Dot(q)/TMath::Sqrt(ptot2);
00329       if(arg >  1.0) arg =  1.0;
00330       if(arg < -1.0) arg = -1.0;
00331       return TMath::ACos(arg);
00332    }
00333 
00334    void SetPhi(Double_t ph) { 
00335       Double_t xy   = Perp();
00336       SetX(xy*TMath::Cos(ph));
00337       SetY(xy*TMath::Sin(ph));
00338    }
00339 
00340   inline void SetTheta(Double_t th) { 
00341      Double_t ma   = Mag3();
00342      Double_t ph   = Phi();
00343      SetX(ma*TMath::Sin(th)*TMath::Cos(ph));
00344      SetY(ma*TMath::Sin(th)*TMath::Sin(ph));
00345      SetZ(ma*TMath::Cos(th));
00346   }
00347 
00348   inline void SetMag3(Double_t ma) { 
00349      Double_t factor = Mag3();
00350      if (factor == 0) {
00351         Warning("SetMag3","zero vector can't be stretched");
00352      }else{
00353         factor = ma/factor;
00354         SetX(fX*factor);
00355         SetY(fY*factor);
00356         SetZ(fZ*factor);
00357      }
00358   }
00359   inline Double_t Perp2() const { 
00361      return fX*fX + fY*fY; }
00362   inline Double_t Pt() const { 
00364      return Perp(); }
00365   inline Double_t Perp() const { 
00367      return TMath::Sqrt(Perp2()); }
00368      
00369   inline void SetPerp(Double_t r) { 
00371      Double_t p = Perp();
00372      if (p != 0.0) {
00373         fX *= r/p;
00374         fY *= r/p;
00375      }
00376   }
00377 
00378   inline Double_t Perp2(const TMBVector3 &p) const { 
00380      Double_t tot = p.Mag32();
00381      Double_t ss  = Dot(p);
00382      return tot > 0.0 ? Mag32()-ss*ss/tot : Mag32();
00383   }
00384 
00385   inline Double_t Pt(const TMBVector3 &p) const {
00387      return Perp(p); }
00388   inline Double_t Perp(const TMBVector3 &p) const {
00390      return TMath::Sqrt(Perp2(p)); }
00391 
00392   inline void SetMag3ThetaPhi(Double_t mag, Double_t theta, Double_t phi) {
00394      Double_t amag = TMath::Abs(mag);
00395      fX = amag * TMath::Sin(theta) * TMath::Cos(phi);
00396      fY = amag * TMath::Sin(theta) * TMath::Sin(phi);
00397      fZ = amag * TMath::Cos(theta);
00398   }
00399 
00400   inline TMBVector3 & operator = (const TMBVector3 &p) {
00402      fX = p.fX; fY = p.fY; fZ = p.fZ;
00403      return *this;
00404   }
00405 
00406   inline Bool_t operator == (const TMBVector3 &v) const {
00408      return (v.fX==fX && v.fY==fY && v.fZ==fZ); }
00409   inline Bool_t operator != (const TMBVector3 &v) const {
00411      return (v.fX!=fX || v.fY!=fY || v.fZ!=fZ); }
00412 
00415   Bool_t is_equal (const TMBVector3 &v) const;
00416 
00417   inline TMBVector3 & operator += (const TMBVector3 &p) { 
00418      fX += p.fX; fY += p.fY; fZ += p.fZ;
00419      return *this;
00420   }
00421 
00422   inline TMBVector3 & operator -= (const TMBVector3 &p) { 
00423      fX -= p.fX; fY -= p.fY; fZ -= p.fZ;
00424      return *this;
00425   }
00426 
00427   inline TMBVector3 operator - () const { 
00428      return TMBVector3(-fX, -fY, -fZ); }
00429 
00430   inline TMBVector3 & operator *= (Double_t a) { 
00432      fX *= a; fY *= a; fZ *= a;
00433      return *this;
00434   }
00435 
00436   inline TMBVector3 Unit() const { 
00437      Double_t  tot = Mag32();
00438      TMBVector3 p(fX,fY,fZ);
00439      return tot>.0 ? p *= (1.0/TMath::Sqrt(tot)) : p;
00440   }
00441 
00442   inline TMBVector3 Orthogonal() const { 
00443      Double_t x = fX < 0.0 ? -fX : fX;
00444      Double_t y = fY < 0.0 ? -fY : fY;
00445      Double_t z = fZ < 0.0 ? -fZ : fZ;
00446      if (x < y)
00447         return x<z ? TMBVector3(0,fZ,-fY) : TMBVector3(fY,-fX,0);
00448      return y<z ? TMBVector3(-fZ,0,fX) : TMBVector3(fY,-fX,0);
00449   }
00450 
00451   Double_t Dot(const TMBVector3 &p) const { 
00452      return fX*p.fX + fY*p.fY + fZ*p.fZ; }
00453 
00454   inline TMBVector3 Cross(const TMBVector3 &p) const { 
00455      return TMBVector3(fY*p.fZ-p.fY*fZ, fZ*p.fX-p.fZ*fX, fX*p.fY-p.fX*fY);
00456   }
00457 
00458   Double_t PseudoRapidity() const {
00460      double cosTheta = CosTheta();
00461      if (cosTheta*cosTheta < 1) 
00462         return (-0.5* TMath::Log( (1.0-cosTheta)/(1.0+cosTheta) ));
00463      Warning("PseudoRapidity","transvers momentum = 0! return +/- 10e10");
00464      if (fZ > 0) return (10e10);
00465      else        return (-10e10);
00466   }
00467   void RotateX(Double_t angle) {
00469      Double_t s = TMath::Sin(angle);
00470      Double_t c = TMath::Cos(angle);
00471      Double_t yy = fY;
00472      fY = c*yy - s*fZ;
00473      fZ = s*yy + c*fZ;
00474   }
00475 
00476   void RotateY(Double_t angle) {
00478      Double_t s = TMath::Sin(angle);
00479      Double_t c = TMath::Cos(angle);
00480      Double_t zz = fZ;
00481      fZ = c*zz - s*fX;
00482      fX = s*zz + c*fX;
00483   }
00484 
00485   void RotateZ(Double_t angle) {
00487      Double_t s = TMath::Sin(angle);
00488      Double_t c = TMath::Cos(angle);
00489      Double_t xx = fX;
00490      fX = c*xx - s*fY;
00491      fY = s*xx + c*fY;
00492   }
00493 
00494   void RotateUz(const TMBVector3& NewUzVector) {
00496      Double_t u1 = NewUzVector.fX;
00497      Double_t u2 = NewUzVector.fY;
00498      Double_t u3 = NewUzVector.fZ;
00499      Double_t up = u1*u1 + u2*u2;
00500 
00501      if (up>0.) {
00502         up = TMath::Sqrt(up);
00503         Double_t px = fX,  py = fY,  pz = fZ;
00504         fX = (u1*u3*px - u2*py + u1*up*pz)/up;
00505         fY = (u2*u3*px + u1*py + u2*up*pz)/up;
00506         fZ = (u3*u3*px -    px + u3*up*pz)/up;
00507      }
00508      else if (u3 < 0.) { fX = -fX; fZ = -fZ; }      // phi=0  teta=pi
00509      else {};
00510   }
00511 
00512   void Rotate(Double_t angle, const TMBVector3 &axis) {
00514      TRotation trans;
00515      trans.Rotate(angle, (TVector3)axis);
00516      operator*=(trans);
00517   }
00518 
00519   TMBVector3 & operator *= (const TRotation &m) {
00521      return *this = m * (*this);
00522   }
00523   TMBVector3 & Transform(const TRotation &m) {
00525      return *this = m * (*this);
00526   }
00527 
00528    TVector2 XYvector() const {
00530      return TVector2(fX,fY); }
00531 
00532    operator TVector3() const { 
00534       return TVector3(X(),Y(),Z()); }
00535 
00536    void Print(Option_t* option="") const {
00538       Printf("%s %s (x,y,z)=(%f,%f,%f) (rho,theta,phi)=(%f,%f,%f)",
00539              GetName(),GetTitle(),X(),Y(),Z(),
00540              Mag3(),Theta()*TMath::RadToDeg(),Phi()*TMath::RadToDeg());
00541    }
00542 
00543    // Helper to test if two FP numbers are equivalent within machine precision.
00544    static bool is_equal (double x1, double x2);
00545 
00546 private:
00547    Double32_t fX; 
00548    Double32_t fY; 
00549    Double32_t fZ; 
00550 
00551    ClassDef(TMBVector3,1) // A 3D physics vector
00552 
00553 };
00554 
00555 inline
00556 TMBVector3 operator + (const TMBVector3 &a, const TMBVector3 &b) {
00558   return TVector3(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z()); }
00559 
00560 inline
00561 TMBVector3 operator - (const TMBVector3 &a, const TMBVector3 &b) {
00563    return TVector3(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z()); }
00564 
00565 inline
00566 Double_t operator * (const TMBVector3 &a, const TMBVector3 &b) {
00568    return a.Dot(b); }
00569 
00570 inline
00571 TMBVector3 operator * (const TMBVector3 &p, Double_t a) {
00573    return TVector3(a*p.X(), a*p.Y(), a*p.Z()); }
00574 
00575 inline
00576 TMBVector3 operator * (Double_t a, const TMBVector3 &p) {
00578    return TVector3(a*p.X(), a*p.Y(), a*p.Z()); }
00579 
00580 inline
00581 TMBVector3 operator * (const TMatrix &m, const TMBVector3 &v) {
00583   return TVector3( m(0,0)*v.X()+m(0,1)*v.Y()+m(0,2)*v.Z(),
00584                    m(1,0)*v.X()+m(1,1)*v.Y()+m(1,2)*v.Z(),
00585                    m(2,0)*v.X()+m(2,1)*v.Y()+m(2,2)*v.Z());
00586 }
00587 
00588 /*************************************************************************
00589  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00590  * All rights reserved.                                                  *
00591  *                                                                       *
00592  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00593  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00594  *************************************************************************/
00595 
00596 // Author: Pasha Murat, Peter Malzacher   12/02/99
00597 
00598 #endif // INCLUDE_TMBVECTOR3

Generated on Tue Mar 28 10:13:05 2006 for CAF by doxygen 1.3.4