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