CLHEP VERSION Reference Documentation
   
CLHEP Home Page     CLHEP Documentation     CLHEP Bug Reports

ThreeVector.h

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 // CLASSDOC OFF
00003 // $Id: ThreeVector.h,v 1.4 2010/06/16 17:15:57 garren Exp $
00004 // ---------------------------------------------------------------------------
00005 // CLASSDOC ON
00006 //
00007 // This file is a part of the CLHEP - a Class Library for High Energy Physics.
00008 //
00009 // Hep3Vector is a general 3-vector class defining vectors in three
00010 // dimension using double components. Rotations of these vectors are
00011 // performed by multiplying with an object of the HepRotation class.
00012 //
00013 // .SS See Also
00014 // LorentzVector.h, Rotation.h, LorentzRotation.h 
00015 //
00016 // .SS Authors
00017 // Leif Lonnblad and Anders Nilsson; Modified by Evgueni Tcherniaev;
00018 // ZOOM additions by Mark Fischler
00019 //
00020 
00021 #ifndef HEP_THREEVECTOR_H
00022 #define HEP_THREEVECTOR_H
00023 
00024 #ifdef GNUPRAGMA
00025 #pragma interface
00026 #endif
00027 
00028 #include <iostream>
00029 #include "CLHEP/Vector/defs.h" 
00030 
00031 namespace CLHEP {
00032 
00033 class HepRotation;
00034 class HepEulerAngles;
00035 class HepAxisAngle;
00036 
00041 class Hep3Vector {
00042 
00043 public:
00044 
00045 // Basic properties and operations on 3-vectors:  
00046 
00047   enum { X=0, Y=1, Z=2, NUM_COORDINATES=3, SIZE=NUM_COORDINATES };
00048   // Safe indexing of the coordinates when using with matrices, arrays, etc.
00049   // (BaBar)
00050 
00051   Hep3Vector();
00052   explicit Hep3Vector(double x);
00053   Hep3Vector(double x, double y);
00054   Hep3Vector(double x, double y, double z);
00055   // The constructor.  
00056 
00057   inline Hep3Vector(const Hep3Vector &);
00058   // The copy constructor.
00059 
00060   inline ~Hep3Vector();
00061   // The destructor.  Not virtual - inheritance from this class is dangerous.
00062 
00063   double operator () (int) const;
00064   // Get components by index -- 0-based (Geant4) 
00065 
00066   inline double operator [] (int) const;
00067   // Get components by index -- 0-based (Geant4) 
00068 
00069   double & operator () (int);
00070   // Set components by index.  0-based.
00071 
00072   inline double & operator [] (int);
00073   // Set components by index.  0-based.
00074 
00075   inline double x() const;
00076   inline double y() const;
00077   inline double z() const;
00078   // The components in cartesian coordinate system.  Same as getX() etc.
00079 
00080   inline void setX(double);
00081   inline void setY(double);
00082   inline void setZ(double);
00083   // Set the components in cartesian coordinate system.
00084 
00085   inline void set( double x, double y, double z); 
00086   // Set all three components in cartesian coordinate system.
00087 
00088   inline double phi() const;
00089   // The azimuth angle.
00090 
00091   inline double theta() const;
00092   // The polar angle.
00093 
00094   inline double cosTheta() const;
00095   // Cosine of the polar angle.
00096 
00097   inline double cos2Theta() const;
00098   // Cosine squared of the polar angle - faster than cosTheta(). (ZOOM)
00099 
00100   inline double mag2() const;
00101   // The magnitude squared (r^2 in spherical coordinate system).
00102 
00103   inline double mag() const;
00104   // The magnitude (r in spherical coordinate system).
00105 
00106   inline void setPhi(double);
00107   // Set phi keeping mag and theta constant (BaBar).
00108 
00109   inline void setTheta(double);
00110   // Set theta keeping mag and phi constant (BaBar).
00111 
00112          void setMag(double);
00113   // Set magnitude keeping theta and phi constant (BaBar).
00114 
00115   inline double perp2() const;
00116   // The transverse component squared (rho^2 in cylindrical coordinate system).
00117 
00118   inline double perp() const;
00119   // The transverse component (rho in cylindrical coordinate system).
00120 
00121   inline void setPerp(double);
00122   // Set the transverse component keeping phi and z constant.
00123 
00124   void setCylTheta(double);
00125   // Set theta while keeping transvers component and phi fixed 
00126 
00127   inline double perp2(const Hep3Vector &) const;
00128   // The transverse component w.r.t. given axis squared.
00129 
00130   inline double perp(const Hep3Vector &) const;
00131   // The transverse component w.r.t. given axis.
00132 
00133   inline Hep3Vector & operator = (const Hep3Vector &);
00134   // Assignment.
00135 
00136   inline bool operator == (const Hep3Vector &) const;
00137   inline bool operator != (const Hep3Vector &) const;
00138   // Comparisons (Geant4). 
00139 
00140   bool isNear (const Hep3Vector &, double epsilon=tolerance) const;
00141   // Check for equality within RELATIVE tolerance (default 2.2E-14). (ZOOM)
00142   // |v1 - v2|**2 <= epsilon**2 * |v1.dot(v2)| 
00143 
00144   double howNear(const Hep3Vector & v ) const;
00145   // sqrt ( |v1-v2|**2 / v1.dot(v2) ) with a maximum of 1.
00146   // If v1.dot(v2) is negative, will return 1.
00147 
00148   double deltaR(const Hep3Vector & v) const;
00149   // sqrt( pseudorapity_difference**2 + deltaPhi **2 )
00150 
00151   inline Hep3Vector & operator += (const Hep3Vector &);
00152   // Addition.
00153 
00154   inline Hep3Vector & operator -= (const Hep3Vector &);
00155   // Subtraction.
00156 
00157   inline Hep3Vector operator - () const;
00158   // Unary minus.
00159 
00160   inline Hep3Vector & operator *= (double);
00161   // Scaling with real numbers.
00162 
00163          Hep3Vector & operator /= (double);
00164   // Division by (non-zero) real number.
00165 
00166   inline Hep3Vector unit() const;
00167   // Vector parallel to this, but of length 1.
00168 
00169   inline Hep3Vector orthogonal() const;
00170   // Vector orthogonal to this (Geant4).
00171 
00172   inline double dot(const Hep3Vector &) const;
00173   // double product.
00174 
00175   inline Hep3Vector cross(const Hep3Vector &) const;
00176   // Cross product.
00177 
00178   double angle(const Hep3Vector &) const;
00179   // The angle w.r.t. another 3-vector.
00180 
00181   double pseudoRapidity() const;
00182   // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
00183 
00184   void setEta  ( double p );
00185   // Set pseudo-rapidity, keeping magnitude and phi fixed.  (ZOOM)
00186 
00187   void setCylEta  ( double p );
00188   // Set pseudo-rapidity, keeping transverse component and phi fixed.  (ZOOM)
00189 
00190   Hep3Vector & rotateX(double);
00191   // Rotates the Hep3Vector around the x-axis.
00192 
00193   Hep3Vector & rotateY(double);
00194   // Rotates the Hep3Vector around the y-axis.
00195 
00196   Hep3Vector & rotateZ(double);
00197   // Rotates the Hep3Vector around the z-axis.
00198 
00199   Hep3Vector & rotateUz(const Hep3Vector&);
00200   // Rotates reference frame from Uz to newUz (unit vector) (Geant4).
00201 
00202     Hep3Vector & rotate(double, const Hep3Vector &);
00203   // Rotates around the axis specified by another Hep3Vector.
00204   // (Uses methods of HepRotation, forcing linking in of Rotation.cc.)
00205 
00206   Hep3Vector & operator *= (const HepRotation &);
00207   Hep3Vector & transform(const HepRotation &);
00208   // Transformation with a Rotation matrix.
00209 
00210 
00211 // = = = = = = = = = = = = = = = = = = = = = = = =
00212 //
00213 // Esoteric properties and operations on 3-vectors:  
00214 //
00215 // 1 - Set vectors in various coordinate systems
00216 // 2 - Synonyms for accessing coordinates and properties
00217 // 3 - Comparisions (dictionary, near-ness, and geometric)
00218 // 4 - Intrinsic properties 
00219 // 5 - Properties releative to z axis and arbitrary directions
00220 // 6 - Polar and azimuthal angle decomposition and deltaPhi
00221 // 7 - Rotations 
00222 //
00223 // = = = = = = = = = = = = = = = = = = = = = = = =
00224 
00225 // 1 - Set vectors in various coordinate systems
00226 
00227   inline void setRThetaPhi  (double r, double theta, double phi);
00228   // Set in spherical coordinates:  Angles are measured in RADIANS
00229 
00230   inline void setREtaPhi  ( double r, double eta,  double phi );
00231   // Set in spherical coordinates, but specify peudorapidiy to determine theta.
00232 
00233   inline void setRhoPhiZ   (double rho, double phi, double z);
00234   // Set in cylindrical coordinates:  Phi angle is measured in RADIANS
00235 
00236   void setRhoPhiTheta ( double rho, double phi, double theta);
00237   // Set in cylindrical coordinates, but specify theta to determine z.
00238 
00239   void setRhoPhiEta ( double rho, double phi, double eta);
00240   // Set in cylindrical coordinates, but specify pseudorapidity to determine z.
00241 
00242 // 2 - Synonyms for accessing coordinates and properties
00243 
00244   inline double getX() const; 
00245   inline double getY() const;
00246   inline double getZ() const; 
00247   // x(), y(), and z()
00248 
00249   inline double getR    () const;
00250   inline double getTheta() const;
00251   inline double getPhi  () const;
00252   // mag(), theta(), and phi()
00253 
00254   inline double r       () const;
00255   // mag()
00256 
00257   inline double rho     () const;
00258   inline double getRho  () const;
00259   // perp()
00260 
00261   double eta     () const;
00262   double getEta  () const;
00263   // pseudoRapidity() 
00264 
00265   inline void setR ( double s );
00266   // setMag()
00267 
00268   inline void setRho ( double s );
00269   // setPerp()
00270 
00271 // 3 - Comparisions (dictionary, near-ness, and geometric)
00272 
00273   int compare (const Hep3Vector & v) const;
00274   bool operator > (const Hep3Vector & v) const;
00275   bool operator < (const Hep3Vector & v) const;
00276   bool operator>= (const Hep3Vector & v) const;
00277   bool operator<= (const Hep3Vector & v) const;
00278   // dictionary ordering according to z, then y, then x component
00279 
00280   inline double diff2 (const Hep3Vector & v) const;
00281   // |v1-v2|**2
00282 
00283   static double setTolerance (double tol);
00284   static inline double getTolerance ();
00285   // Set the tolerance used in isNear() for Hep3Vectors 
00286 
00287   bool isParallel (const Hep3Vector & v, double epsilon=tolerance) const;
00288   // Are the vectors parallel, within the given tolerance?
00289 
00290   bool isOrthogonal (const Hep3Vector & v, double epsilon=tolerance) const;
00291   // Are the vectors orthogonal, within the given tolerance?
00292 
00293   double howParallel   (const Hep3Vector & v) const;
00294   // | v1.cross(v2) / v1.dot(v2) |, to a maximum of 1.
00295 
00296   double howOrthogonal (const Hep3Vector & v) const;
00297   // | v1.dot(v2) / v1.cross(v2) |, to a maximum of 1.
00298 
00299   enum { ToleranceTicks = 100 };
00300 
00301 // 4 - Intrinsic properties 
00302 
00303   double beta    () const;
00304   // relativistic beta (considering v as a velocity vector with c=1)
00305   // Same as mag() but will object if >= 1
00306 
00307   double gamma() const;
00308   // relativistic gamma (considering v as a velocity vector with c=1)
00309 
00310   double coLinearRapidity() const;
00311   // inverse tanh (beta)
00312 
00313 // 5 - Properties relative to Z axis and to an arbitrary direction
00314 
00315           // Note that the non-esoteric CLHEP provides 
00316           // theta(), cosTheta(), cos2Theta, and angle(const Hep3Vector&)
00317 
00318   inline double angle() const;
00319   // angle against the Z axis -- synonym for theta()
00320 
00321   inline double theta(const Hep3Vector & v2) const;  
00322   // synonym for angle(v2)
00323 
00324   double cosTheta (const Hep3Vector & v2) const;
00325   double cos2Theta(const Hep3Vector & v2) const;
00326   // cos and cos^2 of the angle between two vectors
00327 
00328   inline Hep3Vector project () const;
00329          Hep3Vector project (const Hep3Vector & v2) const;
00330   // projection of a vector along a direction.  
00331 
00332   inline Hep3Vector perpPart() const;
00333   inline Hep3Vector perpPart (const Hep3Vector & v2) const;
00334   // vector minus its projection along a direction.
00335 
00336   double rapidity () const;
00337   // inverse tanh(v.z())
00338 
00339   double rapidity (const Hep3Vector & v2) const;
00340   // rapidity with respect to specified direction:  
00341   // inverse tanh (v.dot(u)) where u is a unit in the direction of v2
00342 
00343   double eta(const Hep3Vector & v2) const;
00344   // - ln tan of the angle beween the vector and the ref direction.
00345 
00346 // 6 - Polar and azimuthal angle decomposition and deltaPhi
00347 
00348   // Decomposition of an angle within reference defined by a direction:
00349 
00350   double polarAngle (const Hep3Vector & v2) const;
00351   // The reference direction is Z: the polarAngle is abs(v.theta()-v2.theta()).
00352 
00353   double deltaPhi (const Hep3Vector & v2) const;
00354   // v.phi()-v2.phi(), brought into the range (-PI,PI]
00355 
00356   double azimAngle  (const Hep3Vector & v2) const;
00357   // The reference direction is Z: the azimAngle is the same as deltaPhi
00358 
00359   double polarAngle (const Hep3Vector & v2, 
00360                                         const Hep3Vector & ref) const;
00361   // For arbitrary reference direction, 
00362   //    polarAngle is abs(v.angle(ref) - v2.angle(ref)).
00363 
00364   double azimAngle  (const Hep3Vector & v2, 
00365                                         const Hep3Vector & ref) const;
00366   // To compute azimangle, project v and v2 into the plane normal to
00367   // the reference direction.  Then in that plane take the angle going
00368   // clockwise around the direction from projection of v to that of v2.
00369 
00370 // 7 - Rotations 
00371 
00372 // These mehtods **DO NOT** use anything in the HepRotation class.
00373 // Thus, use of v.rotate(axis,delta) does not force linking in Rotation.cc.
00374 
00375   Hep3Vector & rotate  (const Hep3Vector & axis, double delta);
00376   // Synonym for rotate (delta, axis)
00377 
00378   Hep3Vector & rotate  (const HepAxisAngle & ax);
00379   // HepAxisAngle is a struct holding an axis direction and an angle.
00380 
00381   Hep3Vector & rotate (const HepEulerAngles & e);
00382   Hep3Vector & rotate (double phi,
00383                         double theta,
00384                         double psi);
00385   // Rotate via Euler Angles. Our Euler Angles conventions are 
00386   // those of Goldstein Classical Mechanics page 107.
00387 
00388 protected:
00389   void setSpherical (double r, double theta, double phi);
00390   void setCylindrical (double r, double phi, double z);
00391   double negativeInfinity() const;
00392 
00393 protected:
00394 
00395   double dx;
00396   double dy;
00397   double dz;
00398   // The components.
00399 
00400   static double tolerance;
00401   // default tolerance criterion for isNear() to return true.
00402 };  // Hep3Vector
00403 
00404 // Global Methods
00405 
00406 Hep3Vector rotationXOf (const Hep3Vector & vec, double delta);
00407 Hep3Vector rotationYOf (const Hep3Vector & vec, double delta);
00408 Hep3Vector rotationZOf (const Hep3Vector & vec, double delta);
00409 
00410 Hep3Vector rotationOf (const Hep3Vector & vec, 
00411                                 const Hep3Vector & axis, double delta);
00412 Hep3Vector rotationOf (const Hep3Vector & vec, const HepAxisAngle & ax);
00413 
00414 Hep3Vector rotationOf (const Hep3Vector & vec, 
00415                                 double phi, double theta, double psi);
00416 Hep3Vector rotationOf (const Hep3Vector & vec, const HepEulerAngles & e);
00417 // Return a new vector based on a rotation of the supplied vector
00418 
00419 std::ostream & operator << (std::ostream &, const Hep3Vector &);
00420 // Output to a stream.
00421 
00422 std::istream & operator >> (std::istream &, Hep3Vector &);
00423 // Input from a stream.
00424 
00425 extern const Hep3Vector HepXHat, HepYHat, HepZHat;
00426 
00427 typedef Hep3Vector HepThreeVectorD;
00428 typedef Hep3Vector HepThreeVectorF;
00429 
00430 Hep3Vector operator / (const Hep3Vector &, double a);
00431 // Division of 3-vectors by non-zero real number
00432 
00433 inline Hep3Vector operator + (const Hep3Vector &, const Hep3Vector &);
00434 // Addition of 3-vectors.
00435 
00436 inline Hep3Vector operator - (const Hep3Vector &, const Hep3Vector &);
00437 // Subtraction of 3-vectors.
00438 
00439 inline double operator * (const Hep3Vector &, const Hep3Vector &);
00440 // double product of 3-vectors.
00441 
00442 inline Hep3Vector operator * (const Hep3Vector &, double a);
00443 inline Hep3Vector operator * (double a, const Hep3Vector &);
00444 // Scaling of 3-vectors with a real number
00445 
00446 }  // namespace CLHEP
00447 
00448 #include "CLHEP/Vector/ThreeVector.icc"
00449 
00450 #ifdef ENABLE_BACKWARDS_COMPATIBILITY
00451 //  backwards compatibility will be enabled ONLY in CLHEP 1.9
00452 using namespace CLHEP;
00453 #endif
00454 
00455 #endif /* HEP_THREEVECTOR_H */

Generated on 15 Nov 2012 for CLHEP by  doxygen 1.4.7