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