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 \f$ |\eta| \f$ directly.
00139     double abspseudorapidity() const { return fabs(eta()); }
00140 
00141     /// Get the \f$ |\eta| \f$ directly (alias).
00142     double abseta() const { return fabs(eta()); }
00143 
00144     /// Get the spatial part of the 4-vector as a 3-vector.
00145     Vector3 vector3() const {
00146       return Vector3(get(1), get(2), get(3));
00147     }
00148 
00149 
00150   public:
00151 
00152     /// Contract two 4-vectors, with metric signature (+ - - -).
00153     double contract(const FourVector& v) const {
00154       const double result = t()*v.t() - x()*v.x() - y()*v.y() - z()*v.z();
00155       return result;
00156     }
00157 
00158     /// Contract two 4-vectors, with metric signature (+ - - -).
00159     double dot(const FourVector& v) const {
00160       return contract(v);
00161     }
00162 
00163     /// Contract two 4-vectors, with metric signature (+ - - -).
00164     double operator*(const FourVector& v) const {
00165       return contract(v);
00166     }
00167 
00168     /// Multiply by a scalar.
00169     FourVector& operator*=(double a) {
00170       _vec = multiply(a, *this)._vec;
00171       return *this;
00172     }
00173 
00174     /// Divide by a scalar.
00175     FourVector& operator/=(double a) {
00176       _vec = multiply(1.0/a, *this)._vec;
00177       return *this;
00178     }
00179 
00180     /// Add to this 4-vector.
00181     FourVector& operator+=(const FourVector& v) {
00182       _vec = add(*this, v)._vec;
00183       return *this;
00184     }
00185 
00186     /// Subtract from this 4-vector. NB time as well as space components are subtracted.
00187     FourVector& operator-=(const FourVector& v) {
00188       _vec = add(*this, -v)._vec;
00189       return *this;
00190     }
00191 
00192     /// Multiply all components (space and time) by -1.
00193     FourVector operator-() const {
00194       FourVector result;
00195       result._vec = -_vec;
00196       return result;
00197     }
00198 
00199   };
00200 
00201 
00202   /// Contract two 4-vectors, with metric signature (+ - - -).
00203   inline double contract(const FourVector& a, const FourVector& b) {
00204     return a.contract(b);
00205   }
00206 
00207   /// Contract two 4-vectors, with metric signature (+ - - -).
00208   inline double dot(const FourVector& a, const FourVector& b) {
00209     return contract(a, b);
00210   }
00211 
00212   inline FourVector multiply(const double a, const FourVector& v) {
00213     FourVector result;
00214     result._vec = a * v._vec;
00215     return result;
00216   }
00217 
00218   inline FourVector multiply(const FourVector& v, const double a) {
00219     return multiply(a, v);
00220   }
00221 
00222   inline FourVector operator*(const double a, const FourVector& v) {
00223     return multiply(a, v);
00224   }
00225 
00226   inline FourVector operator*(const FourVector& v, const double a) {
00227     return multiply(a, v);
00228   }
00229 
00230   inline FourVector operator/(const FourVector& v, const double a) {
00231     return multiply(1.0/a, v);
00232   }
00233 
00234   inline FourVector add(const FourVector& a, const FourVector& b) {
00235     FourVector result;
00236     result._vec = a._vec + b._vec;
00237     return result;
00238   }
00239 
00240   inline FourVector operator+(const FourVector& a, const FourVector& b) {
00241     return add(a, b);
00242   }
00243 
00244   inline FourVector operator-(const FourVector& a, const FourVector& b) {
00245     return add(a, -b);
00246   }
00247 
00248   /// Calculate the Lorentz self-invariant of a 4-vector.
00249   /// \f$ v_\mu v^\mu = g_{\mu\nu} x^\mu x^\nu \f$.
00250   inline double invariant(const FourVector& lv) {
00251     return lv.invariant();
00252   }
00253 
00254   /// Angle (in radians) between spatial parts of two Lorentz vectors.
00255   inline double angle(const FourVector& a, const FourVector& b) {
00256     return a.angle(b);
00257   }
00258 
00259   /// Angle (in radians) between spatial parts of two Lorentz vectors.
00260   inline double angle(const Vector3& a, const FourVector& b) {
00261     return angle( a, b.vector3() );
00262   }
00263 
00264   /// Angle (in radians) between spatial parts of two Lorentz vectors.
00265   inline double angle(const FourVector& a, const Vector3& b) {
00266     return a.angle(b);
00267   }
00268 
00269   /// Calculate transverse length sq. \f$ \rho^2 \f$ of a Lorentz vector.
00270   inline double polarRadius2(const FourVector& v) {
00271     return v.polarRadius2();
00272   }
00273   /// Synonym for polarRadius2.
00274   inline double perp2(const FourVector& v) {
00275     return v.perp2();
00276   }
00277   /// Synonym for polarRadius2.
00278   inline double rho2(const FourVector& v) {
00279     return v.rho2();
00280   }
00281 
00282   /// Calculate transverse length \f$ \rho \f$ of a Lorentz vector.
00283   inline double polarRadius(const FourVector& v) {
00284     return v.polarRadius();
00285   }
00286   /// Synonym for polarRadius.
00287   inline double perp(const FourVector& v) {
00288     return v.perp();
00289   }
00290   /// Synonym for polarRadius.
00291   inline double rho(const FourVector& v) {
00292     return v.rho();
00293   }
00294 
00295   /// Calculate azimuthal angle of a Lorentz vector.
00296   inline double azimuthalAngle(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00297     return v.azimuthalAngle(mapping);
00298   }
00299   /// Synonym for azimuthalAngle.
00300   inline double phi(const FourVector& v, const PhiMapping mapping = ZERO_2PI) {
00301     return v.phi(mapping);
00302   }
00303 
00304 
00305   /// Calculate polar angle of a Lorentz vector.
00306   inline double polarAngle(const FourVector& v) {
00307     return v.polarAngle();
00308   }
00309   /// Synonym for polarAngle.
00310   inline double theta(const FourVector& v) {
00311     return v.theta();
00312   }
00313 
00314   /// Calculate pseudorapidity of a Lorentz vector.
00315   inline double pseudorapidity(const FourVector& v) {
00316     return v.pseudorapidity();
00317   }
00318   /// Synonym for pseudorapidity.
00319   inline double eta(const FourVector& v) {
00320     return v.eta();
00321   }
00322 
00323 
00324   /// Calculate absolute pseudorapidity of a Lorentz vector.
00325   inline double abspseudorapidity(const FourVector& v) {
00326     return v.abspseudorapidity();
00327   }
00328   /// Synonym for absolute pseudorapidity.
00329   inline double abseta(const FourVector& v) {
00330     return v.abseta();
00331   }
00332 
00333 
00334 
00335   ////////////////////////////////////////////////
00336 
00337 
00338 
00339   /// Specialized version of the FourVector with momentum/energy functionality.
00340   class FourMomentum : public FourVector {
00341     friend FourMomentum multiply(const double a, const FourMomentum& v);
00342     friend FourMomentum multiply(const FourMomentum& v, const double a);
00343     friend FourMomentum add(const FourMomentum& a, const FourMomentum& b);
00344     friend FourMomentum transform(const LorentzTransform& lt, const FourMomentum& v4);
00345 
00346   public:
00347     FourMomentum() { }
00348 
00349     template<typename V4>
00350     FourMomentum(const V4& other) {
00351       this->setE(other.t());
00352       this->setPx(other.x());
00353       this->setPy(other.y());
00354       this->setPz(other.z());
00355     }
00356 
00357     FourMomentum(const Vector<4>& other)
00358       : FourVector(other) { }
00359 
00360     FourMomentum(const double E, const double px, const double py, const double pz) {
00361       this->setE(E);
00362       this->setPx(px);
00363       this->setPy(py);
00364       this->setPz(pz);
00365     }
00366 
00367     ~FourMomentum() {}
00368 
00369   public:
00370     /// Get energy \f$ E \f$ (time component of momentum).
00371     double E() const { return t(); }
00372 
00373     /// Get 3-momentum part, \f$ p \f$.
00374     Vector3 p() const { return vector3(); }
00375 
00376     /// Get x-component of momentum \f$ p_x \f$.
00377     double px() const { return x(); }
00378 
00379     /// Get y-component of momentum \f$ p_y \f$.
00380     double py() const { return y(); }
00381 
00382     /// Get z-component of momentum \f$ p_z \f$.
00383     double pz() const { return z(); }
00384 
00385     /// Set energy \f$ E \f$ (time component of momentum).
00386     FourMomentum& setE(double E)   { setT(E); return *this; }
00387 
00388     /// Set x-component of momentum \f$ p_x \f$.
00389     FourMomentum& setPx(double px) { setX(px); return *this; }
00390 
00391     /// Set y-component of momentum \f$ p_y \f$.
00392     FourMomentum& setPy(double py) { setY(py); return *this; }
00393 
00394     /// Set z-component of momentum \f$ p_z \f$.
00395     FourMomentum& setPz(double pz) { setZ(pz); return *this; }
00396 
00397     /// @brief Get the mass \f$ m = \sqrt{E^2 - p^2} \f$ (the Lorentz self-invariant).
00398     ///
00399     /// For spacelike momenta, the mass will be -sqrt(|mass2|).
00400     double mass() const {
00401       // assert(Rivet::isZero(mass2()) || mass2() > 0);
00402       // if (Rivet::isZero(mass2())) {
00403       //   return 0.0;
00404       // } else {
00405       //   return sqrt(mass2());
00406       // }
00407       return sign(mass2()) * sqrt(fabs(mass2()));
00408     }
00409 
00410     /// Get the squared mass \f$ m^2 = E^2 - p^2 \f$ (the Lorentz self-invariant).
00411     double mass2() const {
00412       return invariant();
00413     }
00414 
00415     /// Calculate the rapidity.
00416     double rapidity() const {
00417       return 0.5 * std::log( (E() + pz()) / (E() - pz()) );
00418     }
00419 
00420     /// Alias for rapidity.
00421     double rap() const {
00422       return rapidity();
00423     }
00424 
00425     /// Absolute rapidity.
00426     double absrapidity() const {
00427       return fabs(rapidity());
00428     }
00429 
00430     /// Absolute rapidity.
00431     double absrap() const {
00432       return fabs(rap());
00433     }
00434 
00435     /// Calculate the squared transverse momentum \f$ p_T^2 \f$.
00436     double pT2() const {
00437       return vector3().polarRadius2();
00438     }
00439 
00440     /// Calculate the squared transverse momentum \f$ p_T^2 \f$.
00441     double pt2() const {
00442       return vector3().polarRadius2();
00443     }
00444 
00445     /// Calculate the transverse momentum \f$ p_T \f$.
00446     double pT() const {
00447       return sqrt(pT2());
00448     }
00449 
00450     /// Calculate the transverse momentum \f$ p_T \f$.
00451     double pt() const {
00452       return sqrt(pT2());
00453     }
00454 
00455     /// Calculate the transverse energy \f$ E_T^2 = E^2 \sin^2{\theta} \f$.
00456     double Et2() const {
00457       return Et() * Et();
00458     }
00459 
00460     /// Calculate the transverse energy \f$ E_T = E \sin{\theta} \f$.
00461     double Et() const {
00462       return E() * sin(polarAngle());
00463     }
00464 
00465     /// Calculate the boost vector (in units of \f$ \beta \f$).
00466     Vector3 boostVector() const {
00467       // const Vector3 p3 = vector3();
00468       // const double m2 = mass2();
00469       // if (Rivet::isZero(m2)) return p3.unit();
00470       // else {
00471       //   // Could also do this via beta = tanh(rapidity), but that's
00472       //   // probably more messy from a numerical hygiene point of view.
00473       //   const double p2 = p3.mod2();
00474       //   const double beta = sqrt( p2 / (m2 + p2) );
00475       //   return beta * p3.unit();
00476       // }
00477       /// @todo Be careful about c=1 convention...
00478       return Vector3(px()/E(), py()/E(), pz()/E());
00479     }
00480 
00481     /// Struct for sorting by increasing energy
00482     struct byEAscending {
00483       bool operator()(const FourMomentum& left, const FourMomentum& right) const{
00484         double pt2left = left.E();
00485         double pt2right = right.E();
00486         return pt2left < pt2right;
00487       }
00488 
00489       bool operator()(const FourMomentum* left, const FourMomentum* right) const{
00490         return (*this)(left, right);
00491       }
00492     };
00493 
00494     /// Struct for sorting by decreasing energy
00495     struct byEDescending {
00496       bool operator()(const FourMomentum& left, const FourMomentum& right) const{
00497         return byEAscending()(right, left);
00498       }
00499 
00500       bool operator()(const FourMomentum* left, const FourVector* right) const{
00501         return (*this)(left, right);
00502       }
00503     };
00504 
00505 
00506     /// Multiply by a scalar
00507     FourMomentum& operator*=(double a) {
00508       _vec = multiply(a, *this)._vec;
00509       return *this;
00510     }
00511 
00512     /// Divide by a scalar
00513     FourMomentum& operator/=(double a) {
00514       _vec = multiply(1.0/a, *this)._vec;
00515       return *this;
00516     }
00517 
00518     /// Add to this 4-vector. NB time as well as space components are added.
00519     FourMomentum& operator+=(const FourMomentum& v) {
00520       _vec = add(*this, v)._vec;
00521       return *this;
00522     }
00523 
00524     /// Subtract from this 4-vector. NB time as well as space components are subtracted.
00525     FourMomentum& operator-=(const FourMomentum& v) {
00526       _vec = add(*this, -v)._vec;
00527       return *this;
00528     }
00529 
00530     /// Multiply all components (time and space) by -1.
00531     FourMomentum operator-() const {
00532       FourMomentum result;
00533       result._vec = -_vec;
00534       return result;
00535     }
00536 
00537 
00538   };
00539 
00540 
00541   inline FourMomentum multiply(const double a, const FourMomentum& v) {
00542     FourMomentum result;
00543     result._vec = a * v._vec;
00544     return result;
00545   }
00546 
00547   inline FourMomentum multiply(const FourMomentum& v, const double a) {
00548     return multiply(a, v);
00549   }
00550 
00551   inline FourMomentum operator*(const double a, const FourMomentum& v) {
00552     return multiply(a, v);
00553   }
00554 
00555   inline FourMomentum operator*(const FourMomentum& v, const double a) {
00556     return multiply(a, v);
00557   }
00558 
00559   inline FourMomentum operator/(const FourMomentum& v, const double a) {
00560     return multiply(1.0/a, v);
00561   }
00562 
00563   inline FourMomentum add(const FourMomentum& a, const FourMomentum& b) {
00564     FourMomentum result;
00565     result._vec = a._vec + b._vec;
00566     return result;
00567   }
00568 
00569   inline FourMomentum operator+(const FourMomentum& a, const FourMomentum& b) {
00570     return add(a, b);
00571   }
00572 
00573   inline FourMomentum operator-(const FourMomentum& a, const FourMomentum& b) {
00574     return add(a, -b);
00575   }
00576 
00577 
00578 
00579   /// Get the mass \f$ m = \sqrt{E^2 - p^2} \f$ (the Lorentz self-invariant) of a momentum 4-vector.
00580   inline double mass(const FourMomentum& v) {
00581     return v.mass();
00582   }
00583 
00584   /// Get the squared mass \f$ m^2 = E^2 - p^2 \f$ (the Lorentz self-invariant) of a momentum 4-vector.
00585   inline double mass2(const FourMomentum& v) {
00586     return v.mass2();
00587   }
00588 
00589   /// Calculate the rapidity of a momentum 4-vector.
00590   inline double rapidity(const FourMomentum& v) {
00591     return v.rapidity();
00592   }
00593 
00594   /// Calculate the rapidity of a momentum 4-vector.
00595   inline double rap(const FourMomentum& v) {
00596     return v.rap();
00597   }
00598 
00599   /// Calculate the absolute rapidity of a momentum 4-vector.
00600   inline double absrapidity(const FourMomentum& v) {
00601     return v.absrapidity();
00602   }
00603 
00604   /// Calculate the absolute rapidity of a momentum 4-vector.
00605   inline double absrap(const FourMomentum& v) {
00606     return v.absrap();
00607   }
00608 
00609   /// Calculate the squared transverse momentum \f$ p_T^2 \f$ of a momentum 4-vector.
00610   inline double pT2(const FourMomentum& v) {
00611     return v.pT2();
00612   }
00613 
00614   /// Calculate the transverse momentum \f$ p_T \f$ of a momentum 4-vector.
00615   inline double pT(const FourMomentum& v) {
00616     return v.pT();
00617   }
00618 
00619   /// Calculate the transverse energy squared, \f$ E_T^2 = E^2 \sin^2{\theta} \f$ of a momentum 4-vector.
00620   inline double Et2(const FourMomentum& v) {
00621     return v.Et2();}
00622 
00623   /// Calculate the transverse energy \f$ E_T = E \sin{\theta} \f$ of a momentum 4-vector.
00624   inline double Et(const FourMomentum& v) {
00625     return v.Et();
00626   }
00627 
00628   /// Calculate the velocity boost vector of a momentum 4-vector.
00629   inline Vector3 boostVector(const FourMomentum& v) {
00630     return v.boostVector();
00631   }
00632 
00633 
00634   //////////////////////////////////////////////////////
00635 
00636 
00637   /// @name \f$ \Delta R \f$ calculations from 4-vectors
00638   //@{
00639 
00640   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00641   /// There is a scheme ambiguity for momentum-type four vectors as to whether
00642   /// the pseudorapidity (a purely geometric concept) or the rapidity (a
00643   /// relativistic energy-momentum quantity) is to be used: this can be chosen
00644   /// via the optional scheme parameter. Use of this scheme option is
00645   /// discouraged in this case since @c RAPIDITY is only a valid option for
00646   /// vectors whose type is really the FourMomentum derived class.
00647   inline double deltaR(const FourVector& a, const FourVector& b,
00648                        RapScheme scheme = PSEUDORAPIDITY) {
00649     switch (scheme) {
00650     case PSEUDORAPIDITY :
00651       return deltaR(a.vector3(), b.vector3());
00652     case RAPIDITY:
00653       {
00654         const FourMomentum* ma = dynamic_cast<const FourMomentum*>(&a);
00655         const FourMomentum* mb = dynamic_cast<const FourMomentum*>(&b);
00656         if (!ma || !mb) {
00657           string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00658           throw std::runtime_error(err);
00659         }
00660         return deltaR(*ma, *mb, scheme);
00661       }
00662     default:
00663       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00664     }
00665   }
00666 
00667 
00668   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00669   /// There is a scheme ambiguity for momentum-type four vectors
00670   /// as to whether the pseudorapidity (a purely geometric concept) or the
00671   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00672   /// be chosen via the optional scheme parameter.
00673   inline double deltaR(const FourVector& v,
00674                        double eta2, double phi2,
00675                        RapScheme scheme = PSEUDORAPIDITY) {
00676     switch (scheme) {
00677     case PSEUDORAPIDITY :
00678       return deltaR(v.vector3(), eta2, phi2);
00679     case RAPIDITY:
00680       {
00681         const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00682         if (!mv) {
00683           string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00684           throw std::runtime_error(err);
00685         }
00686         return deltaR(*mv, eta2, phi2, scheme);
00687       }
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 FourVector& v,
00701                        RapScheme scheme = PSEUDORAPIDITY) {
00702     switch (scheme) {
00703     case PSEUDORAPIDITY :
00704       return deltaR(eta1, phi1, v.vector3());
00705     case RAPIDITY:
00706       {
00707         const FourMomentum* mv = dynamic_cast<const FourMomentum*>(&v);
00708         if (!mv) {
00709           string err = "deltaR with scheme RAPIDITY can only be called with FourMomentum objects, not FourVectors";
00710           throw std::runtime_error(err);
00711         }
00712         return deltaR(eta1, phi1, *mv, scheme);
00713       }
00714     default:
00715       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00716     }
00717   }
00718 
00719 
00720   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00721   /// There is a scheme ambiguity for momentum-type four vectors
00722   /// as to whether the pseudorapidity (a purely geometric concept) or the
00723   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00724   /// be chosen via the optional scheme parameter.
00725   inline double deltaR(const FourMomentum& a, const FourMomentum& b,
00726                        RapScheme scheme = PSEUDORAPIDITY) {
00727     switch (scheme) {
00728     case PSEUDORAPIDITY:
00729       return deltaR(a.vector3(), b.vector3());
00730     case RAPIDITY:
00731       return deltaR(a.rapidity(), a.azimuthalAngle(), b.rapidity(), b.azimuthalAngle());
00732     default:
00733       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00734     }
00735   }
00736 
00737   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00738   /// There is a scheme ambiguity for momentum-type four vectors
00739   /// as to whether the pseudorapidity (a purely geometric concept) or the
00740   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00741   /// be chosen via the optional scheme parameter.
00742   inline double deltaR(const FourMomentum& v,
00743                        double eta2, double phi2,
00744                        RapScheme scheme = PSEUDORAPIDITY) {
00745     switch (scheme) {
00746     case PSEUDORAPIDITY:
00747       return deltaR(v.vector3(), eta2, phi2);
00748     case RAPIDITY:
00749       return deltaR(v.rapidity(), v.azimuthalAngle(), eta2, phi2);
00750     default:
00751       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00752     }
00753   }
00754 
00755 
00756   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00757   /// There is a scheme ambiguity for momentum-type four vectors
00758   /// as to whether the pseudorapidity (a purely geometric concept) or the
00759   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00760   /// be chosen via the optional scheme parameter.
00761   inline double deltaR(double eta1, double phi1,
00762                        const FourMomentum& v,
00763                        RapScheme scheme = PSEUDORAPIDITY) {
00764     switch (scheme) {
00765     case PSEUDORAPIDITY:
00766       return deltaR(eta1, phi1, v.vector3());
00767     case RAPIDITY:
00768       return deltaR(eta1, phi1, v.rapidity(), v.azimuthalAngle());
00769     default:
00770       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00771     }
00772   }
00773 
00774   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00775   /// There is a scheme ambiguity for momentum-type four vectors
00776   /// as to whether the pseudorapidity (a purely geometric concept) or the
00777   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00778   /// be chosen via the optional scheme parameter.
00779   inline double deltaR(const FourMomentum& a, const FourVector& b,
00780                        RapScheme scheme = PSEUDORAPIDITY) {
00781     switch (scheme) {
00782     case PSEUDORAPIDITY:
00783       return deltaR(a.vector3(), b.vector3());
00784     case RAPIDITY:
00785       return deltaR(a.rapidity(), a.azimuthalAngle(), FourMomentum(b).rapidity(), b.azimuthalAngle());
00786     default:
00787       throw std::runtime_error("The specified deltaR scheme is not yet implemented");
00788     }
00789   }
00790 
00791   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between two four-vectors.
00792   /// There is a scheme ambiguity for momentum-type four vectors
00793   /// as to whether the pseudorapidity (a purely geometric concept) or the
00794   /// rapidity (a relativistic energy-momentum quantity) is to be used: this can
00795   /// be chosen via the optional scheme parameter.
00796   inline double deltaR(const FourVector& a, const FourMomentum& b,
00797                        RapScheme scheme = PSEUDORAPIDITY) {
00798     return deltaR(b, a, scheme);
00799   }
00800 
00801   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00802   /// three-vector and a four-vector.
00803   inline double deltaR(const FourMomentum& a, const Vector3& b) {
00804     return deltaR(a.vector3(), b);
00805   }
00806 
00807   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00808   /// three-vector and a four-vector.
00809   inline double deltaR(const Vector3& a, const FourMomentum& b) {
00810     return deltaR(a, b.vector3());
00811   }
00812 
00813   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00814   /// three-vector and a four-vector.
00815   inline double deltaR(const FourVector& a, const Vector3& b) {
00816     return deltaR(a.vector3(), b);
00817   }
00818 
00819   /// @brief Calculate the 2D rapidity-azimuthal ("eta-phi") distance between a
00820   /// three-vector and a four-vector.
00821   inline double deltaR(const Vector3& a, const FourVector& b) {
00822     return deltaR(a, b.vector3());
00823   }
00824 
00825   //@}
00826 
00827 
00828   //////////////////////////////////////////////////////
00829 
00830 
00831   /// @name \f$ \Delta phi \f$ calculations from 4-vectors
00832   //@{
00833 
00834   /// Calculate the difference in azimuthal angle between two vectors.
00835   inline double deltaPhi(const FourMomentum& a, const FourMomentum& b) {
00836     return deltaPhi(a.vector3(), b.vector3());
00837   }
00838 
00839   /// Calculate the difference in azimuthal angle between two vectors.
00840   inline double deltaPhi(const FourMomentum& v, double phi2) {
00841     return deltaPhi(v.vector3(), phi2);
00842   }
00843 
00844   /// Calculate the difference in azimuthal angle between two vectors.
00845   inline double deltaPhi(double phi1, const FourMomentum& v) {
00846     return deltaPhi(phi1, v.vector3());
00847   }
00848 
00849   /// Calculate the difference in azimuthal angle between two vectors.
00850   inline double deltaPhi(const FourVector& a, const FourVector& b) {
00851     return deltaPhi(a.vector3(), b.vector3());
00852   }
00853 
00854   /// Calculate the difference in azimuthal angle between two vectors.
00855   inline double deltaPhi(const FourVector& v, double phi2) {
00856     return deltaPhi(v.vector3(), phi2);
00857   }
00858 
00859   /// Calculate the difference in azimuthal angle between two vectors.
00860   inline double deltaPhi(double phi1, const FourVector& v) {
00861     return deltaPhi(phi1, v.vector3());
00862   }
00863 
00864   /// Calculate the difference in azimuthal angle between two vectors.
00865   inline double deltaPhi(const FourVector& a, const FourMomentum& b) {
00866     return deltaPhi(a.vector3(), b.vector3());
00867   }
00868 
00869   /// Calculate the difference in azimuthal angle between two vectors.
00870   inline double deltaPhi(const FourMomentum& a, const FourVector& b) {
00871     return deltaPhi(a.vector3(), b.vector3());
00872   }
00873 
00874   /// Calculate the difference in azimuthal angle between two vectors.
00875   inline double deltaPhi(const FourVector& a, const Vector3& b) {
00876     return deltaPhi(a.vector3(), b);
00877   }
00878 
00879   /// Calculate the difference in azimuthal angle between two vectors.
00880   inline double deltaPhi(const Vector3& a, const FourVector& b) {
00881     return deltaPhi(a, b.vector3());
00882   }
00883 
00884   /// Calculate the difference in azimuthal angle between two vectors.
00885   inline double deltaPhi(const FourMomentum& a, const Vector3& b) {
00886     return deltaPhi(a.vector3(), b);
00887   }
00888 
00889   /// Calculate the difference in azimuthal angle between two vectors.
00890   inline double deltaPhi(const Vector3& a, const FourMomentum& b) {
00891     return deltaPhi(a, b.vector3());
00892   }
00893 
00894   //@}
00895 
00896 
00897   //////////////////////////////////////////////////////
00898 
00899 
00900   /// @name \f$ |\Delta eta| \f$ calculations from 4-vectors
00901   //@{
00902 
00903   /// Calculate the difference in pseudorapidity between two vectors.
00904   inline double deltaEta(const FourMomentum& a, const FourMomentum& b) {
00905     return deltaEta(a.vector3(), b.vector3());
00906   }
00907 
00908   /// Calculate the difference in pseudorapidity between two vectors.
00909   inline double deltaEta(const FourMomentum& v, double eta2) {
00910     return deltaEta(v.vector3(), eta2);
00911   }
00912 
00913   /// Calculate the difference in pseudorapidity between two vectors.
00914   inline double deltaEta(double eta1, const FourMomentum& v) {
00915     return deltaEta(eta1, v.vector3());
00916   }
00917 
00918   /// Calculate the difference in pseudorapidity between two vectors.
00919   inline double deltaEta(const FourVector& a, const FourVector& b) {
00920     return deltaEta(a.vector3(), b.vector3());
00921   }
00922 
00923   /// Calculate the difference in pseudorapidity between two vectors.
00924   inline double deltaEta(const FourVector& v, double eta2) {
00925     return deltaEta(v.vector3(), eta2);
00926   }
00927 
00928   /// Calculate the difference in pseudorapidity between two vectors.
00929   inline double deltaEta(double eta1, const FourVector& v) {
00930     return deltaEta(eta1, v.vector3());
00931   }
00932 
00933   /// Calculate the difference in pseudorapidity between two vectors.
00934   inline double deltaEta(const FourVector& a, const FourMomentum& b) {
00935     return deltaEta(a.vector3(), b.vector3());
00936   }
00937 
00938   /// Calculate the difference in pseudorapidity between two vectors.
00939   inline double deltaEta(const FourMomentum& a, const FourVector& b) {
00940     return deltaEta(a.vector3(), b.vector3());
00941   }
00942 
00943   /// Calculate the difference in pseudorapidity between two vectors.
00944   inline double deltaEta(const FourVector& a, const Vector3& b) {
00945     return deltaEta(a.vector3(), b);
00946   }
00947 
00948   /// Calculate the difference in pseudorapidity between two vectors.
00949   inline double deltaEta(const Vector3& a, const FourVector& b) {
00950     return deltaEta(a, b.vector3());
00951   }
00952 
00953   /// Calculate the difference in pseudorapidity between two vectors.
00954   inline double deltaEta(const FourMomentum& a, const Vector3& b) {
00955     return deltaEta(a.vector3(), b);
00956   }
00957 
00958   /// Calculate the difference in pseudorapidity between two vectors.
00959   inline double deltaEta(const Vector3& a, const FourMomentum& b) {
00960     return deltaEta(a, b.vector3());
00961   }
00962 
00963   //@}
00964 
00965 
00966   /// @name \f$ |\Delta y| \f$ calculations from 4-momentum vectors
00967   //@{
00968 
00969   /// Calculate the difference in rapidity between two 4-momentum vectors.
00970   inline double deltaRap(const FourMomentum& a, const FourMomentum& b) {
00971     return deltaRap(a.rapidity(), b.rapidity());
00972   }
00973 
00974   /// Calculate the difference in rapidity between two 4-momentum vectors.
00975   inline double deltaRap(const FourMomentum& v, double y2) {
00976     return deltaRap(v.rapidity(), y2);
00977   }
00978 
00979   /// Calculate the difference in rapidity between two 4-momentum vectors.
00980   inline double deltaRap(double y1, const FourMomentum& v) {
00981     return deltaRap(y1, v.rapidity());
00982   }
00983 
00984   //@}
00985 
00986 
00987   //////////////////////////////////////////////////////
00988 
00989 
00990   /// @name 4-vector comparison functions (for sorting)
00991   //@{
00992 
00993   inline bool cmpMomByPt(const FourMomentum& a, const FourMomentum& b) {
00994     return a.pt() > b.pt();
00995   }
00996   inline bool cmpMomByAscPt(const FourMomentum& a, const FourMomentum& b) {
00997     return a.pt() < b.pt();
00998   }
00999   inline bool cmpMomByP(const FourMomentum& a, const FourMomentum& b) {
01000     return a.vector3().mod() > b.vector3().mod();
01001   }
01002   inline bool cmpMomByAscP(const FourMomentum& a, const FourMomentum& b) {
01003     return a.vector3().mod() < b.vector3().mod();
01004   }
01005   inline bool cmpMomByEt(const FourMomentum& a, const FourMomentum& b) {
01006     return a.Et() > b.Et();
01007   }
01008   inline bool cmpMomByAscEt(const FourMomentum& a, const FourMomentum& b) {
01009     return a.Et() < b.Et();
01010   }
01011   inline bool cmpMomByE(const FourMomentum& a, const FourMomentum& b) {
01012     return a.E() > b.E();
01013   }
01014   inline bool cmpMomByAscE(const FourMomentum& a, const FourMomentum& b) {
01015     return a.E() < b.E();
01016   }
01017   inline bool cmpMomByDescPseudorapidity(const FourMomentum& a, const FourMomentum& b) {
01018     return a.pseudorapidity() > b.pseudorapidity();
01019   }
01020   inline bool cmpMomByAscPseudorapidity(const FourMomentum& a, const FourMomentum& b) {
01021     return a.pseudorapidity() < b.pseudorapidity();
01022   }
01023   inline bool cmpMomByDescAbsPseudorapidity(const FourMomentum& a, const FourMomentum& b) {
01024     return fabs(a.pseudorapidity()) > fabs(b.pseudorapidity());
01025   }
01026   inline bool cmpMomByAscAbsPseudorapidity(const FourMomentum& a, const FourMomentum& b) {
01027     return fabs(a.pseudorapidity()) < fabs(b.pseudorapidity());
01028   }
01029   inline bool cmpMomByDescRapidity(const FourMomentum& a, const FourMomentum& b) {
01030     return a.rapidity() > b.rapidity();
01031   }
01032   inline bool cmpMomByAscRapidity(const FourMomentum& a, const FourMomentum& b) {
01033     return a.rapidity() < b.rapidity();
01034   }
01035   inline bool cmpMomByDescAbsRapidity(const FourMomentum& a, const FourMomentum& b) {
01036     return fabs(a.rapidity()) > fabs(b.rapidity());
01037   }
01038   inline bool cmpMomByAscAbsRapidity(const FourMomentum& a, const FourMomentum& b) {
01039     return fabs(a.rapidity()) < fabs(b.rapidity());
01040   }
01041 
01042   //@}
01043 
01044 
01045   //////////////////////////////////////////////////////
01046 
01047 
01048   /// @name 4-vector string representations
01049   //@{
01050 
01051   /// Render a 4-vector as a string.
01052   inline std::string toString(const FourVector& lv) {
01053     ostringstream out;
01054     out << "("  << (fabs(lv.t()) < 1E-30 ? 0.0 : lv.t())
01055         << "; " << (fabs(lv.x()) < 1E-30 ? 0.0 : lv.x())
01056         << ", " << (fabs(lv.y()) < 1E-30 ? 0.0 : lv.y())
01057         << ", " << (fabs(lv.z()) < 1E-30 ? 0.0 : lv.z())
01058         << ")";
01059     return out.str();
01060   }
01061 
01062   /// Write a 4-vector to an ostream.
01063   inline std::ostream& operator<<(std::ostream& out, const FourVector& lv) {
01064     out << toString(lv);
01065     return out;
01066   }
01067 
01068   //@}
01069 
01070 
01071 }
01072 
01073 #endif