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