Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BasicVector3D.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 // $Id:$
3 // ---------------------------------------------------------------------------
4 //
5 // This file is a part of the CLHEP - a Class Library for High Energy Physics.
6 //
7 // History:
8 // 12.06.01 E.Chernyaev - CLHEP-1.7: initial version
9 // 14.03.03 E.Chernyaev - CLHEP-1.9: template version
10 //
11 
12 #ifndef BASIC_VECTOR3D_H
13 #define BASIC_VECTOR3D_H
14 
15 #include <iosfwd>
17 
18 namespace HepGeom {
19  /**
20  * Base class for Point3D<T>, Vector3D<T> and Normal3D<T>.
21  * It defines only common functionality for those classes and
22  * should not be used as separate class.
23  *
24  * @author Evgeni Chernyaev <Evgueni.Tcherniaev@cern.ch>
25  * @ingroup geometry
26  */
27  template<class T> class BasicVector3D {
28  protected:
29  T v_[3];
30 
31  /**
32  * Default constructor.
33  * It is protected - this class should not be instantiated directly.
34  */
35  BasicVector3D() { v_[0] = 0; v_[1] = 0; v_[2] = 0; }
36 
37  public:
38  /**
39  * Safe indexing of the coordinates when using with matrices, arrays, etc.
40  */
41  enum {
42  X = 0, /**< index for x-component */
43  Y = 1, /**< index for y-component */
44  Z = 2, /**< index for z-component */
45  NUM_COORDINATES = 3, /**< number of components */
46  SIZE = NUM_COORDINATES /**< number of components */
47  };
48 
49  /**
50  * Constructor from three numbers. */
51  BasicVector3D(T x1, T y1, T z1) { v_[0] = x1; v_[1] = y1; v_[2] = z1; }
52 
53  /**
54  * Copy constructor.
55  * Note: BasicVector3D<double> has constructors
56  * from BasicVector3D<double> (provided by compiler) and
57  * from BasicVector3D<float> (defined in this file);
58  * BasicVector3D<float> has only the last one.
59  */
61  v_[0] = v.x(); v_[1] = v.y(); v_[2] = v.z();
62  }
63 
64  /**
65  * Destructor. */
66  virtual ~BasicVector3D() {}
67 
68  // -------------------------
69  // Interface to "good old C"
70  // -------------------------
71 
72  /**
73  * Conversion (cast) to ordinary array. */
74  operator T * () { return v_; }
75 
76  /**
77  * Conversion (cast) to ordinary const array. */
78  operator const T * () const { return v_; }
79 
80  /**
81  * Conversion (cast) to CLHEP::Hep3Vector.
82  * This operator is needed only for backward compatibility and
83  * in principle should not exit.
84  */
85  operator CLHEP::Hep3Vector () const { return CLHEP::Hep3Vector(x(),y(),z()); }
86 
87  // -----------------------------
88  // General arithmetic operations
89  // -----------------------------
90 
91  /**
92  * Assignment. */
94  v_[0] = v.v_[0]; v_[1] = v.v_[1]; v_[2] = v.v_[2]; return *this;
95  }
96  /**
97  * Addition. */
99  v_[0] += v.v_[0]; v_[1] += v.v_[1]; v_[2] += v.v_[2]; return *this;
100  }
101  /**
102  * Subtraction. */
104  v_[0] -= v.v_[0]; v_[1] -= v.v_[1]; v_[2] -= v.v_[2]; return *this;
105  }
106  /**
107  * Multiplication by scalar. */
109  v_[0] *= a; v_[1] *= a; v_[2] *= a; return *this;
110  }
111  /**
112  * Division by scalar. */
114  v_[0] /= a; v_[1] /= a; v_[2] /= a; return *this;
115  }
116 
117  // ------------
118  // Subscripting
119  // ------------
120 
121  /**
122  * Gets components by index. */
123  T operator()(int i) const { return v_[i]; }
124  /**
125  * Gets components by index. */
126  T operator[](int i) const { return v_[i]; }
127 
128  /**
129  * Sets components by index. */
130  T & operator()(int i) { return v_[i]; }
131  /**
132  * Sets components by index. */
133  T & operator[](int i) { return v_[i]; }
134 
135  // ------------------------------------
136  // Cartesian coordinate system: x, y, z
137  // ------------------------------------
138 
139  /**
140  * Gets x-component in cartesian coordinate system. */
141  T x() const { return v_[0]; }
142  /**
143  * Gets y-component in cartesian coordinate system. */
144  T y() const { return v_[1]; }
145  /**
146  * Gets z-component in cartesian coordinate system. */
147  T z() const { return v_[2]; }
148 
149  /**
150  * Sets x-component in cartesian coordinate system. */
151  void setX(T a) { v_[0] = a; }
152  /**
153  * Sets y-component in cartesian coordinate system. */
154  void setY(T a) { v_[1] = a; }
155  /**
156  * Sets z-component in cartesian coordinate system. */
157  void setZ(T a) { v_[2] = a; }
158 
159  /**
160  * Sets components in cartesian coordinate system. */
161  void set(T x1, T y1, T z1) { v_[0] = x1; v_[1] = y1; v_[2] = z1; }
162 
163  // ------------------------------------------
164  // Cylindrical coordinate system: rho, phi, z
165  // ------------------------------------------
166 
167  /**
168  * Gets transverse component squared. */
169  T perp2() const { return x()*x()+y()*y(); }
170  /**
171  * Gets transverse component. */
172  T perp() const { return std::sqrt(perp2()); }
173  /**
174  * Gets rho-component in cylindrical coordinate system */
175  T rho() const { return perp(); }
176 
177  /**
178  * Sets transverse component keeping phi and z constant. */
179  void setPerp(T rh) {
180  T factor = perp();
181  if (factor > 0) {
182  factor = rh/factor; v_[0] *= factor; v_[1] *= factor;
183  }
184  }
185 
186  // ------------------------------------------
187  // Spherical coordinate system: r, phi, theta
188  // ------------------------------------------
189 
190  /**
191  * Gets magnitude squared of the vector. */
192  T mag2() const { return x()*x()+y()*y()+z()*z(); }
193  /**
194  * Gets magnitude of the vector. */
195  T mag() const { return std::sqrt(mag2()); }
196  /**
197  * Gets r-component in spherical coordinate system */
198  T r() const { return mag(); }
199  /**
200  * Gets azimuth angle. */
201  T phi() const {
202  return x() == 0 && y() == 0 ? 0 : std::atan2(y(),x());
203  }
204  /**
205  * Gets polar angle. */
206  T theta() const {
207  return x() == 0 && y() == 0 && z() == 0 ? 0 : std::atan2(perp(),z());
208  }
209  /**
210  * Gets cosine of polar angle. */
211  T cosTheta() const { T ma = mag(); return ma == 0 ? 1 : z()/ma; }
212 
213  /**
214  * Gets r-component in spherical coordinate system */
215  T getR() const { return r(); }
216  /**
217  * Gets phi-component in spherical coordinate system */
218  T getPhi() const { return phi(); }
219  /**
220  * Gets theta-component in spherical coordinate system */
221  T getTheta() const { return theta(); }
222 
223  /**
224  * Sets magnitude. */
225  void setMag(T ma) {
226  T factor = mag();
227  if (factor > 0) {
228  factor = ma/factor; v_[0] *= factor; v_[1] *= factor; v_[2] *= factor;
229  }
230  }
231  /**
232  * Sets r-component in spherical coordinate system. */
233  void setR(T ma) { setMag(ma); }
234  /**
235  * Sets phi-component in spherical coordinate system. */
236  void setPhi(T ph) { T xy = perp(); setX(xy*std::cos(ph)); setY(xy*std::sin(ph)); }
237  /**
238  * Sets theta-component in spherical coordinate system. */
239  void setTheta(T th) {
240  T ma = mag();
241  T ph = phi();
242  set(ma*std::sin(th)*std::cos(ph), ma*std::sin(th)*std::sin(ph), ma*std::cos(th));
243  }
244 
245  // ---------------
246  // Pseudo rapidity
247  // ---------------
248 
249  /**
250  * Gets pseudo-rapidity: -ln(tan(theta/2)) */
251  T pseudoRapidity() const;
252  /**
253  * Gets pseudo-rapidity. */
254  T eta() const { return pseudoRapidity(); }
255  /**
256  * Gets pseudo-rapidity. */
257  T getEta() const { return pseudoRapidity(); }
258 
259  /**
260  * Sets pseudo-rapidity, keeping magnitude and phi fixed. */
261  void setEta(T a);
262 
263  // -------------------
264  // Combine two vectors
265  // -------------------
266 
267  /**
268  * Scalar product. */
269  T dot(const BasicVector3D<T> & v) const {
270  return x()*v.x()+y()*v.y()+z()*v.z();
271  }
272 
273  /**
274  * Vector product. */
276  return BasicVector3D<T>(y()*v.z()-v.y()*z(),
277  z()*v.x()-v.z()*x(),
278  x()*v.y()-v.x()*y());
279  }
280 
281  /**
282  * Returns transverse component w.r.t. given axis squared. */
283  T perp2(const BasicVector3D<T> & v) const {
284  T tot = v.mag2(), s = dot(v);
285  return tot > 0 ? mag2()-s*s/tot : mag2();
286  }
287 
288  /**
289  * Returns transverse component w.r.t. given axis. */
290  T perp(const BasicVector3D<T> & v) const {
291  return std::sqrt(perp2(v));
292  }
293 
294  /**
295  * Returns angle w.r.t. another vector. */
296  T angle(const BasicVector3D<T> & v) const;
297 
298  // ---------------
299  // Related vectors
300  // ---------------
301 
302  /**
303  * Returns unit vector parallel to this. */
305  T len = mag();
306  return (len > 0) ?
308  }
309 
310  /**
311  * Returns orthogonal vector. */
313  T dx = x() < 0 ? -x() : x();
314  T dy = y() < 0 ? -y() : y();
315  T dz = z() < 0 ? -z() : z();
316  if (dx < dy) {
317  return dx < dz ?
318  BasicVector3D<T>(0,z(),-y()) : BasicVector3D<T>(y(),-x(),0);
319  }else{
320  return dy < dz ?
321  BasicVector3D<T>(-z(),0,x()) : BasicVector3D<T>(y(),-x(),0);
322  }
323  }
324 
325  // ---------
326  // Rotations
327  // ---------
328 
329  /**
330  * Rotates around x-axis. */
332  /**
333  * Rotates around y-axis. */
335  /**
336  * Rotates around z-axis. */
338  /**
339  * Rotates around the axis specified by another vector. */
341  };
342 
343  /*************************************************************************
344  * *
345  * Non-member functions for BasicVector3D<float> *
346  * *
347  *************************************************************************/
348 
349  /**
350  * Output to stream.
351  * @relates BasicVector3D
352  */
353  std::ostream &
354  operator<<(std::ostream &, const BasicVector3D<float> &);
355 
356  /**
357  * Input from stream.
358  * @relates BasicVector3D
359  */
360  std::istream &
361  operator>>(std::istream &, BasicVector3D<float> &);
362 
363  /**
364  * Unary plus.
365  * @relates BasicVector3D
366  */
367  inline BasicVector3D<float>
368  operator+(const BasicVector3D<float> & v) { return v; }
369 
370  /**
371  * Addition of two vectors.
372  * @relates BasicVector3D
373  */
374  inline BasicVector3D<float>
376  return BasicVector3D<float>(a.x()+b.x(), a.y()+b.y(), a.z()+b.z());
377  }
378 
379  /**
380  * Unary minus.
381  * @relates BasicVector3D
382  */
383  inline BasicVector3D<float>
385  return BasicVector3D<float>(-v.x(), -v.y(), -v.z());
386  }
387 
388  /**
389  * Subtraction of two vectors.
390  * @relates BasicVector3D
391  */
392  inline BasicVector3D<float>
394  return BasicVector3D<float>(a.x()-b.x(), a.y()-b.y(), a.z()-b.z());
395  }
396 
397  /**
398  * Multiplication vector by scalar.
399  * @relates BasicVector3D
400  */
401  inline BasicVector3D<float>
402  operator*(const BasicVector3D<float> & v, double a) {
403  return BasicVector3D<float>(v.x()*static_cast<float>(a), v.y()*static_cast<float>(a), v.z()*static_cast<float>(a));
404  }
405 
406  /**
407  * Scalar product of two vectors.
408  * @relates BasicVector3D
409  */
410  inline float
412  return a.dot(b);
413  }
414 
415  /**
416  * Multiplication scalar by vector.
417  * @relates BasicVector3D
418  */
419  inline BasicVector3D<float>
420  operator*(double a, const BasicVector3D<float> & v) {
421  return BasicVector3D<float>(static_cast<float>(a)*v.x(), static_cast<float>(a)*v.y(), static_cast<float>(a)*v.z());
422  }
423 
424  /**
425  * Division vector by scalar.
426  * @relates BasicVector3D
427  */
428  inline BasicVector3D<float>
429  operator/(const BasicVector3D<float> & v, double a) {
430  return BasicVector3D<float>(v.x()/static_cast<float>(a), v.y()/static_cast<float>(a), v.z()/static_cast<float>(a));
431  }
432 
433  /**
434  * Comparison of two vectors for equality.
435  * @relates BasicVector3D
436  */
437  inline bool
439  return (a.x()==b.x() && a.y()==b.y() && a.z()==b.z());
440  }
441 
442  /**
443  * Comparison of two vectors for inequality.
444  * @relates BasicVector3D
445  */
446  inline bool
448  return (a.x()!=b.x() || a.y()!=b.y() || a.z()!=b.z());
449  }
450 
451  /*************************************************************************
452  * *
453  * Non-member functions for BasicVector3D<double> *
454  * *
455  *************************************************************************/
456 
457  /**
458  * Output to stream.
459  * @relates BasicVector3D
460  */
461  std::ostream &
462  operator<<(std::ostream &, const BasicVector3D<double> &);
463 
464  /**
465  * Input from stream.
466  * @relates BasicVector3D
467  */
468  std::istream &
469  operator>>(std::istream &, BasicVector3D<double> &);
470 
471  /**
472  * Unary plus.
473  * @relates BasicVector3D
474  */
475  inline BasicVector3D<double>
476  operator+(const BasicVector3D<double> & v) { return v; }
477 
478  /**
479  * Addition of two vectors.
480  * @relates BasicVector3D
481  */
482  inline BasicVector3D<double>
484  return BasicVector3D<double>(a.x()+b.x(), a.y()+b.y(), a.z()+b.z());
485  }
486 
487  /**
488  * Unary minus.
489  * @relates BasicVector3D
490  */
491  inline BasicVector3D<double>
493  return BasicVector3D<double>(-v.x(), -v.y(), -v.z());
494  }
495 
496  /**
497  * Subtraction of two vectors.
498  * @relates BasicVector3D
499  */
500  inline BasicVector3D<double>
502  return BasicVector3D<double>(a.x()-b.x(), a.y()-b.y(), a.z()-b.z());
503  }
504 
505  /**
506  * Multiplication vector by scalar.
507  * @relates BasicVector3D
508  */
509  inline BasicVector3D<double>
510  operator*(const BasicVector3D<double> & v, double a) {
511  return BasicVector3D<double>(v.x()*a, v.y()*a, v.z()*a);
512  }
513 
514  /**
515  * Scalar product of two vectors.
516  * @relates BasicVector3D
517  */
518  inline double
520  return a.dot(b);
521  }
522 
523  /**
524  * Multiplication scalar by vector.
525  * @relates BasicVector3D
526  */
527  inline BasicVector3D<double>
528  operator*(double a, const BasicVector3D<double> & v) {
529  return BasicVector3D<double>(a*v.x(), a*v.y(), a*v.z());
530  }
531 
532  /**
533  * Division vector by scalar.
534  * @relates BasicVector3D
535  */
536  inline BasicVector3D<double>
537  operator/(const BasicVector3D<double> & v, double a) {
538  return BasicVector3D<double>(v.x()/a, v.y()/a, v.z()/a);
539  }
540 
541  /**
542  * Comparison of two vectors for equality.
543  * @relates BasicVector3D
544  */
545  inline bool
547  {
548  return (a.x()==b.x() && a.y()==b.y() && a.z()==b.z());
549  }
550 
551  /**
552  * Comparison of two vectors for inequality.
553  * @relates BasicVector3D
554  */
555  inline bool
557  {
558  return (a.x()!=b.x() || a.y()!=b.y() || a.z()!=b.z());
559  }
560 } /* namespace HepGeom */
561 
562 #endif /* BASIC_VECTOR3D_H */
bool operator!=(const BasicVector3D< float > &a, const BasicVector3D< float > &b)
BasicVector3D< double > operator*(const BasicVector3D< double > &v, double a)
BasicVector3D< T > & rotate(T a, const BasicVector3D< T > &v)
BasicVector3D< float > operator*(double a, const BasicVector3D< float > &v)
BasicVector3D< double > operator/(const BasicVector3D< double > &v, double a)
void set(T x1, T y1, T z1)
BasicVector3D< double > operator+(const BasicVector3D< double > &v)
bool operator==(const BasicVector3D< double > &a, const BasicVector3D< double > &b)
BasicVector3D(const BasicVector3D< float > &v)
Definition: BasicVector3D.h:60
const XML_Char * s
T operator()(int i) const
BasicVector3D< double > operator*(double a, const BasicVector3D< double > &v)
BasicVector3D< T > unit() const
BasicVector3D< T > & operator=(const BasicVector3D< T > &v)
Definition: BasicVector3D.h:93
BasicVector3D< double > operator-(const BasicVector3D< double > &a, const BasicVector3D< double > &b)
BasicVector3D< float > operator-(const BasicVector3D< float > &v)
BasicVector3D< float > operator+(const BasicVector3D< float > &v)
BasicVector3D< double > operator-(const BasicVector3D< double > &v)
BasicVector3D< T > & operator+=(const BasicVector3D< T > &v)
Definition: BasicVector3D.h:98
double operator*(const BasicVector3D< double > &a, const BasicVector3D< double > &b)
T operator[](int i) const
BasicVector3D< double > operator+(const BasicVector3D< double > &a, const BasicVector3D< double > &b)
bool operator!=(const BasicVector3D< double > &a, const BasicVector3D< double > &b)
T perp(const BasicVector3D< T > &v) const
float operator*(const BasicVector3D< float > &a, const BasicVector3D< float > &b)
BasicVector3D< T > & operator-=(const BasicVector3D< T > &v)
bool operator==(const BasicVector3D< float > &a, const BasicVector3D< float > &b)
BasicVector3D< float > operator-(const BasicVector3D< float > &a, const BasicVector3D< float > &b)
T angle(const BasicVector3D< T > &v) const
BasicVector3D< T > & rotateX(T a)
BasicVector3D< T > orthogonal() const
BasicVector3D< float > operator*(const BasicVector3D< float > &v, double a)
BasicVector3D< T > & operator*=(double a)
std::istream & operator>>(std::istream &is, BasicVector3D< float > &a)
BasicVector3D< float > operator/(const BasicVector3D< float > &v, double a)
BasicVector3D(T x1, T y1, T z1)
Definition: BasicVector3D.h:51
const XML_Char int len
BasicVector3D< float > operator+(const BasicVector3D< float > &a, const BasicVector3D< float > &b)
T perp2(const BasicVector3D< T > &v) const
T dot(const BasicVector3D< T > &v) const
BasicVector3D< T > & rotateY(T a)
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > & operator/=(double a)
BasicVector3D< T > & rotateZ(T a)