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