Geant4 Cross Reference |
1 // -*- C++ -*- 2 // CLASSDOC OFF 3 // --------------------------------------------------------------------------- 4 // CLASSDOC ON 5 // 6 // This file is a part of the CLHEP - a Class Library for High Energy Physics. 7 // 8 // Hep3Vector is a general 3-vector class defining vectors in three 9 // dimension using double components. Rotations of these vectors are 10 // performed by multiplying with an object of the HepRotation class. 11 // 12 // .SS See Also 13 // LorentzVector.h, Rotation.h, LorentzRotation.h 14 // 15 // .SS Authors 16 // Leif Lonnblad and Anders Nilsson; Modified by Evgueni Tcherniaev; 17 // ZOOM additions by Mark Fischler 18 // 19 20 #ifndef HEP_THREEVECTOR_H 21 #define HEP_THREEVECTOR_H 22 23 #include <iostream> 24 #include "CLHEP/Utility/defs.h" 25 26 namespace CLHEP { 27 28 class HepRotation; 29 class HepEulerAngles; 30 class HepAxisAngle; 31 32 /** 33 * @author 34 * @ingroup vector 35 */ 36 class Hep3Vector { 37 38 public: 39 40 // Basic properties and operations on 3-vectors: 41 42 enum { X=0, Y=1, Z=2, NUM_COORDINATES=3, SIZE=NUM_COORDINATES }; 43 // Safe indexing of the coordinates when using with matrices, arrays, etc. 44 // (BaBar) 45 46 Hep3Vector(); 47 explicit Hep3Vector(double x); 48 Hep3Vector(double x, double y); 49 Hep3Vector(double x, double y, double z); 50 // The constructor. 51 52 inline Hep3Vector(const Hep3Vector &); 53 inline Hep3Vector(Hep3Vector &&) = default; 54 // The copy and move constructors. 55 56 inline ~Hep3Vector(); 57 // The destructor. Not virtual - inheritance from this class is dangerous. 58 59 inline double operator () (int) const; 60 // Get components by index -- 0-based (Geant4) 61 62 inline double operator [] (int) const; 63 // Get components by index -- 0-based (Geant4) 64 65 inline double & operator () (int); 66 // Set components by index. 0-based. 67 68 inline double & operator [] (int); 69 // Set components by index. 0-based. 70 71 inline double x() const; 72 inline double y() const; 73 inline double z() const; 74 // The components in cartesian coordinate system. Same as getX() etc. 75 76 inline void setX(double); 77 inline void setY(double); 78 inline void setZ(double); 79 // Set the components in cartesian coordinate system. 80 81 inline void set( double x, double y, double z); 82 // Set all three components in cartesian coordinate system. 83 84 inline double phi() const; 85 // The azimuth angle. 86 87 inline double theta() const; 88 // The polar angle. 89 90 inline double cosTheta() const; 91 // Cosine of the polar angle. 92 93 inline double cos2Theta() const; 94 // Cosine squared of the polar angle - faster than cosTheta(). (ZOOM) 95 96 inline double mag2() const; 97 // The magnitude squared (r^2 in spherical coordinate system). 98 99 inline double mag() const; 100 // The magnitude (r in spherical coordinate system). 101 102 inline void setPhi(double); 103 // Set phi keeping mag and theta constant (BaBar). 104 105 inline void setTheta(double); 106 // Set theta keeping mag and phi constant (BaBar). 107 108 void setMag(double); 109 // Set magnitude keeping theta and phi constant (BaBar). 110 111 inline double perp2() const; 112 // The transverse component squared (rho^2 in cylindrical coordinate system). 113 114 inline double perp() const; 115 // The transverse component (rho in cylindrical coordinate system). 116 117 inline void setPerp(double); 118 // Set the transverse component keeping phi and z constant. 119 120 void setCylTheta(double); 121 // Set theta while keeping transvers component and phi fixed 122 123 inline double perp2(const Hep3Vector &) const; 124 // The transverse component w.r.t. given axis squared. 125 126 inline double perp(const Hep3Vector &) const; 127 // The transverse component w.r.t. given axis. 128 129 inline Hep3Vector & operator = (const Hep3Vector &); 130 inline Hep3Vector & operator = (Hep3Vector &&) = default; 131 // The copy and move assignment operators. 132 133 inline bool operator == (const Hep3Vector &) const; 134 inline bool operator != (const Hep3Vector &) const; 135 // Comparisons (Geant4). 136 137 bool isNear (const Hep3Vector &, double epsilon=tolerance) const; 138 // Check for equality within RELATIVE tolerance (default 2.2E-14). (ZOOM) 139 // |v1 - v2|**2 <= epsilon**2 * |v1.dot(v2)| 140 141 double howNear(const Hep3Vector & v ) const; 142 // sqrt ( |v1-v2|**2 / v1.dot(v2) ) with a maximum of 1. 143 // If v1.dot(v2) is negative, will return 1. 144 145 double deltaR(const Hep3Vector & v) const; 146 // sqrt( pseudorapity_difference**2 + deltaPhi **2 ) 147 148 inline Hep3Vector & operator += (const Hep3Vector &); 149 // Addition. 150 151 inline Hep3Vector & operator -= (const Hep3Vector &); 152 // Subtraction. 153 154 inline Hep3Vector operator - () const; 155 // Unary minus. 156 157 inline Hep3Vector & operator *= (double); 158 // Scaling with real numbers. 159 160 Hep3Vector & operator /= (double); 161 // Division by (non-zero) real number. 162 163 inline Hep3Vector unit() const; 164 // Vector parallel to this, but of length 1. 165 166 inline Hep3Vector orthogonal() const; 167 // Vector orthogonal to this (Geant4). 168 169 inline double dot(const Hep3Vector &) const; 170 // double product. 171 172 inline Hep3Vector cross(const Hep3Vector &) const; 173 // Cross product. 174 175 double angle(const Hep3Vector &) const; 176 // The angle w.r.t. another 3-vector. 177 178 double pseudoRapidity() const; 179 // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2)) 180 181 void setEta ( double p ); 182 // Set pseudo-rapidity, keeping magnitude and phi fixed. (ZOOM) 183 184 void setCylEta ( double p ); 185 // Set pseudo-rapidity, keeping transverse component and phi fixed. (ZOOM) 186 187 Hep3Vector & rotateX(double); 188 // Rotates the Hep3Vector around the x-axis. 189 190 Hep3Vector & rotateY(double); 191 // Rotates the Hep3Vector around the y-axis. 192 193 Hep3Vector & rotateZ(double); 194 // Rotates the Hep3Vector around the z-axis. 195 196 Hep3Vector & rotateUz(const Hep3Vector&); 197 // Rotates reference frame from Uz to newUz (unit vector) (Geant4). 198 199 Hep3Vector & rotate(double, const Hep3Vector &); 200 // Rotates around the axis specified by another Hep3Vector. 201 // (Uses methods of HepRotation, forcing linking in of Rotation.cc.) 202 203 Hep3Vector & operator *= (const HepRotation &); 204 Hep3Vector & transform(const HepRotation &); 205 // Transformation with a Rotation matrix. 206 207 // = = = = = = = = = = = = = = = = = = = = = = = = 208 // 209 // Esoteric properties and operations on 3-vectors: 210 // 211 // 1 - Set vectors in various coordinate systems 212 // 2 - Synonyms for accessing coordinates and properties 213 // 3 - Comparisions (dictionary, near-ness, and geometric) 214 // 4 - Intrinsic properties 215 // 5 - Properties releative to z axis and arbitrary directions 216 // 6 - Polar and azimuthal angle decomposition and deltaPhi 217 // 7 - Rotations 218 // 219 // = = = = = = = = = = = = = = = = = = = = = = = = 220 221 // 1 - Set vectors in various coordinate systems 222 223 inline void setRThetaPhi (double r, double theta, double phi); 224 // Set in spherical coordinates: Angles are measured in RADIANS 225 226 inline void setREtaPhi ( double r, double eta, double phi ); 227 // Set in spherical coordinates, but specify peudorapidiy to determine theta. 228 229 inline void setRhoPhiZ (double rho, double phi, double z); 230 // Set in cylindrical coordinates: Phi angle is measured in RADIANS 231 232 void setRhoPhiTheta ( double rho, double phi, double theta); 233 // Set in cylindrical coordinates, but specify theta to determine z. 234 235 void setRhoPhiEta ( double rho, double phi, double eta); 236 // Set in cylindrical coordinates, but specify pseudorapidity to determine z. 237 238 // 2 - Synonyms for accessing coordinates and properties 239 240 inline double getX() const; 241 inline double getY() const; 242 inline double getZ() const; 243 // x(), y(), and z() 244 245 inline double getR () const; 246 inline double getTheta() const; 247 inline double getPhi () const; 248 // mag(), theta(), and phi() 249 250 inline double r () const; 251 // mag() 252 253 inline double rho () const; 254 inline double getRho () const; 255 // perp() 256 257 double eta () const; 258 double getEta () const; 259 // pseudoRapidity() 260 261 inline void setR ( double s ); 262 // setMag() 263 264 inline void setRho ( double s ); 265 // setPerp() 266 267 // 3 - Comparisions (dictionary, near-ness, and geometric) 268 269 int compare (const Hep3Vector & v) const; 270 bool operator > (const Hep3Vector & v) const; 271 bool operator < (const Hep3Vector & v) const; 272 bool operator>= (const Hep3Vector & v) const; 273 bool operator<= (const Hep3Vector & v) const; 274 // dictionary ordering according to z, then y, then x component 275 276 inline double diff2 (const Hep3Vector & v) const; 277 // |v1-v2|**2 278 279 static double setTolerance (double tol); 280 static inline double getTolerance (); 281 // Set the tolerance used in isNear() for Hep3Vectors 282 283 bool isParallel (const Hep3Vector & v, double epsilon=tolerance) const; 284 // Are the vectors parallel, within the given tolerance? 285 286 bool isOrthogonal (const Hep3Vector & v, double epsilon=tolerance) const; 287 // Are the vectors orthogonal, within the given tolerance? 288 289 double howParallel (const Hep3Vector & v) const; 290 // | v1.cross(v2) / v1.dot(v2) |, to a maximum of 1. 291 292 double howOrthogonal (const Hep3Vector & v) const; 293 // | v1.dot(v2) / v1.cross(v2) |, to a maximum of 1. 294 295 static const int ToleranceTicks = 100; 296 297 // 4 - Intrinsic properties 298 299 double beta () const; 300 // relativistic beta (considering v as a velocity vector with c=1) 301 // Same as mag() but will object if >= 1 302 303 double gamma() const; 304 // relativistic gamma (considering v as a velocity vector with c=1) 305 306 double coLinearRapidity() const; 307 // inverse tanh (beta) 308 309 // 5 - Properties relative to Z axis and to an arbitrary direction 310 311 // Note that the non-esoteric CLHEP provides 312 // theta(), cosTheta(), cos2Theta, and angle(const Hep3Vector&) 313 314 inline double angle() const; 315 // angle against the Z axis -- synonym for theta() 316 317 inline double theta(const Hep3Vector & v2) const; 318 // synonym for angle(v2) 319 320 double cosTheta (const Hep3Vector & v2) const; 321 double cos2Theta(const Hep3Vector & v2) const; 322 // cos and cos^2 of the angle between two vectors 323 324 inline Hep3Vector project () const; 325 Hep3Vector project (const Hep3Vector & v2) const; 326 // projection of a vector along a direction. 327 328 inline Hep3Vector perpPart() const; 329 inline Hep3Vector perpPart (const Hep3Vector & v2) const; 330 // vector minus its projection along a direction. 331 332 double rapidity () const; 333 // inverse tanh(v.z()) 334 335 double rapidity (const Hep3Vector & v2) const; 336 // rapidity with respect to specified direction: 337 // inverse tanh (v.dot(u)) where u is a unit in the direction of v2 338 339 double eta(const Hep3Vector & v2) const; 340 // - ln tan of the angle beween the vector and the ref direction. 341 342 // 6 - Polar and azimuthal angle decomposition and deltaPhi 343 344 // Decomposition of an angle within reference defined by a direction: 345 346 double polarAngle (const Hep3Vector & v2) const; 347 // The reference direction is Z: the polarAngle is abs(v.theta()-v2.theta()). 348 349 double deltaPhi (const Hep3Vector & v2) const; 350 // v.phi()-v2.phi(), brought into the range (-PI,PI] 351 352 double azimAngle (const Hep3Vector & v2) const; 353 // The reference direction is Z: the azimAngle is the same as deltaPhi 354 355 double polarAngle (const Hep3Vector & v2, 356 const Hep3Vector & ref) const; 357 // For arbitrary reference direction, 358 // polarAngle is abs(v.angle(ref) - v2.angle(ref)). 359 360 double azimAngle (const Hep3Vector & v2, 361 const Hep3Vector & ref) const; 362 // To compute azimangle, project v and v2 into the plane normal to 363 // the reference direction. Then in that plane take the angle going 364 // clockwise around the direction from projection of v to that of v2. 365 366 // 7 - Rotations 367 368 // These mehtods **DO NOT** use anything in the HepRotation class. 369 // Thus, use of v.rotate(axis,delta) does not force linking in Rotation.cc. 370 371 Hep3Vector & rotate (const Hep3Vector & axis, double delta); 372 // Synonym for rotate (delta, axis) 373 374 Hep3Vector & rotate (const HepAxisAngle & ax); 375 // HepAxisAngle is a struct holding an axis direction and an angle. 376 377 Hep3Vector & rotate (const HepEulerAngles & e); 378 Hep3Vector & rotate (double phi, 379 double theta, 380 double psi); 381 // Rotate via Euler Angles. Our Euler Angles conventions are 382 // those of Goldstein Classical Mechanics page 107. 383 384 protected: 385 void setSpherical (double r, double theta, double phi); 386 void setCylindrical (double r, double phi, double z); 387 double negativeInfinity() const; 388 389 protected: 390 391 double data[3]; 392 // The components. 393 394 DLL_API static double tolerance; 395 // default tolerance criterion for isNear() to return true. 396 }; // Hep3Vector 397 398 // Global Methods 399 400 Hep3Vector rotationXOf (const Hep3Vector & vec, double delta); 401 Hep3Vector rotationYOf (const Hep3Vector & vec, double delta); 402 Hep3Vector rotationZOf (const Hep3Vector & vec, double delta); 403 404 Hep3Vector rotationOf (const Hep3Vector & vec, 405 const Hep3Vector & axis, double delta); 406 Hep3Vector rotationOf (const Hep3Vector & vec, const HepAxisAngle & ax); 407 408 Hep3Vector rotationOf (const Hep3Vector & vec, 409 double phi, double theta, double psi); 410 Hep3Vector rotationOf (const Hep3Vector & vec, const HepEulerAngles & e); 411 // Return a new vector based on a rotation of the supplied vector 412 413 std::ostream & operator << (std::ostream &, const Hep3Vector &); 414 // Output to a stream. 415 416 std::istream & operator >> (std::istream &, Hep3Vector &); 417 // Input from a stream. 418 419 extern DLL_API const Hep3Vector HepXHat, HepYHat, HepZHat; 420 421 typedef Hep3Vector HepThreeVectorD; 422 typedef Hep3Vector HepThreeVectorF; 423 424 Hep3Vector operator / (const Hep3Vector &, double a); 425 // Division of 3-vectors by non-zero real number 426 427 inline Hep3Vector operator + (const Hep3Vector &, const Hep3Vector &); 428 // Addition of 3-vectors. 429 430 inline Hep3Vector operator - (const Hep3Vector &, const Hep3Vector &); 431 // Subtraction of 3-vectors. 432 433 inline double operator * (const Hep3Vector &, const Hep3Vector &); 434 // double product of 3-vectors. 435 436 inline Hep3Vector operator * (const Hep3Vector &, double a); 437 inline Hep3Vector operator * (double a, const Hep3Vector &); 438 // Scaling of 3-vectors with a real number 439 440 } // namespace CLHEP 441 442 #include "CLHEP/Vector/ThreeVector.icc" 443 444 #endif /* HEP_THREEVECTOR_H */ 445