rivet is hosted by Hepforge, IPPP Durham
Vector4.hh
Go to the documentation of this file.
00001 #ifndef RIVET_MATH_VECTOR4
00002 #define RIVET_MATH_VECTOR4
00003 
00004 #include "Rivet/Math/MathHeader.hh"
00005 #include "Rivet/Math/MathUtils.hh"
00006 #include "Rivet/Math/VectorN.hh"
00007 #include "Rivet/Math/Vector3.hh"
00008 
00009 namespace Rivet {
00010 
00011 
00012   class FourVector;
00013   class FourMomentum;
00014   class LorentzTransform;
00015   typedef FourVector Vector4;
00016   FourVector transform(const LorentzTransform& lt, const FourVector& v4);
00017 
00018 
00019   /// @brief Specialisation of VectorN to a general (non-momentum) Lorentz 4-vector.
00020   class FourVector : public Vector<4> {
00021     friend FourVector multiply(const double a, const FourVector& v);
00022     friend FourVector multiply(const FourVector& v, const double a);
00023     friend FourVector add(const FourVector& a, const FourVector& b);
00024     friend FourVector transform(const LorentzTransform& lt, const FourVector& v4);
00025 
00026   public:
00027 
00028     FourVector() : Vector<4>() { }
00029 
00030     template<typename V4>
00031     FourVector(const V4& other) {
00032       this->setT(other.t());
00033       this->setX(other.x());
00034       this->setY(other.y());
00035       this->setZ(other.z());
00036     }
00037 
00038     FourVector(const Vector<4>& other)
00039     : Vector<4>(other) { }
00040 
00041     FourVector(const double t, const double x, const double y, const double z) {
00042       this->setT(t);
00043       this->setX(x);
00044       this->setY(y);
00045       this->setZ(z);
00046     }
00047 
00048     virtual ~FourVector() { }
00049 
00050   public:
00051 
00052     double t() const { return get(0); }
00053     double x() const { return get(1); }
00054     double y() const { return get(2); }
00055     double z() const { return get(3); }
00056     FourVector& setT(const double t) { set(0, t); return *this; }
00057     FourVector& setX(const double x) { set(1, x); return *this; }
00058     FourVector& setY(const double y) { set(2, y); return *this; }
00059     FourVector& setZ(const double z) { set(3, z); return *this; }
00060 
00061     double invariant() const {
00062       // Done this way for numerical precision
00063       return (t() + z())*(t() - z()) - x()*x() - y()*y();
00064     }
00065 
00066     /// Angle between this vector and another
00067     double angle(const FourVector& v) const {
00068       return vector3().angle( v.vector3() );
00069     }
00070 
00071     /// Angle between this vector and another (3-vector)
00072     double angle(const Vector3& v3) const {
00073       return vector3().angle(v3);
00074     }
00075 
00076     /// @brief Square of the projection of the 3-vector on to the \f$ x-y \f$ plane
00077     /// This is a more efficient function than @c polarRadius, as it avoids the square root.
00078     /// Use it if you only need the squared value, or e.g. an ordering by magnitude.
00079     double polarRadius2() const {
00080       return vector3().polarRadius2();
00081     }
00082 
00083     /// Synonym for polarRadius2
00084     double perp2() const {
00085       return vector3().perp2();
00086     }
00087 
00088     /// Synonym for polarRadius2
00089     double rho2() const {
00090       return vector3().rho2();
00091     }
00092 
00093     /// Projection of 3-vector on to the \f$ x-y \f$ plane
00094     double polarRadius() const {
00095       return vector3().polarRadius();
00096     }
00097 
00098     /// Synonym for polarRadius
00099     double perp() const {
00100       return vector3().perp();
00101     }
00102 
00103     /// Synonym for polarRadius
00104     double rho() const {
00105       return vector3().rho();
00106     }
00107 
00108     /// Angle subtended by the 3-vector's projection in x-y and the x-axis.
00109     double azimuthalAngle(const PhiMapping mapping = ZERO_2PI) const {
00110       return vector3().azimuthalAngle(mapping);
00111     }
00112 
00113     /// Synonym for azimuthalAngle.
00114     double phi(const PhiMapping mapping = ZERO_2PI) const {
00115       return vector3().phi(mapping);
00116     }
00117 
00118     /// Angle subtended by the 3-vector and the z-axis.
00119     double polarAngle() const {
00120       return vector3().polarAngle();
00121     }
00122 
00123     /// Synonym for polarAngle.
00124     double theta() const {
00125       return vector3().theta();
00126     }
00127 
00128     /// Pseudorapidity (defined purely by the 3-vector components)
00129     double pseudorapidity() const {
00130       return vector3().pseudorapidity();
00131     }
00132 
00133     /// Synonym for pseudorapidity.
00134     double eta() const {
00135       return vector3().eta();
00136     }
00137 
00138     /// Get the spatial part of the 4-vector as a 3-vector.
00139     Vector3 vector3() const {
00140       return Vector3(get(1), get(2), get(3));
00141     }
00142 
00143 
00144   public:
00145 
00146     /// Contract two 4-vectors, with metric signature (+ - - -).
00147     double contract(const FourVector& v) const {
00148       const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
00149       return result;
00150     }
00151 
00152     /// Contract two 4-vectors, with metric signature (+ - - -).
00153     double dot(const FourVector& v) const {
00154       return contract(v);
00155     }
00156 
00157     /// Contract two 4-vectors, with metric signature (+ - - -).
00158     double operator*(const FourVector& v) const {
00159       return contract(v);
00160     }
00161 
00162     /// Multiply by a scalar
00163     FourVector& operator*=(double a) {
00164       _vec = multiply(a, *this)._vec;
00165       return *this;
00166     }
00167 
00168     /// Divide by a scalar
00169     FourVector& operator/=(double a) {
00170       _vec = multiply(1.0/a, *this)._vec;
00171       return *this;
00172     }
00173 
00174     FourVector& operator+=(const FourVector& v) {
00175       _vec = add(*this, v)._vec;
00176       return *this;
00177     }
00178 
00179     FourVector& operator-=(const FourVector& v) {
00180       _vec = add(*this, -v)._vec;
00181       return *this;
00182     }
00183 
00184     FourVector operator-() const {
00185       FourVector result;
00186       result._vec = -_vec;
00187       return result;
00188     }
00189 
00190   };
00191 
00192 
00193   /// Contract two 4-vectors, with metric signature (+ - - -).
00194   inline double contract(const FourVector& a, const FourVector& b) {
00195     return a.contract(b);
00196   }
00197 
00198   /// Contract two 4-vectors, with metric signature (+ - - -).
00199   inline double dot(const FourVector& a, const FourVector& b) {
00200     return contract(a, b);
00201   }
00202 
00203   inline FourVector multiply(const double a, const FourVector& v) {
00204     FourVector result;
00205     result._vec = a * v._vec;
00206     return result;
00207   }
00208 
00209   inline FourVector multiply(const FourVector& v, const double a) {
00210     return multiply(a, v);
00211   }
00212 
00213   inline FourVector operator*(const double a, const FourVector& v) {
00214     return multiply(a, v);
00215   }
00216 
00217   inline FourVector operator*(const FourVector& v, const double a) {
00218     return multiply(a, v);
00219   }
00220 
00221   inline FourVector operator/(const FourVector& v, const double a) {
00222     return multiply(1.0/a, v);
00223   }
00224 
00225   inline FourVector add(const FourVector& a, const FourVector& b) {
00226     FourVector result;
00227     result._vec = a._vec + b._vec;
00228     return result;
00229   }
00230 
00231   inline FourVector operator+(const FourVector& a, const FourVector& b) {
00232     return add(a, b);
00233   }
00234 
00235   inline FourVector operator-(const FourVector& a, const FourVector& b) {
00236     return add(a, -b);
00237   }
00238 
00239   /// Calculate the Lorentz self-invariant of a 4-vector.
00240   /// \f$ v_\mu v^\mu = g_{\mu\nu} x^\mu x^\nu \f$.
00241   inline double invariant(const FourVector& lv) {
00242     return lv.invariant();
00243   }
00244 
00245   /// Angle (in radians) between spatial parts of two Lorentz vectors.
00246   inline double angle(const FourVector& a, const FourVector& b) {
00247     return a.angle(b);
00248   }
00249 
00250   /// Angle (in radians) between spatial parts of two Lorentz vectors.
00251   inline double angle(const Vector3& a, const FourVector& b) {
00252     return angle( a, b.vector3() );
00253   }
00254 
00255   /// Angle (in radians) between spatial parts of two Lorentz vectors.
00256   inline double angle(const FourVector& a, const Vector3& b) {
00257     return a.angle(b);
00258   }
00259 
00260   /// Calculate transverse length sq. \f$ \rho^2 \f$ of a Lorentz vector.
00261   inline double polarRadius2(const FourVector& v) {
00262     return v.polarRadius2();
00263   }
00264   /// Synonym for polarRadius2.
00265   inline double perp2(const FourVector& v) {
00266     return v.perp2();
00267   }
00268   /// Synonym for polarRadius2.
00269   inline double rho2(const FourVector& v) {
00270     return v.rho2();
00271   }
00272 
00273   /// Calculate transverse length \f$ \rho \f$ of a Lorentz vector.
00274   inline double polarRadius(const FourVector& v) {
00275     return v.polarRadius();
00276   }
00277   /// Synonym for polarRadius.
00278   inline double perp(const FourVector& v) {
00279     return v.perp();
00280   }
00281   /// Synonym for polarRadius.
00282   inline double rho(const FourVector& v) {
00283     return v.rho();
00284   }
00285 
00286   /// Calculate azimuthal angle of a Lorentz vector.
00287   inline double azimuthalAngle(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00288     return v.azimuthalAngle(mapping);
00289   }
00290   /// Synonym for azimuthalAngle.
00291   inline double phi(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00292     return v.phi(mapping);
00293   }
00294 
00295 
00296   /// Calculate polar angle of a Lorentz vector.
00297   inline double polarAngle(const FourVector& v) {
00298     return v.polarAngle();
00299   }
00300   /// Synonym for polarAngle.
00301   inline double theta(const FourVector& v) {
00302     return v.theta();
00303   }
00304 
00305   /// Calculate pseudorapidity of a Lorentz vector.
00306   inline double pseudorapidity(const FourVector& v) {
00307     return v.pseudorapidity();
00308   }
00309   /// Synonym for pseudorapidity.
00310   inline double eta(const FourVector& v) {
00311     return v.eta();
00312   }
00313 
00314 
00315 
00316   ////////////////////////////////////////////////
00317 
00318 
00319 
00320   /// Specialized version of the FourVector with momentum/energy functionality.
00321   class FourMomentum : public FourVector {
00322     friend FourMomentum multiply(const double a, const FourMomentum& v);
00323     friend FourMomentum multiply(const FourMomentum& v, const double a);
00324     friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
00325     friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
00326 
00327   public:
00328     FourMomentum() { }
00329 
00330     template<typename V4>
00331     FourMomentum(const V4& other) {
00332       this->setE(other.t());
00333       this->setPx(other.x());
00334       this->setPy(other.y());
00335       this->setPz(other.z());
00336     }
00337 
00338     FourMomentum(const Vector<4>& other)
00339       : FourVector(other) { }
00340 
00341     FourMomentum(const double E, const double px, const double py, const double pz) {
00342       this->setE(E);
00343       this->setPx(px);
00344       this->setPy(py);
00345       this->setPz(pz);
00346     }
00347 
00348     ~FourMomentum() {}
00349 
00350   public:
00351     /// Get energy \f$ E \f$ (time component of momentum).
00352     double E() const { return t(); }
00353 
00354     /// Get 3-momentum part, \f$ p \f$.
00355     Vector3 p() const { return vector3(); }
00356 
00357     /// Get x-component of momentum \f$ p_x \f$.
00358     double px() const { return x(); }
00359 
00360     /// Get y-component of momentum \f$ p_y \f$.
00361     double py() const { return y(); }
00362 
00363     /// Get z-component of momentum \f$ p_z \f$.
00364     double pz() const { return z(); }
00365 
00366     /// Set energy \f$ E \f$ (time component of momentum).
00367     FourMomentum& setE(double E)   { setT(E); return *this; }
00368 
00369     /// Set x-component of momentum \f$ p_x \f$.
00370     FourMomentum& setPx(double px) { setX(px); return *this; }
00371 
00372     /// Set y-component of momentum \f$ p_y \f$.
00373     FourMomentum& setPy(double py) { setY(py); return *this; }
00374 
00375     /// Set z-component of momentum \f$ p_z \f$.
00376     FourMomentum& setPz(double pz) { setZ(pz); return *this; }
00377 
00378     /// @brief Get the mass \f$ m = \sqrt{E^2 - p^2} \f$ (the Lorentz self-invariant).
00379     ///
00380     /// For spacelike momenta, the mass will be -sqrt(|mass2|).
00381     double mass() const {
00382       // assert(Rivet::isZero(mass2()) || mass2() > 0);
00383       // if (Rivet::isZero(mass2())) {
00384       //   return 0.0;
00385       // } else {
00386       //   return sqrt(mass2());
00387       // }
00388       return sign(mass2()) * sqrt(fabs(mass2()));
00389     }
00390 
00391     /// Get the squared mass \f$ m^2 = E^2 - p^2 \f$ (the Lorentz self-invariant).
00392     double mass2() const {
00393       return invariant();
00394     }
00395 
00396     /// Calculate the rapidity.
00397     double rapidity() const {
00398       return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
00399     }
00400 
00401     /// Calculate the squared transverse momentum \f$ p_T^2 \f$.
00402     double pT2() const {
00403       return vector3().polarRadius2();
00404     }
00405 
00406     /// Calculate the transverse momentum \f$ p_T \f$.
00407     double pT() const {
00408       return sqrt(pT2());
00409     }
00410 
00411     /// Calculate the transverse energy \f$ E_T^2 = E^2 \sin^2{\theta} \f$.
00412     double Et2() const {
00413       return Et() * Et();
00414     }
00415 
00416     /// Calculate the transverse energy \f$ E_T = E \sin{\theta} \f$.
00417     double Et() const {
00418       return E() * sin(polarAngle());
00419     }
00420 
00421     /// Calculate the boost vector (in units of \f$ \beta \f$).
00422     Vector3 boostVector() const {
00423       // const Vector3 p3 = vector3();
00424       // const double m2 = mass2();
00425       // if (Rivet::isZero(m2)) return p3.unit();
00426       // else {
00427       //   // Could also do this via beta = tanh(rapidity), but that's
00428       //   // probably more messy from a numerical hygiene point of view.
00429       //   const double p2 = p3.mod2();
00430       //   const double beta = sqrt( p2 / (m2 + p2) );
00431       //   return beta * p3.unit();
00432       // }
00433       /// @todo Be careful about c=1 convention...
00434       return Vector3(px()/E(), py()/E(), pz()/E());
00435     }
00436 
00437     /// Struct for sorting by increasing energy
00438     struct byEAscending {
00439       bool operator()(const FourMomentum& left, const FourMomentum& right) const{
00440         double pt2left = left.E();
00441         double pt2right = right.E();
00442         return pt2left < pt2right;
00443       }
00444 
00445       bool operator()(const FourMomentum* left, const FourMomentum* right) const{
00446         return (*this)(left, right);
00447       }
00448     };
00449 
00450     /// Struct for sorting by decreasing energy
00451     struct byEDescending {
00452       bool operator()(const FourMomentum& left, const FourMomentum& right) const{
00453         return byEAscending()(right, left);
00454       }
00455 
00456       bool operator()(const FourMomentum* left, const FourVector* right) const{
00457         return (*this)(left, right);
00458       }
00459     };
00460 
00461 
00462     /// Multiply by a scalar
00463     FourMomentum& operator*=(double a) {
00464       _vec = multiply(a, *this)._vec;
00465       return *this;
00466     }
00467 
00468     /// Divide by a scalar
00469     FourMomentum& operator/=(double a) {
00470       _vec = multiply(1.0/a, *this)._vec;
00471       return *this;
00472     }
00473 
00474     FourMomentum& operator+=(const FourMomentum& v) {
00475       _vec = add(*this, v)._vec;
00476       return *this;
00477     }
00478 
00479     FourMomentum& operator-=(const FourMomentum& v) {
00480       _vec = add(*this, -v)._vec;
00481       return *this;
00482     }
00483 
00484     FourMomentum operator-() const {
00485       FourMomentum result;
00486       result._vec = -_vec;
00487       return result;
00488     }
00489 
00490 
00491   };
00492 
00493 
00494   inline FourMomentum multiply(const double a, const FourMomentum& v) {
00495     FourMomentum result;
00496     result._vec = a * v._vec;
00497     return result;
00498   }
00499 
00500   inline FourMomentum multiply(const FourMomentum& v, const double a) {
00501     return multiply(a, v);
00502   }
00503 
00504   inline FourMomentum operator*(const double a, const FourMomentum& v) {
00505     return multiply(a, v);
00506   }
00507 
00508   inline FourMomentum operator*(const FourMomentum& v, const double a) {
00509     return multiply(a, v);
00510   }
00511 
00512   inline FourMomentum operator/(const FourMomentum& v, const double a) {
00513     return multiply(1.0/a, v);
00514   }
00515 
00516   inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
00517     FourMomentum result;
00518     result._vec = a._vec + b._vec;
00519     return result;
00520   }
00521 
00522   inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
00523     return add(a, b);
00524   }
00525 
00526   inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
00527     return add(a, -b);
00528   }
00529 
00530 
00531 
00532   /// Get the mass \f$ m = \sqrt{E^2 - p^2} \f$ (the Lorentz self-invariant) of a momentum 4-vector.
00533   inline double mass(const FourMomentum& v) {
00534     return v.mass();
00535   }
00536 
00537   /// Get the squared mass \f$ m^2 = E^2 - p^2 \f$ (the Lorentz self-invariant) of a momentum 4-vector.
00538   inline double mass2(const FourMomentum& v) {
00539     return v.mass2();
00540   }
00541 
00542   /// Calculate the rapidity of a momentum 4-vector.
00543   inline double rapidity(const FourMomentum& v) {
00544     return v.rapidity();
00545   }
00546 
00547   /// Calculate the squared transverse momentum \f$ p_T^2 \f$ of a momentum 4-vector.
00548   inline double pT2(const FourMomentum& v) {
00549     return v.pT2();
00550   }
00551 
00552   /// Calculate the transverse momentum \f$ p_T \f$ of a momentum 4-vector.
00553   inline double pT(const FourMomentum& v) {
00554     return v.pT();
00555   }
00556 
00557   /// Calculate the transverse energy squared, \f$ E_T^2 = E^2 \sin^2{\theta} \f$ of a momentum 4-vector.
00558   inline double Et2(const FourMomentum& v) {
00559     return v.Et2();}
00560 
00561   /// Calculate the transverse energy \f$ E_T = E \sin{\theta} \f$ of a momentum 4-vector.
00562   inline double Et(const FourMomentum& v) {
00563     return v.Et();
00564   }
00565 
00566   /// Calculate the velocity boost vector of a momentum 4-vector.
00567   inline Vector3 boostVector(const FourMomentum& v) {
00568     return v.boostVector();
00569   }
00570 
00571 
00572   //////////////////////////////////////////////////////
00573 
00574 
00575   /// @name \f$ \Delta R \f$ calculations from 4-vectors
00576   //@{
00577 
00578   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00579   /// There is a scheme ambiguity for momentum-type four vectors as to whether
00580   /// the pseudorapidity (a purely geometric concept) or the rapidity (a
00581   /// relativistic energy-momentum quantity) is to be used: this can be chosen
00582   /// via the optional scheme parameter. Use of this scheme option is
00583   /// discouraged in this case since @c RAPIDITY is only a valid option for
00584   /// vectors whose type is really the FourMomentum derived class.
00585   inline double deltaR(const FourVector& a, const FourVector& b,
00586                        RapScheme scheme = PSEUDORAPIDITY) {
00587     switch (scheme) {
00588     case PSEUDORAPIDITY :
00589       return deltaR(a.vector3(), b.vector3());
00590     case RAPIDITY:
00591       {
00592         const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
00593         const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
00594         if (!ma || !mb) {
00595           string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00596           throw std::runtime_error(err);
00597         }
00598         return deltaR(*ma, *mb, scheme);
00599       }
00600     default:
00601       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00602     }
00603   }
00604 
00605 
00606   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00607   /// There is a scheme ambiguity for momentum-type four vectors
00608   /// as to whether the pseudorapidity (a purely geometric concept) or the
00609   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00610   /// be chosen via the optional scheme parameter.
00611   inline double deltaR(const FourVector& v,
00612                        double eta2, double phi2,
00613                        RapScheme scheme = PSEUDORAPIDITY) {
00614     switch (scheme) {
00615     case PSEUDORAPIDITY :
00616       return deltaR(v.vector3(), eta2, phi2);
00617     case RAPIDITY:
00618       {
00619         const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00620         if (!mv) {
00621           string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00622           throw std::runtime_error(err);
00623         }
00624         return deltaR(*mv, eta2, phi2, scheme);
00625       }
00626     default:
00627       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00628     }
00629   }
00630 
00631 
00632   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00633   /// There is a scheme ambiguity for momentum-type four vectors
00634   /// as to whether the pseudorapidity (a purely geometric concept) or the
00635   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00636   /// be chosen via the optional scheme parameter.
00637   inline double deltaR(double eta1, double phi1,
00638                        const FourVector& v,
00639                        RapScheme scheme = PSEUDORAPIDITY) {
00640     switch (scheme) {
00641     case PSEUDORAPIDITY :
00642       return deltaR(eta1, phi1, v.vector3());
00643     case RAPIDITY:
00644       {
00645         const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00646         if (!mv) {
00647           string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00648           throw std::runtime_error(err);
00649         }
00650         return deltaR(eta1, phi1, *mv, scheme);
00651       }
00652     default:
00653       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00654     }
00655   }
00656 
00657 
00658   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00659   /// There is a scheme ambiguity for momentum-type four vectors
00660   /// as to whether the pseudorapidity (a purely geometric concept) or the
00661   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00662   /// be chosen via the optional scheme parameter.
00663   inline double deltaR(const FourMomentum& a, const FourMomentum& b,
00664                        RapScheme scheme = PSEUDORAPIDITY) {
00665     switch (scheme) {
00666     case PSEUDORAPIDITY:
00667       return deltaR(a.vector3(), b.vector3());
00668     case RAPIDITY:
00669       return deltaR(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
00670     default:
00671       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00672     }
00673   }
00674 
00675   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00676   /// There is a scheme ambiguity for momentum-type four vectors
00677   /// as to whether the pseudorapidity (a purely geometric concept) or the
00678   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00679   /// be chosen via the optional scheme parameter.
00680   inline double deltaR(const FourMomentum& v,
00681                        double eta2, double phi2,
00682                        RapScheme scheme = PSEUDORAPIDITY) {
00683     switch (scheme) {
00684     case PSEUDORAPIDITY:
00685       return deltaR(v.vector3(), eta2, phi2);
00686     case RAPIDITY:
00687       return deltaR(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
00688     default:
00689       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00690     }
00691   }
00692 
00693 
00694   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00695   /// There is a scheme ambiguity for momentum-type four vectors
00696   /// as to whether the pseudorapidity (a purely geometric concept) or the
00697   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00698   /// be chosen via the optional scheme parameter.
00699   inline double deltaR(double eta1, double phi1,
00700                        const FourMomentum& v,
00701                        RapScheme scheme = PSEUDORAPIDITY) {
00702     switch (scheme) {
00703     case PSEUDORAPIDITY:
00704       return deltaR(eta1, phi1, v.vector3());
00705     case RAPIDITY:
00706       return deltaR(eta1, phi1, v.rapidity(), v.azimuthalAngle());
00707     default:
00708       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00709     }
00710   }
00711 
00712   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00713   /// There is a scheme ambiguity for momentum-type four vectors
00714   /// as to whether the pseudorapidity (a purely geometric concept) or the
00715   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00716   /// be chosen via the optional scheme parameter.
00717   inline double deltaR(const FourMomentum& a, const FourVector& b,
00718                        RapScheme scheme = PSEUDORAPIDITY) {
00719     switch (scheme) {
00720     case PSEUDORAPIDITY:
00721       return deltaR(a.vector3(), b.vector3());
00722     case RAPIDITY:
00723       return deltaR(a.rapidity(), a.azimuthalAngle(), FourMomentum(b).rapidity(), b.azimuthalAngle());
00724     default:
00725       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00726     }
00727   }
00728 
00729   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00730   /// There is a scheme ambiguity for momentum-type four vectors
00731   /// as to whether the pseudorapidity (a purely geometric concept) or the
00732   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00733   /// be chosen via the optional scheme parameter.
00734   inline double deltaR(const FourVector& a, const FourMomentum& b,
00735                        RapScheme scheme = PSEUDORAPIDITY) {
00736     switch (scheme) {
00737     case PSEUDORAPIDITY:
00738       return deltaR(a.vector3(), b.vector3());
00739     case RAPIDITY:
00740       return deltaR(FourMomentum(a).rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
00741     default:
00742       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00743     }
00744   }
00745 
00746   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00747   /// three-vector and a four-vector.
00748   inline double deltaR(const FourMomentum& a, const Vector3& b) {
00749     return deltaR(a.vector3(), b);
00750   }
00751 
00752   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00753   /// three-vector and a four-vector.
00754   inline double deltaR(const Vector3& a, const FourMomentum& b) {
00755     return deltaR(a, b.vector3());
00756   }
00757 
00758   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00759   /// three-vector and a four-vector.
00760   inline double deltaR(const FourVector& a, const Vector3& b) {
00761     return deltaR(a.vector3(), b);
00762   }
00763 
00764   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00765   /// three-vector and a four-vector.
00766   inline double deltaR(const Vector3& a, const FourVector& b) {
00767     return deltaR(a, b.vector3());
00768   }
00769 
00770   //@}
00771 
00772   /// @name \f$ \Delta phi \f$ calculations from 4-vectors
00773   //@{
00774 
00775   /// Calculate the difference in azimuthal angle between two spatial vectors.
00776   inline double deltaPhi(const FourMomentum& a, const FourMomentum& b) {
00777     return deltaPhi(a.vector3(), b.vector3());
00778   }
00779 
00780   /// Calculate the difference in azimuthal angle between two spatial vectors.
00781   inline double deltaPhi(const FourMomentum& v, double phi2) {
00782     return deltaPhi(v.vector3(), phi2);
00783   }
00784 
00785   /// Calculate the difference in azimuthal angle between two spatial vectors.
00786   inline double deltaPhi(double phi1, const FourMomentum& v) {
00787     return deltaPhi(phi1, v.vector3());
00788   }
00789 
00790   /// Calculate the difference in azimuthal angle between two spatial vectors.
00791   inline double deltaPhi(const FourVector& a, const FourVector& b) {
00792     return deltaPhi(a.vector3(), b.vector3());
00793   }
00794 
00795   /// Calculate the difference in azimuthal angle between two spatial vectors.
00796   inline double deltaPhi(const FourVector& v, double phi2) {
00797     return deltaPhi(v.vector3(), phi2);
00798   }
00799 
00800   /// Calculate the difference in azimuthal angle between two spatial vectors.
00801   inline double deltaPhi(double phi1, const FourVector& v) {
00802     return deltaPhi(phi1, v.vector3());
00803   }
00804 
00805   /// Calculate the difference in azimuthal angle between two spatial vectors.
00806   inline double deltaPhi(const FourVector& a, const FourMomentum& b) {
00807     return deltaPhi(a.vector3(), b.vector3());
00808   }
00809 
00810   /// Calculate the difference in azimuthal angle between two spatial vectors.
00811   inline double deltaPhi(const FourMomentum& a, const FourVector& b) {
00812     return deltaPhi(a.vector3(), b.vector3());
00813   }
00814 
00815   /// Calculate the difference in azimuthal angle between two spatial vectors.
00816   inline double deltaPhi(const FourVector& a, const Vector3& b) {
00817     return deltaPhi(a.vector3(), b);
00818   }
00819 
00820   /// Calculate the difference in azimuthal angle between two spatial vectors.
00821   inline double deltaPhi(const Vector3& a, const FourVector& b) {
00822     return deltaPhi(a, b.vector3());
00823   }
00824 
00825   /// Calculate the difference in azimuthal angle between two spatial vectors.
00826   inline double deltaPhi(const FourMomentum& a, const Vector3& b) {
00827     return deltaPhi(a.vector3(), b);
00828   }
00829 
00830   /// Calculate the difference in azimuthal angle between two spatial vectors.
00831   inline double deltaPhi(const Vector3& a, const FourMomentum& b) {
00832     return deltaPhi(a, b.vector3());
00833   }
00834 
00835   //@}
00836 
00837   /// @name \f$ |\Delta eta| \f$ calculations from 4-vectors
00838   //@{
00839 
00840   /// Calculate the difference in azimuthal angle between two spatial vectors.
00841   inline double deltaEta(const FourMomentum& a, const FourMomentum& b) {
00842     return deltaEta(a.vector3(), b.vector3());
00843   }
00844 
00845   /// Calculate the difference in azimuthal angle between two spatial vectors.
00846   inline double deltaEta(const FourMomentum& v, double eta2) {
00847     return deltaEta(v.vector3(), eta2);
00848   }
00849 
00850   /// Calculate the difference in azimuthal angle between two spatial vectors.
00851   inline double deltaEta(double eta1, const FourMomentum& v) {
00852     return deltaEta(eta1, v.vector3());
00853   }
00854 
00855   /// Calculate the difference in azimuthal angle between two spatial vectors.
00856   inline double deltaEta(const FourVector& a, const FourVector& b) {
00857     return deltaEta(a.vector3(), b.vector3());
00858   }
00859 
00860   /// Calculate the difference in azimuthal angle between two spatial vectors.
00861   inline double deltaEta(const FourVector& v, double eta2) {
00862     return deltaEta(v.vector3(), eta2);
00863   }
00864 
00865   /// Calculate the difference in azimuthal angle between two spatial vectors.
00866   inline double deltaEta(double eta1, const FourVector& v) {
00867     return deltaEta(eta1, v.vector3());
00868   }
00869 
00870   /// Calculate the difference in azimuthal angle between two spatial vectors.
00871   inline double deltaEta(const FourVector& a, const FourMomentum& b) {
00872     return deltaEta(a.vector3(), b.vector3());
00873   }
00874 
00875   /// Calculate the difference in azimuthal angle between two spatial vectors.
00876   inline double deltaEta(const FourMomentum& a, const FourVector& b) {
00877     return deltaEta(a.vector3(), b.vector3());
00878   }
00879 
00880   /// Calculate the difference in azimuthal angle between two spatial vectors.
00881   inline double deltaEta(const FourVector& a, const Vector3& b) {
00882     return deltaEta(a.vector3(), b);
00883   }
00884 
00885   /// Calculate the difference in azimuthal angle between two spatial vectors.
00886   inline double deltaEta(const Vector3& a, const FourVector& b) {
00887     return deltaEta(a, b.vector3());
00888   }
00889 
00890   /// Calculate the difference in azimuthal angle between two spatial vectors.
00891   inline double deltaEta(const FourMomentum& a, const Vector3& b) {
00892     return deltaEta(a.vector3(), b);
00893   }
00894 
00895   /// Calculate the difference in azimuthal angle between two spatial vectors.
00896   inline double deltaEta(const Vector3& a, const FourMomentum& b) {
00897     return deltaEta(a, b.vector3());
00898   }
00899 
00900   //@}
00901 
00902   //////////////////////////////////////////////////////
00903 
00904 
00905   /// @name 4-vector string representations
00906   //@{
00907 
00908   /// Render a 4-vector as a string.
00909   inline std::string toString(const FourVector& lv) {
00910     ostringstream out;
00911     out << "("  << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
00912         << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
00913         << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
00914         << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
00915         << ")";
00916     return out.str();
00917   }
00918 
00919   /// Write a 4-vector to an ostream.
00920   inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
00921     out << toString(lv);
00922     return out;
00923   }
00924 
00925   //@}
00926 
00927 
00928 }
00929 
00930 #endif