Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:54:36

0001 // -*- C++ -*-
0002 // $Id: BasicVector3D.h,v 1.5 2010/06/16 16:21:27 garren Exp $
0003 // ---------------------------------------------------------------------------
0004 //
0005 // This file is a part of the CLHEP - a Class Library for High Energy Physics.
0006 //
0007 // History:
0008 // 12.06.01 E.Chernyaev - CLHEP-1.7: initial  version
0009 // 14.03.03 E.Chernyaev - CLHEP-1.9: template version
0010 //
0011 
0012 #ifndef BASIC_VECTOR3D_H
0013 #define BASIC_VECTOR3D_H
0014 
0015 #include <iosfwd>
0016 #include <type_traits>
0017 #include "CLHEP/Geometry/defs.h"
0018 #include "CLHEP/Vector/ThreeVector.h"
0019 
0020 namespace HepGeom {
0021   /**
0022    * Base class for Point3D<T>, Vector3D<T> and Normal3D<T>.
0023    * It defines only common functionality for those classes and
0024    * should not be used as separate class.
0025    *
0026    * @author Evgeni Chernyaev <Evgueni.Tcherniaev@cern.ch>
0027    * @ingroup geometry
0028    */
0029   template<class T> class BasicVector3D {
0030   protected:
0031     T v_[3];
0032 
0033     /**
0034      * Default constructor.
0035      * It is protected - this class should not be instantiated directly.
0036      */
0037     BasicVector3D() { v_[0] = 0; v_[1] = 0; v_[2] = 0; }
0038 
0039   public:
0040     /**
0041      * Safe indexing of the coordinates when using with matrices, arrays, etc.
0042      */
0043     enum {
0044       X = 0,                 /**< index for x-component */
0045       Y = 1,                 /**< index for y-component */
0046       Z = 2,                 /**< index for z-component */
0047       NUM_COORDINATES = 3,   /**< number of components  */
0048       SIZE = NUM_COORDINATES /**< number of components  */
0049     };
0050 
0051     /**
0052      * Constructor from three numbers. */
0053     BasicVector3D(T x1, T y1, T z1) { v_[0] = x1; v_[1] = y1; v_[2] = z1; }
0054 
0055     /**
0056      * Copy constructor. */
0057     BasicVector3D(const BasicVector3D<T> &) = default;
0058 
0059     /**
0060      * Constructor for BasicVector3D<double> from BasicVector3D<float>. */
0061     template<typename U = T,
0062              typename = typename std::enable_if<!std::is_same<U,float>::value >::type>
0063     BasicVector3D(const BasicVector3D<float> & v) {
0064       v_[0] = v.x(); v_[1] = v.y(); v_[2] = v.z();
0065     }
0066 
0067     /**
0068      * Move constructor. */
0069     BasicVector3D(BasicVector3D<T> &&) = default;
0070 
0071     /**
0072      * Destructor. */
0073     virtual ~BasicVector3D() = default;
0074 
0075     // -------------------------
0076     // Interface to "good old C"
0077     // -------------------------
0078 
0079     /**
0080      * Conversion (cast) to ordinary array. */
0081     operator T * () { return v_; }
0082 
0083     /**
0084      * Conversion (cast) to ordinary const array. */
0085     operator const T * () const { return v_; }
0086 
0087     /**
0088      * Conversion (cast) to CLHEP::Hep3Vector.
0089      * This operator is needed only for backward compatibility and
0090      * in principle should not exit.
0091      */
0092     operator CLHEP::Hep3Vector () const { return CLHEP::Hep3Vector(x(),y(),z()); }
0093 
0094     // -----------------------------
0095     // General arithmetic operations
0096     // -----------------------------
0097 
0098     /**
0099      * Assignment. */
0100     BasicVector3D<T> & operator= (const BasicVector3D<T> &) = default;
0101     /**
0102      * Move assignment. */
0103     BasicVector3D<T> & operator= (BasicVector3D<T> &&) = default;
0104     /**
0105      * Addition. */
0106     BasicVector3D<T> & operator+=(const BasicVector3D<T> & v) {
0107       v_[0] += v.v_[0]; v_[1] += v.v_[1]; v_[2] += v.v_[2]; return *this;
0108     }
0109     /**
0110      * Subtraction. */
0111     BasicVector3D<T> & operator-=(const BasicVector3D<T> & v) {
0112       v_[0] -= v.v_[0]; v_[1] -= v.v_[1]; v_[2] -= v.v_[2]; return *this;
0113     }
0114     /**
0115      * Multiplication by scalar. */
0116     BasicVector3D<T> & operator*=(double a) {
0117       v_[0] *= a; v_[1] *= a; v_[2] *= a; return *this;
0118     }
0119     /**
0120      * Division by scalar. */
0121     BasicVector3D<T> & operator/=(double a) {
0122       v_[0] /= a; v_[1] /= a; v_[2] /= a; return *this;
0123     }
0124 
0125     // ------------
0126     // Subscripting
0127     // ------------
0128 
0129     /**
0130      * Gets components by index. */
0131     T operator()(int i) const { return v_[i]; }
0132     /**
0133      * Gets components by index. */
0134     T operator[](int i) const { return v_[i]; }
0135 
0136     /**
0137      * Sets components by index. */
0138     T & operator()(int i) { return v_[i]; }
0139     /**
0140      * Sets components by index. */
0141     T & operator[](int i) { return v_[i]; }
0142 
0143     // ------------------------------------
0144     // Cartesian coordinate system: x, y, z
0145     // ------------------------------------
0146 
0147     /**
0148      * Gets x-component in cartesian coordinate system. */
0149     T x() const { return v_[0]; }
0150     /**
0151      * Gets y-component in cartesian coordinate system. */
0152     T y() const { return v_[1]; }
0153     /**
0154      * Gets z-component in cartesian coordinate system. */
0155     T z() const { return v_[2]; }
0156 
0157     /**
0158      * Sets x-component in cartesian coordinate system. */
0159     void setX(T a) { v_[0] = a; }
0160     /**
0161      * Sets y-component in cartesian coordinate system. */
0162     void setY(T a) { v_[1] = a; }
0163     /**
0164      * Sets z-component in cartesian coordinate system. */
0165     void setZ(T a) { v_[2] = a; }
0166 
0167     /**
0168      * Sets components in cartesian coordinate system.  */
0169     void set(T x1, T y1, T z1) { v_[0] = x1; v_[1] = y1; v_[2] = z1; }
0170 
0171     // ------------------------------------------
0172     // Cylindrical coordinate system: rho, phi, z
0173     // ------------------------------------------
0174 
0175     /**
0176      * Gets transverse component squared. */
0177     T perp2() const { return x()*x()+y()*y(); }
0178     /**
0179      * Gets transverse component. */
0180     T perp() const { return std::sqrt(perp2()); }
0181     /**
0182      * Gets rho-component in cylindrical coordinate system */
0183     T rho() const { return perp(); }
0184 
0185     /**
0186      * Sets transverse component keeping phi and z constant. */
0187     void setPerp(T rh) {
0188       T factor = perp();
0189       if (factor > 0) {
0190     factor = rh/factor; v_[0] *= factor; v_[1] *= factor;
0191       }
0192     }
0193 
0194     // ------------------------------------------
0195     // Spherical coordinate system: r, phi, theta
0196     // ------------------------------------------
0197 
0198     /**
0199      * Gets magnitude squared of the vector. */
0200     T mag2() const { return x()*x()+y()*y()+z()*z(); }
0201     /**
0202      * Gets magnitude of the vector. */
0203     T mag() const { return std::sqrt(mag2()); }
0204     /**
0205      * Gets r-component in spherical coordinate system */
0206     T r() const { return mag(); }
0207     /**
0208      * Gets azimuth angle. */
0209     T phi() const {
0210       return x() == 0 && y() == 0 ? 0 : std::atan2(y(),x());
0211     }
0212     /**
0213      * Gets polar angle. */
0214     T theta() const {
0215       return x() == 0 && y() == 0 && z() == 0 ? 0 : std::atan2(perp(),z());
0216     }
0217     /**
0218      * Gets cosine of polar angle. */
0219     T cosTheta() const { T ma = mag(); return ma == 0 ? 1 : z()/ma; }
0220 
0221     /**
0222      * Gets r-component in spherical coordinate system */
0223     T getR() const { return r(); }
0224     /**
0225      * Gets phi-component in spherical coordinate system */
0226     T getPhi() const { return phi(); }
0227     /**
0228      * Gets theta-component in spherical coordinate system */
0229     T getTheta() const { return theta(); }
0230 
0231     /**
0232      * Sets magnitude. */
0233     void setMag(T ma) {
0234       T factor = mag();
0235       if (factor > 0) {
0236     factor = ma/factor; v_[0] *= factor; v_[1] *= factor; v_[2] *= factor;
0237       }
0238     }
0239     /**
0240      * Sets r-component in spherical coordinate system. */
0241     void setR(T ma) { setMag(ma); }
0242     /**
0243      * Sets phi-component in spherical coordinate system. */
0244     void setPhi(T ph) { T xy = perp(); setX(xy*std::cos(ph)); setY(xy*std::sin(ph)); }
0245     /**
0246      * Sets theta-component in spherical coordinate system. */
0247     void setTheta(T th) {
0248       T ma = mag();
0249       T ph = phi();
0250       set(ma*std::sin(th)*std::cos(ph), ma*std::sin(th)*std::sin(ph), ma*std::cos(th));
0251     }
0252 
0253     // ---------------
0254     // Pseudo rapidity
0255     // ---------------
0256 
0257     /**
0258      * Gets pseudo-rapidity: -ln(tan(theta/2)) */
0259     T pseudoRapidity() const;
0260     /**
0261      * Gets pseudo-rapidity. */
0262     T eta() const { return pseudoRapidity(); }
0263     /**
0264      * Gets pseudo-rapidity. */
0265     T getEta() const { return pseudoRapidity(); }
0266 
0267     /**
0268      * Sets pseudo-rapidity, keeping magnitude and phi fixed. */
0269     void setEta(T a);
0270 
0271     // -------------------
0272     // Combine two vectors
0273     // -------------------
0274 
0275     /**
0276      * Scalar product. */
0277     T dot(const BasicVector3D<T> & v) const {
0278       return x()*v.x()+y()*v.y()+z()*v.z();
0279     }
0280 
0281     /**
0282      * Vector product. */
0283     BasicVector3D<T> cross(const BasicVector3D<T> & v) const {
0284       return BasicVector3D<T>(y()*v.z()-v.y()*z(),
0285                   z()*v.x()-v.z()*x(),
0286                   x()*v.y()-v.x()*y());
0287     }
0288 
0289     /**
0290      * Returns transverse component w.r.t. given axis squared. */
0291     T perp2(const BasicVector3D<T> & v) const {
0292       T tot = v.mag2(), s = dot(v);
0293       return tot > 0 ? mag2()-s*s/tot : mag2();
0294     }
0295 
0296     /**
0297      * Returns transverse component w.r.t. given axis. */
0298     T perp(const BasicVector3D<T> & v) const {
0299       return std::sqrt(perp2(v));
0300     }
0301 
0302     /**
0303      * Returns angle w.r.t. another vector. */
0304     T angle(const BasicVector3D<T> & v) const;
0305 
0306     // ---------------
0307     // Related vectors
0308     // ---------------
0309 
0310     /**
0311      * Returns unit vector parallel to this. */
0312     BasicVector3D<T> unit() const {
0313       T len = mag();
0314       return (len > 0) ?
0315     BasicVector3D<T>(x()/len, y()/len, z()/len) : BasicVector3D<T>();
0316     }
0317 
0318     /**
0319      * Returns orthogonal vector. */
0320     BasicVector3D<T> orthogonal() const {
0321       T dx = x() < 0 ? -x() : x();
0322       T dy = y() < 0 ? -y() : y();
0323       T dz = z() < 0 ? -z() : z();
0324       if (dx < dy) {
0325     return dx < dz ?
0326       BasicVector3D<T>(0,z(),-y()) : BasicVector3D<T>(y(),-x(),0);
0327       }else{
0328     return dy < dz ?
0329       BasicVector3D<T>(-z(),0,x()) : BasicVector3D<T>(y(),-x(),0);
0330       }
0331     }
0332 
0333     // ---------
0334     // Rotations
0335     // ---------
0336 
0337     /**
0338      * Rotates around x-axis. */
0339     BasicVector3D<T> & rotateX(T a);
0340     /**
0341      * Rotates around y-axis. */
0342     BasicVector3D<T> & rotateY(T a);
0343     /**
0344      * Rotates around z-axis. */
0345     BasicVector3D<T> & rotateZ(T a);
0346     /**
0347      * Rotates around the axis specified by another vector. */
0348     BasicVector3D<T> & rotate(T a, const BasicVector3D<T> & v);
0349   };
0350 
0351   /*************************************************************************
0352    *                                                                       *
0353    * Non-member functions for BasicVector3D<float>                         *
0354    *                                                                       *
0355    *************************************************************************/
0356 
0357   /**
0358    * Output to stream.
0359    * @relates BasicVector3D
0360    */
0361   std::ostream &
0362   operator<<(std::ostream &, const BasicVector3D<float> &);
0363 
0364   /**
0365    * Input from stream.
0366    * @relates BasicVector3D
0367    */
0368   std::istream &
0369   operator>>(std::istream &, BasicVector3D<float> &);
0370 
0371   /**
0372    * Unary plus.
0373    * @relates BasicVector3D
0374    */
0375   inline BasicVector3D<float>
0376   operator+(const BasicVector3D<float> & v) { return v; }
0377 
0378   /**
0379    * Addition of two vectors.
0380    * @relates BasicVector3D
0381    */
0382   inline BasicVector3D<float>
0383   operator+(const BasicVector3D<float> & a, const BasicVector3D<float> & b) {
0384     return BasicVector3D<float>(a.x()+b.x(), a.y()+b.y(), a.z()+b.z());
0385   }
0386 
0387   /**
0388    * Unary minus.
0389    * @relates BasicVector3D
0390    */
0391   inline BasicVector3D<float>
0392   operator-(const BasicVector3D<float> & v) {
0393     return BasicVector3D<float>(-v.x(), -v.y(), -v.z());
0394   }
0395 
0396   /**
0397    * Subtraction of two vectors.
0398    * @relates BasicVector3D
0399    */
0400   inline BasicVector3D<float>
0401   operator-(const BasicVector3D<float> & a, const BasicVector3D<float> & b) {
0402     return BasicVector3D<float>(a.x()-b.x(), a.y()-b.y(), a.z()-b.z());
0403   }
0404 
0405   /**
0406    * Multiplication vector by scalar.
0407    * @relates BasicVector3D
0408    */
0409   inline BasicVector3D<float>
0410   operator*(const BasicVector3D<float> & v, double a) {
0411     return BasicVector3D<float>(v.x()*static_cast<float>(a), v.y()*static_cast<float>(a), v.z()*static_cast<float>(a));
0412   }
0413 
0414   /**
0415    * Scalar product of two vectors.
0416    * @relates BasicVector3D
0417    */
0418   inline float
0419   operator*(const BasicVector3D<float> & a, const BasicVector3D<float> & b) {
0420     return a.dot(b);
0421   }
0422 
0423   /**
0424    * Multiplication scalar by vector.
0425    * @relates BasicVector3D
0426    */
0427   inline BasicVector3D<float>
0428   operator*(double a, const BasicVector3D<float> & v) {
0429     return BasicVector3D<float>(static_cast<float>(a)*v.x(), static_cast<float>(a)*v.y(), static_cast<float>(a)*v.z());
0430   }
0431 
0432   /**
0433    * Division vector by scalar.
0434    * @relates BasicVector3D
0435    */
0436   inline BasicVector3D<float>
0437   operator/(const BasicVector3D<float> & v, double a) {
0438     return BasicVector3D<float>(v.x()/static_cast<float>(a), v.y()/static_cast<float>(a), v.z()/static_cast<float>(a));
0439   }
0440 
0441   /**
0442    * Comparison of two vectors for equality.
0443    * @relates BasicVector3D
0444    */
0445   inline bool
0446   operator==(const BasicVector3D<float> & a, const BasicVector3D<float> & b) {
0447     return (a.x()==b.x() && a.y()==b.y() && a.z()==b.z());
0448   }
0449 
0450   /**
0451    * Comparison of two vectors for inequality.
0452    * @relates BasicVector3D
0453    */
0454   inline bool
0455   operator!=(const BasicVector3D<float> & a, const BasicVector3D<float> & b) {
0456     return (a.x()!=b.x() || a.y()!=b.y() || a.z()!=b.z());
0457   }
0458 
0459   /*************************************************************************
0460    *                                                                       *
0461    * Non-member functions for BasicVector3D<double>                        *
0462    *                                                                       *
0463    *************************************************************************/
0464 
0465   /**
0466    * Output to stream.
0467    * @relates BasicVector3D
0468    */
0469   std::ostream &
0470   operator<<(std::ostream &, const BasicVector3D<double> &);
0471 
0472   /**
0473    * Input from stream.
0474    * @relates BasicVector3D
0475    */
0476   std::istream &
0477   operator>>(std::istream &, BasicVector3D<double> &);
0478 
0479   /**
0480    * Unary plus.
0481    * @relates BasicVector3D
0482    */
0483   inline BasicVector3D<double>
0484   operator+(const BasicVector3D<double> & v) { return v; }
0485 
0486   /**
0487    * Addition of two vectors.
0488    * @relates BasicVector3D
0489    */
0490   inline BasicVector3D<double>
0491   operator+(const BasicVector3D<double> & a,const BasicVector3D<double> & b) {
0492     return BasicVector3D<double>(a.x()+b.x(), a.y()+b.y(), a.z()+b.z());
0493   }
0494 
0495   /**
0496    * Unary minus.
0497    * @relates BasicVector3D
0498    */
0499   inline BasicVector3D<double>
0500   operator-(const BasicVector3D<double> & v) {
0501     return BasicVector3D<double>(-v.x(), -v.y(), -v.z());
0502   }
0503 
0504   /**
0505    * Subtraction of two vectors.
0506    * @relates BasicVector3D
0507    */
0508   inline BasicVector3D<double>
0509   operator-(const BasicVector3D<double> & a,const BasicVector3D<double> & b) {
0510     return BasicVector3D<double>(a.x()-b.x(), a.y()-b.y(), a.z()-b.z());
0511   }
0512 
0513   /**
0514    * Multiplication vector by scalar.
0515    * @relates BasicVector3D
0516    */
0517   inline BasicVector3D<double>
0518   operator*(const BasicVector3D<double> & v, double a) {
0519     return BasicVector3D<double>(v.x()*a, v.y()*a, v.z()*a);
0520   }
0521 
0522   /**
0523    * Scalar product of two vectors.
0524    * @relates BasicVector3D
0525    */
0526   inline double
0527   operator*(const BasicVector3D<double> & a,const BasicVector3D<double> & b) {
0528     return a.dot(b);
0529   }
0530 
0531   /**
0532    * Multiplication scalar by vector.
0533    * @relates BasicVector3D
0534    */
0535   inline BasicVector3D<double>
0536   operator*(double a, const BasicVector3D<double> & v) {
0537     return BasicVector3D<double>(a*v.x(), a*v.y(), a*v.z());
0538   }
0539 
0540   /**
0541    * Division vector by scalar.
0542    * @relates BasicVector3D
0543    */
0544   inline BasicVector3D<double>
0545   operator/(const BasicVector3D<double> & v, double a) {
0546     return BasicVector3D<double>(v.x()/a, v.y()/a, v.z()/a);
0547   }
0548 
0549   /**
0550    * Comparison of two vectors for equality.
0551    * @relates BasicVector3D
0552    */
0553   inline bool
0554   operator==(const BasicVector3D<double> & a, const BasicVector3D<double> & b)
0555   {
0556     return (a.x()==b.x() && a.y()==b.y() && a.z()==b.z());
0557   }
0558 
0559   /**
0560    * Comparison of two vectors for inequality.
0561    * @relates BasicVector3D
0562    */
0563   inline bool
0564   operator!=(const BasicVector3D<double> & a, const BasicVector3D<double> & b)
0565   {
0566     return (a.x()!=b.x() || a.y()!=b.y() || a.z()!=b.z());
0567   }
0568 } /* namespace HepGeom */
0569 
0570 #ifdef ENABLE_BACKWARDS_COMPATIBILITY
0571 //  backwards compatibility will be enabled ONLY in CLHEP 1.9
0572 using namespace HepGeom;
0573 #endif
0574 
0575 #endif /* BASIC_VECTOR3D_H */