CDF II
CDF KITS
source navigation ]
diff markup ]
identifier search ]
freetext search ]
file search ]
 
Architecture: i386 ]
Version: 4.10.4 ] [ 4.10.5 ] [ 4.8.4 ] [ 4.8.4l3s ] [ 4.8.5 ] [ 4.9.0 ] [ 4.9.1 ] [ 4.9.1.hpt3 ] [ 4.9.1hpt3 ] [ 4.9.1top1 ] [ 5.0.0 ] [ 5.1.0 ] [ 5.1.0beamonly ] [ 5.1.1 ] [ 5.2.0 ] [ 5.3.0 ] [ 5.3.1 ] [ 5.3.1dsp ] [ 5.3.3 ] [ 5.3.3_nt ] [ 5.3.4 ] [ 6.1.1 ] [ 6.1.1b ] [ 6.1.2 ] [ 6.1.3 ] [ 6.1.4 ] [ 6.1.4int3 ] [ 6.1.4mc ] [ 6.1.4mc_a ] [ 6.1.6 ] [ development ]

001 // -*- C++ -*-
002 // CLASSDOC OFF
003 // $Id: ThreeVector.h,v 1.12 2002/05/09 15:14:32 mf Exp $
004 // ---------------------------------------------------------------------------
005 // CLASSDOC ON
006 //
007 // This file is a part of the CLHEP - a Class Library for High Energy Physics.
008 //
009 // Hep3Vector is a general 3-vector class defining vectors in three
010 // dimension using double components. Rotations of these vectors are
011 // performed by multiplying with an object of the HepRotation class.
012 //
013 // .SS See Also
014 // LorentzVector.h, Rotation.h, LorentzRotation.h 
015 //
016 // .SS Authors
017 // Leif Lonnblad and Anders Nilsson; Modified by Evgueni Tcherniaev;
018 // ZOOM additions by Mark Fischler
019 //
020 
021 #ifndef HEP_THREEVECTOR_H
022 #define HEP_THREEVECTOR_H
023 
024 #ifdef GNUPRAGMA
025 #pragma interface
026 #endif
027 
028 #include "CLHEP/config/CLHEP.h"
029 #include "CLHEP/config/iostream.h"
030 
031 #ifdef HEP_NO_INLINE_IN_DECLARATION
032 #define inline
033 #endif
034 
035 class HepRotation;
036 class HepEulerAngles;
037 class HepAxisAngle;
038 
039 /**
040  * @author
041  * @ingroup vector
042  */
043 class Hep3Vector {
044 
045 public:
046 
047 // Basic properties and operations on 3-vectors:  
048 
049   enum { X=0, Y=1, Z=2, NUM_COORDINATES=3, SIZE=NUM_COORDINATES };
050   // Safe indexing of the coordinates when using with matrices, arrays, etc.
051   // (BaBar)
052 
053   inline Hep3Vector(double x = 0.0, double y = 0.0, double z = 0.0);
054   // The constructor.  
055 
056   inline Hep3Vector(const Hep3Vector &);
057   // The copy constructor.
058 
059   inline ~Hep3Vector();
060   // The destructor.  Not virtual - inheritance from this class is dangerous.
061 
062   double operator () (int) const;
063   // Get components by index -- 0-based (Geant4) 
064 
065   inline double operator [] (int) const;
066   // Get components by index -- 0-based (Geant4) 
067 
068   double & operator () (int);
069   // Set components by index.  0-based.
070 
071   inline double & operator [] (int);
072   // Set components by index.  0-based.
073 
074   inline double x() const;
075   inline double y() const;
076   inline double z() const;
077   // The components in cartesian coordinate system.  Same as getX() etc.
078 
079   inline void setX(double);
080   inline void setY(double);
081   inline void setZ(double);
082   // Set the components in cartesian coordinate system.
083 
084   inline void set( double x, double y, double z); 
085   // Set all three components in cartesian coordinate system.
086 
087   inline double phi() const;
088   // The azimuth angle.
089 
090   inline double theta() const;
091   // The polar angle.
092 
093   inline double cosTheta() const;
094   // Cosine of the polar angle.
095 
096   inline double cos2Theta() const;
097   // Cosine squared of the polar angle - faster than cosTheta(). (ZOOM)
098 
099   inline double mag2() const;
100   // The magnitude squared (r^2 in spherical coordinate system).
101 
102   inline double mag() const;
103   // The magnitude (r in spherical coordinate system).
104 
105   inline void setPhi(double);
106   // Set phi keeping mag and theta constant (BaBar).
107 
108   inline void setTheta(double);
109   // Set theta keeping mag and phi constant (BaBar).
110 
111          void setMag(double);
112   // Set magnitude keeping theta and phi constant (BaBar).
113 
114   inline double perp2() const;
115   // The transverse component squared (rho^2 in cylindrical coordinate system).
116 
117   inline double perp() const;
118   // The transverse component (rho in cylindrical coordinate system).
119 
120   inline void setPerp(double);
121   // Set the transverse component keeping phi and z constant.
122 
123   void setCylTheta(double);
124   // Set theta while keeping transvers component and phi fixed 
125 
126   inline double perp2(const Hep3Vector &) const;
127   // The transverse component w.r.t. given axis squared.
128 
129   inline double perp(const Hep3Vector &) const;
130   // The transverse component w.r.t. given axis.
131 
132   inline Hep3Vector & operator = (const Hep3Vector &);
133   // Assignment.
134 
135   inline bool operator == (const Hep3Vector &) const;
136   inline bool operator != (const Hep3Vector &) const;
137   // Comparisons (Geant4). 
138 
139   bool isNear (const Hep3Vector &, double epsilon=tolerance) const;
140   // Check for equality within RELATIVE tolerance (default 2.2E-14). (ZOOM)
141   // |v1 - v2|**2 <= epsilon**2 * |v1.dot(v2)| 
142 
143   double howNear(const Hep3Vector & v ) const;
144   // sqrt ( |v1-v2|**2 / v1.dot(v2) ) with a maximum of 1.
145   // If v1.dot(v2) is negative, will return 1.
146 
147   double deltaR(const Hep3Vector & v) const;
148   // sqrt( pseudorapity_difference**2 + deltaPhi **2 )
149 
150   inline Hep3Vector & operator += (const Hep3Vector &);
151   // Addition.
152 
153   inline Hep3Vector & operator -= (const Hep3Vector &);
154   // Subtraction.
155 
156   inline Hep3Vector operator - () const;
157   // Unary minus.
158 
159   inline Hep3Vector & operator *= (double);
160   // Scaling with real numbers.
161 
162          Hep3Vector & operator /= (double);
163   // Division by (non-zero) real number.
164 
165   inline Hep3Vector unit() const;
166   // Vector parallel to this, but of length 1.
167 
168   inline Hep3Vector orthogonal() const;
169   // Vector orthogonal to this (Geant4).
170 
171   inline double dot(const Hep3Vector &) const;
172   // double product.
173 
174   inline Hep3Vector cross(const Hep3Vector &) const;
175   // Cross product.
176 
177   double angle(const Hep3Vector &) const;
178   // The angle w.r.t. another 3-vector.
179 
180   double pseudoRapidity() const;
181   // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
182 
183   void setEta  ( double p );
184   // Set pseudo-rapidity, keeping magnitude and phi fixed.  (ZOOM)
185 
186   void setCylEta  ( double p );
187   // Set pseudo-rapidity, keeping transverse component and phi fixed.  (ZOOM)
188 
189   Hep3Vector & rotateX(double);
190   // Rotates the Hep3Vector around the x-axis.
191 
192   Hep3Vector & rotateY(double);
193   // Rotates the Hep3Vector around the y-axis.
194 
195   Hep3Vector & rotateZ(double);
196   // Rotates the Hep3Vector around the z-axis.
197 
198   Hep3Vector & rotateUz(const Hep3Vector&);
199   // Rotates reference frame from Uz to newUz (unit vector) (Geant4).
200 
201     Hep3Vector & rotate(double, const Hep3Vector &);
202   // Rotates around the axis specified by another Hep3Vector.
203   // (Uses methods of HepRotation, forcing linking in of Rotation.cc.)
204 
205   Hep3Vector & operator *= (const HepRotation &);
206   Hep3Vector & transform(const HepRotation &);
207   // Transformation with a Rotation matrix.
208 
209 
210 // = = = = = = = = = = = = = = = = = = = = = = = =
211 //
212 // Esoteric properties and operations on 3-vectors:  
213 //
214 // 1 - Set vectors in various coordinate systems
215 // 2 - Synonyms for accessing coordinates and properties
216 // 3 - Comparisions (dictionary, near-ness, and geometric)
217 // 4 - Intrinsic properties 
218 // 5 - Properties releative to z axis and arbitrary directions
219 // 6 - Polar and azimuthal angle decomposition and deltaPhi
220 // 7 - Rotations 
221 //
222 // = = = = = = = = = = = = = = = = = = = = = = = =
223 
224 // 1 - Set vectors in various coordinate systems
225 
226   inline void setRThetaPhi  (double r, double theta, double phi);
227   // Set in spherical coordinates:  Angles are measured in RADIANS
228 
229   inline void setREtaPhi  ( double r, double eta,  double phi );
230   // Set in spherical coordinates, but specify peudorapidiy to determine theta.
231 
232   inline void setRhoPhiZ   (double rho, double phi, double z);
233   // Set in cylindrical coordinates:  Phi angle is measured in RADIANS
234 
235   void setRhoPhiTheta ( double rho, double phi, double theta);
236   // Set in cylindrical coordinates, but specify theta to determine z.
237 
238   void setRhoPhiEta ( double rho, double phi, double eta);
239   // Set in cylindrical coordinates, but specify pseudorapidity to determine z.
240 
241 // 2 - Synonyms for accessing coordinates and properties
242 
243   inline double getX() const; 
244   inline double getY() const;
245   inline double getZ() const; 
246   // x(), y(), and z()
247 
248   inline double getR    () const;
249   inline double getTheta() const;
250   inline double getPhi  () const;
251   // mag(), theta(), and phi()
252 
253   inline double r       () const;
254   // mag()
255 
256   inline double rho     () const;
257   inline double getRho  () const;
258   // perp()
259 
260   double eta     () const;
261   double getEta  () const;
262   // pseudoRapidity() 
263 
264   inline void setR ( double s );
265   // setMag()
266 
267   inline void setRho ( double s );
268   // setPerp()
269 
270 // 3 - Comparisions (dictionary, near-ness, and geometric)
271 
272   int compare (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   bool operator<= (const Hep3Vector & v) const;
277   // dictionary ordering according to z, then y, then x component
278 
279   inline double diff2 (const Hep3Vector & v) const;
280   // |v1-v2|**2
281 
282   static double setTolerance (double tol);
283   static inline double getTolerance ();
284   // Set the tolerance used in isNear() for Hep3Vectors 
285 
286   bool isParallel (const Hep3Vector & v, double epsilon=tolerance) const;
287   // Are the vectors parallel, within the given tolerance?
288 
289   bool isOrthogonal (const Hep3Vector & v, double epsilon=tolerance) const;
290   // Are the vectors orthogonal, within the given tolerance?
291 
292   double howParallel   (const Hep3Vector & v) const;
293   // | v1.cross(v2) / v1.dot(v2) |, to a maximum of 1.
294 
295   double howOrthogonal (const Hep3Vector & v) const;
296   // | v1.dot(v2) / v1.cross(v2) |, to a maximum of 1.
297 
298   enum { ToleranceTicks = 100 };
299 
300 // 4 - Intrinsic properties 
301 
302   double beta    () const;
303   // relativistic beta (considering v as a velocity vector with c=1)
304   // Same as mag() but will object if >= 1
305 
306   double gamma() const;
307   // relativistic gamma (considering v as a velocity vector with c=1)
308 
309   double coLinearRapidity() const;
310   // inverse tanh (beta)
311 
312 // 5 - Properties relative to Z axis and to an arbitrary direction
313 
314           // Note that the non-esoteric CLHEP provides 
315           // theta(), cosTheta(), cos2Theta, and angle(const Hep3Vector&)
316 
317   inline double angle() const;
318   // angle against the Z axis -- synonym for theta()
319 
320   inline double theta(const Hep3Vector & v2) const;  
321   // synonym for angle(v2)
322 
323   double cosTheta (const Hep3Vector & v2) const;
324   double cos2Theta(const Hep3Vector & v2) const;
325   // cos and cos^2 of the angle between two vectors
326 
327   inline Hep3Vector project () const;
328          Hep3Vector project (const Hep3Vector & v2) const;
329   // projection of a vector along a direction.  
330 
331   inline Hep3Vector perpPart() const;
332   inline Hep3Vector perpPart (const Hep3Vector & v2) const;
333   // vector minus its projection along a direction.
334 
335   double rapidity () const;
336   // inverse tanh(v.z())
337 
338   double rapidity (const Hep3Vector & v2) const;
339   // rapidity with respect to specified direction:  
340   // inverse tanh (v.dot(u)) where u is a unit in the direction of v2
341 
342   double eta(const Hep3Vector & v2) const;
343   // - ln tan of the angle beween the vector and the ref direction.
344 
345 // 6 - Polar and azimuthal angle decomposition and deltaPhi
346 
347   // Decomposition of an angle within reference defined by a direction:
348 
349   double polarAngle (const Hep3Vector & v2) const;
350   // The reference direction is Z: the polarAngle is abs(v.theta()-v2.theta()).
351 
352   double deltaPhi (const Hep3Vector & v2) const;
353   // v.phi()-v2.phi(), brought into the range (-PI,PI]
354 
355   double azimAngle  (const Hep3Vector & v2) const;
356   // The reference direction is Z: the azimAngle is the same as deltaPhi
357 
358   double polarAngle (const Hep3Vector & v2, 
359                                         const Hep3Vector & ref) const;
360   // For arbitrary reference direction, 
361   //    polarAngle is abs(v.angle(ref) - v2.angle(ref)).
362 
363   double azimAngle  (const Hep3Vector & v2, 
364                                         const Hep3Vector & ref) const;
365   // To compute azimangle, project v and v2 into the plane normal to
366   // the reference direction.  Then in that plane take the angle going
367   // clockwise around the direction from projection of v to that of v2.
368 
369 // 7 - Rotations 
370 
371 // These mehtods **DO NOT** use anything in the HepRotation class.
372 // Thus, use of v.rotate(axis,delta) does not force linking in Rotation.cc.
373 
374   Hep3Vector & rotate  (const Hep3Vector & axis, double delta);
375   // Synonym for rotate (delta, axis)
376 
377   Hep3Vector & rotate  (const HepAxisAngle & ax);
378   // HepAxisAngle is a struct holding an axis direction and an angle.
379 
380   Hep3Vector & rotate (const HepEulerAngles & e);
381   Hep3Vector & rotate (double phi,
382                         double theta,
383                         double psi);
384   // Rotate via Euler Angles. Our Euler Angles conventions are 
385   // those of Goldstein Classical Mechanics page 107.
386 
387 protected:
388   void setSpherical (double r, double theta, double phi);
389   void setCylindrical (double r, double phi, double z);
390   double negativeInfinity() const;
391 
392 protected:
393 
394   double dx;
395   double dy;
396   double dz;
397   // The components.
398 
399   static double tolerance;
400   // default tolerance criterion for isNear() to return true.
401 };
402 
403 // Global Methods
404 
405 Hep3Vector rotationXOf (const Hep3Vector & vec, double delta);
406 Hep3Vector rotationYOf (const Hep3Vector & vec, double delta);
407 Hep3Vector rotationZOf (const Hep3Vector & vec, double delta);
408 
409 Hep3Vector rotationOf (const Hep3Vector & vec, 
410                                 const Hep3Vector & axis, double delta);
411 Hep3Vector rotationOf (const Hep3Vector & vec, const HepAxisAngle & ax);
412 
413 Hep3Vector rotationOf (const Hep3Vector & vec, 
414                                 double phi, double theta, double psi);
415 Hep3Vector rotationOf (const Hep3Vector & vec, const HepEulerAngles & e);
416 // Return a new vector based on a rotation of the supplied vector
417 
418 HepStd::ostream & operator << (HepStd::ostream &, const Hep3Vector &);
419 // Output to a stream.
420 
421 HepStd::istream & operator >> (HepStd::istream &, Hep3Vector &);
422 // Input from a stream.
423 
424 extern const Hep3Vector HepXHat, HepYHat, HepZHat;
425 
426 #ifdef HEP_SHORT_NAMES
427 typedef Hep3Vector Vector3;
428 typedef Hep3Vector DVector3;
429 typedef Hep3Vector FVector3;
430 static const Hep3Vector & xhat = HepXHat;
431 static const Hep3Vector & yhat = HepYHat;
432 static const Hep3Vector & zhat = HepZHat;
433 #endif
434 typedef Hep3Vector HepThreeVectorD;
435 typedef Hep3Vector HepThreeVectorF;
436 
437 Hep3Vector operator / (const Hep3Vector &, double a);
438 // Division of 3-vectors by non-zero real number
439 
440 #ifdef HEP_NO_INLINE_IN_DECLARATION
441 #undef inline
442 #endif
443 
444 #ifdef HEP_DEBUG_INLINE
445 
446 // The following methods are at global scope, and are inline in the .icc file.
447 // Likely, for some early compilers one could not declare such a method and 
448 // then define it - so one could not merely declare operator+ here.  Instead,
449 // if inline is enabled, this header you relies on the definition in the .icc
450 // file to act as a delcaration as well.  But if inline is supressed, then
451 // we do need declarations for these.
452 
453 Hep3Vector operator + (const Hep3Vector &, const Hep3Vector &);
454 // Addition of 3-vectors.
455 
456 Hep3Vector operator - (const Hep3Vector &, const Hep3Vector &);
457 // Subtraction of 3-vectors.
458 
459 double operator * (const Hep3Vector &, const Hep3Vector &);
460 // double product of 3-vectors.
461 
462 Hep3Vector operator * (const Hep3Vector &, double a);
463 Hep3Vector operator * (double a, const Hep3Vector &);
464 // Scaling of 3-vectors with a real number
465 
466 #else
467 #include "CLHEP/Vector/ThreeVector.icc"
468 #endif
469 
470 #endif /* HEP_THREEVECTOR_H */

source navigation ] diff markup ] identifier search ] freetext search ] file search ]

This page was automatically generated by the LXR engine.
The LXR team
Valid HTML 4.01!

Send problems or questions to cdfcode@fnal.gov