Specialized version of the FourVector with momentum/energy functionality. More...
#include <Vector4.hh>
Classes | |
struct | byEAscending |
struct for sorting by increasing energy More... | |
struct | byEDescending |
struct for sorting by decreasing energy More... | |
Public Member Functions | |
FourMomentum () | |
template<typename V4 > | |
FourMomentum (const V4 &other) | |
FourMomentum (const Vector< 4 > &other) | |
FourMomentum (const double E, const double px, const double py, const double pz) | |
~FourMomentum () | |
double | E () const |
Get energy ![]() | |
Vector3 | p () const |
Get 3-momentum part, ![]() | |
double | px () const |
Get x-component of momentum ![]() | |
double | py () const |
Get y-component of momentum ![]() | |
double | pz () const |
Get z-component of momentum ![]() | |
FourMomentum & | setE (double E) |
Set energy ![]() | |
FourMomentum & | setPx (double px) |
Set x-component of momentum ![]() | |
FourMomentum & | setPy (double py) |
Set y-component of momentum ![]() | |
FourMomentum & | setPz (double pz) |
Set z-component of momentum ![]() | |
double | mass2 () const |
Get squared mass ![]() | |
double | mass () const |
Get mass ![]() | |
double | rapidity () const |
Calculate rapidity. | |
double | pT2 () const |
Calculate squared transverse momentum ![]() | |
double | pT () const |
Calculate transverse momentum ![]() | |
double | Et2 () const |
Calculate transverse energy ![]() | |
double | Et () const |
Calculate transverse energy ![]() | |
Vector3 | boostVector () const |
Calculate boost vector (in units of ![]() | |
double | t () const |
double | x () const |
double | y () const |
double | z () const |
FourVector & | setT (const double t) |
FourVector & | setX (const double x) |
FourVector & | setY (const double y) |
FourVector & | setZ (const double z) |
double | invariant () const |
double | angle (const FourVector &v) const |
double | angle (const Vector3 &v3) const |
double | polarRadius2 () const |
double | perp2 () const |
Synonym for polarRadius2. | |
double | rho2 () const |
Synonym for polarRadius2. | |
double | polarRadius () const |
double | perp () const |
Synonym for polarRadius. | |
double | rho () const |
Synonym for polarRadius. | |
double | azimuthalAngle (const PhiMapping mapping=ZERO_2PI) const |
Angle subtended by the 3-vector's projection in x-y and the x-axis. | |
double | phi (const PhiMapping mapping=ZERO_2PI) const |
Synonym for azimuthalAngle. | |
double | polarAngle () const |
Angle subtended by the 3-vector and the z-axis. | |
double | theta () const |
Synonym for polarAngle. | |
double | pseudorapidity () const |
double | eta () const |
Synonym for pseudorapidity. | |
Vector3 | vector3 () const |
Get the spatial part of the 4-vector as a 3-vector. | |
double | contract (const FourVector &v) const |
Contract two 4-vectors, with metric signature (+ - - -). | |
double | dot (const FourVector &v) const |
Contract two 4-vectors, with metric signature (+ - - -). | |
double | operator* (const FourVector &v) const |
Contract two 4-vectors, with metric signature (+ - - -). | |
FourVector & | operator*= (double a) |
Multiply by a scalar. | |
FourVector & | operator/= (double a) |
Divide by a scalar. | |
FourVector & | operator+= (const FourVector &v) |
FourVector & | operator-= (const FourVector &v) |
FourVector | operator- () const |
Invert the vector. | |
const double & | get (const size_t index) const |
const double & | operator[] (const size_t index) const |
Direct access to vector elements by index. | |
double & | operator[] (const size_t index) |
Direct access to vector elements by index. | |
Vector< N > & | set (const size_t index, const double value) |
Set indexed value. | |
size_t | size () const |
Vector dimensionality. | |
bool | isZero (double tolerance=1E-5) const |
Check for nullness, allowing for numerical precision. | |
double | mod2 () const |
Calculate the modulus-squared of a vector. ![]() | |
double | mod () const |
Calculate the modulus of a vector. ![]() | |
bool | operator== (const Vector< N > &a) const |
bool | operator!= (const Vector< N > &a) const |
bool | operator< (const Vector< N > &a) const |
bool | operator<= (const Vector< N > &a) const |
bool | operator> (const Vector< N > &a) const |
bool | operator>= (const Vector< N > &a) const |
Protected Member Functions | |
double & | get (const size_t index) |
Protected Attributes | |
Eigen::Vector< double, N > | _vec |
Vector. |
Specialized version of the FourVector with momentum/energy functionality.
Definition at line 310 of file Vector4.hh.
FourMomentum | ( | ) | [inline] |
Definition at line 312 of file Vector4.hh.
FourMomentum | ( | const V4 & | other | ) | [inline] |
Definition at line 315 of file Vector4.hh.
References FourMomentum::setE(), FourMomentum::setPx(), FourMomentum::setPy(), and FourMomentum::setPz().
FourMomentum | ( | const Vector< 4 > & | other | ) | [inline] |
Definition at line 322 of file Vector4.hh.
00323 : FourVector(other) { }
FourMomentum | ( | const double | E, | |
const double | px, | |||
const double | py, | |||
const double | pz | |||
) | [inline] |
Definition at line 325 of file Vector4.hh.
References FourMomentum::setE(), FourMomentum::setPx(), FourMomentum::setPy(), and FourMomentum::setPz().
~FourMomentum | ( | ) | [inline] |
Definition at line 332 of file Vector4.hh.
double angle | ( | const Vector3 & | v3 | ) | const [inline, inherited] |
Definition at line 68 of file Vector4.hh.
References Vector3::angle(), and FourVector::vector3().
00068 { 00069 return vector3().angle(v3); 00070 }
double angle | ( | const FourVector & | v | ) | const [inline, inherited] |
Definition at line 64 of file Vector4.hh.
References Vector3::angle(), and FourVector::vector3().
Referenced by D0_1996_S3214044::_fourJetAnalysis(), H1_2000_S4129130::analyze(), H1_1994_S2919893::analyze(), and Rivet::angle().
00064 { 00065 return vector3().angle( v.vector3() ); 00066 }
double azimuthalAngle | ( | const PhiMapping | mapping = ZERO_2PI |
) | const [inline, inherited] |
Angle subtended by the 3-vector's projection in x-y and the x-axis.
Definition at line 101 of file Vector4.hh.
References Vector3::azimuthalAngle(), and FourVector::vector3().
Referenced by CDF_2004_S5839831::_calcTransCones(), MC_LEADINGJETS::analyze(), H1_1994_S2919893::analyze(), D0_2008_S7863608::analyze(), D0_2008_S7719523::analyze(), D0_2008_S6879055::analyze(), D0_2006_S6438750::analyze(), CDF_2008_S8095620::analyze(), CDF_2008_S7540469::analyze(), CDF_2006_S6653332::analyze(), CDF_2004_S5839831::analyze(), Rivet::azimuthalAngle(), Rivet::deltaR(), DISKinematics::project(), and ClusteredPhotons::project().
00101 { 00102 return vector3().azimuthalAngle(mapping); 00103 }
Vector3 boostVector | ( | ) | const [inline] |
Calculate boost vector (in units of ).
Definition at line 403 of file Vector4.hh.
References FourMomentum::E(), FourMomentum::px(), FourMomentum::py(), and FourMomentum::pz().
Referenced by CDF_1997_S3541940::_avg_beam_in_lab(), CDF_1996_S3349578::_avg_beam_in_lab(), CDF_1996_S3349578::_fiveJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), D0_1996_S3214044::_threeJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), CDF_1997_S3541940::analyze(), CDF_1996_S3108457::analyze(), BELLE_2006_S6265367::analyze(), Rivet::boostVector(), and DISKinematics::project().
00403 { 00404 // const Vector3 p3 = vector3(); 00405 // const double m2 = mass2(); 00406 // if (Rivet::isZero(m2)) return p3.unit(); 00407 // else { 00408 // // Could also do this via beta = tanh(rapidity), but that's 00409 // // probably more messy from a numerical hygiene point of view. 00410 // const double p2 = p3.mod2(); 00411 // const double beta = sqrt( p2 / (m2 + p2) ); 00412 // return beta * p3.unit(); 00413 // } 00414 /// @todo Be careful about c=1 convention... 00415 return Vector3(px()/E(), py()/E(), pz()/E()); 00416 }
double contract | ( | const FourVector & | v | ) | const [inline, inherited] |
Contract two 4-vectors, with metric signature (+ - - -).
Definition at line 136 of file Vector4.hh.
References FourVector::t(), FourVector::x(), FourVector::y(), and FourVector::z().
Referenced by Rivet::contract(), FourVector::dot(), and FourVector::operator*().
double dot | ( | const FourVector & | v | ) | const [inline, inherited] |
Contract two 4-vectors, with metric signature (+ - - -).
Definition at line 142 of file Vector4.hh.
References FourVector::contract().
00142 { 00143 return contract(v); 00144 }
double E | ( | ) | const [inline] |
Get energy (time component of momentum).
Definition at line 336 of file Vector4.hh.
References FourVector::t().
Referenced by CDF_1997_S3541940::_avg_beam_in_lab(), CDF_1996_S3349578::_avg_beam_in_lab(), CDF_1996_S3349578::_fiveJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), D0_1996_S3214044::_threeJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), OPAL_1998_S3780481::analyze(), OPAL_1993_S2692198::analyze(), MC_PHOTONJETS::analyze(), MC_GENERIC::analyze(), H1_2000_S4129130::analyze(), H1_1994_S2919893::analyze(), D0_2008_S7719523::analyze(), D0_2006_S6438750::analyze(), CDF_1997_S3541940::analyze(), ALEPH_2004_S5765862::analyze(), ALEPH_1996_S3486095::analyze(), ALEPH_1996_S3196992::analyze(), FourMomentum::boostVector(), FastJets::calc(), Rivet::cmpMomByAscE(), Rivet::cmpMomByE(), Rivet::cmpParticleByAscE(), Rivet::cmpParticleByE(), FourMomentum::Et(), byEAscending::operator()(), ALEPH_1996_S3196992::particleInJet(), Hemispheres::project(), FourMomentum::rapidity(), and Rivet::sqrtS().
00336 { return t(); }
double Et | ( | ) | const [inline] |
Calculate transverse energy .
Definition at line 398 of file Vector4.hh.
References FourMomentum::E(), and FourVector::polarAngle().
Referenced by MC_WJETS::analyze(), MC_DIPHOTON::analyze(), D0_2010_S8570965::analyze(), D0_2008_S7837160::analyze(), D0_1998_S3711838::analyze(), D0_1996_S3324664::analyze(), CDF_2009_S8436959::analyze(), CDF_2009_S8233977::analyze(), CDF_2008_S7541902::analyze(), CDF_2008_LEADINGJETS::analyze(), CDF_2005_S6080774::analyze(), CDF_2001_S4517016::analyze(), CDF_1993_S2742446::analyze(), CDF_1991_S2313472::analyze(), Rivet::cmpMomByAscEt(), Rivet::cmpMomByEt(), Rivet::cmpParticleByAscEt(), Rivet::cmpParticleByEt(), Rivet::Et(), FourMomentum::Et2(), TotalVisibleMomentum::project(), NeutralFinalState::project(), and CentralEtHCM::project().
00398 { 00399 return E() * sin(polarAngle()); 00400 }
double Et2 | ( | ) | const [inline] |
Calculate transverse energy .
Definition at line 393 of file Vector4.hh.
References FourMomentum::Et().
Referenced by Rivet::Et2(), and byETAscending::operator()().
double eta | ( | ) | const [inline, inherited] |
Synonym for pseudorapidity.
Definition at line 125 of file Vector4.hh.
References Vector3::eta(), and FourVector::vector3().
Referenced by FinalState::accept(), UA1_1990_S2044935::analyze(), STAR_2008_S7993412::analyze(), MC_ZZJETS::analyze(), MC_ZJETS::analyze(), MC_WWJETS::analyze(), MC_WJETS::analyze(), MC_SUSY::analyze(), MC_PHOTONJETUE::analyze(), MC_PHOTONJETS::analyze(), MC_HJETS::analyze(), MC_GENERIC::analyze(), MC_DIPHOTON::analyze(), D0_2010_S8570965::analyze(), D0_2008_S7837160::analyze(), D0_2008_S7719523::analyze(), D0_1996_S3324664::analyze(), CDF_2009_S8436959::analyze(), CDF_2008_S7541902::analyze(), CDF_2008_NOTE_9351::analyze(), CDF_2005_S6080774::analyze(), CDF_2001_S4517016::analyze(), CDF_2000_S4266730::analyze(), CDF_1996_S3418421::analyze(), CDF_1994_S2952106::analyze(), CDF_1993_S2742446::analyze(), ATLAS_2010_S8591806::analyze(), ATLAS_2010_CONF_2010_031::analyze(), JetShape::calc(), Rivet::eta(), JetAlg::jets(), and NeutralFinalState::project().
00125 { 00126 return vector3().eta(); 00127 }
double& get | ( | const size_t | index | ) | [inline, protected, inherited] |
Definition at line 126 of file VectorN.hh.
References Vector< N >::_vec.
00126 { 00127 if (index >= N) { 00128 throw std::runtime_error("Tried to access an invalid vector index."); 00129 } else { 00130 return _vec(index); 00131 } 00132 }
const double& get | ( | const size_t | index | ) | const [inline, inherited] |
Definition at line 33 of file VectorN.hh.
References Vector< N >::_vec.
00033 { 00034 if (index >= N) { 00035 throw std::runtime_error("Tried to access an invalid vector index."); 00036 } else { 00037 return _vec(index); 00038 } 00039 }
double invariant | ( | ) | const [inline, inherited] |
Definition at line 59 of file Vector4.hh.
References FourVector::t(), FourVector::x(), FourVector::y(), and FourVector::z().
Referenced by D0_2001_S4674421::analyze(), Rivet::invariant(), and FourMomentum::mass2().
bool isZero | ( | double | tolerance = 1E-5 |
) | const [inline, inherited] |
Check for nullness, allowing for numerical precision.
Definition at line 67 of file VectorN.hh.
References Vector< N >::_vec, and Rivet::isZero().
00067 { 00068 for (size_t i=0; i < N; ++i) { 00069 if (! Rivet::isZero(_vec[i], tolerance) ) return false; 00070 } 00071 return true; 00072 }
double mass | ( | ) | const [inline] |
Get mass (the Lorentz self-invariant).
Definition at line 368 of file Vector4.hh.
References Rivet::isZero(), and FourMomentum::mass2().
Referenced by OPAL_1993_S2692198::analyze(), MC_DIPHOTON::analyze(), D0_2010_S8570965::analyze(), D0_1998_S3711838::analyze(), CDF_2008_NOTE_9351::analyze(), CDF_2005_S6080774::analyze(), CDF_2000_S4155203::analyze(), CDF_1996_S3108457::analyze(), CDF_1991_S2313472::analyze(), ATLAS_2010_S8817804::analyze(), Rivet::mass(), Particle::mass(), WFinder::project(), and InvMassFinalState::project().
00368 { 00369 assert(Rivet::isZero(mass2()) || mass2() > 0); 00370 if (Rivet::isZero(mass2())) { 00371 return 0.0; 00372 } else { 00373 return sqrt(mass2()); 00374 } 00375 }
double mass2 | ( | ) | const [inline] |
Get squared mass (the Lorentz self-invariant).
Definition at line 363 of file Vector4.hh.
References FourVector::invariant().
Referenced by D0_1996_S3214044::_safeMass(), CDF_1997_S3541940::_safeMass(), CDF_1996_S3349578::_safeMass(), CDF_2000_S4155203::analyze(), FourMomentum::mass(), Rivet::mass2(), InvMassFinalState::project(), Hemispheres::project(), and DISKinematics::project().
00363 { 00364 return invariant(); 00365 }
double mod | ( | ) | const [inline, inherited] |
Calculate the modulus of a vector. .
Definition at line 87 of file VectorN.hh.
References Vector< N >::mod2().
00087 { 00088 const double norm = mod2(); 00089 assert(norm >= 0); 00090 return sqrt(norm); 00091 }
double mod2 | ( | ) | const [inline, inherited] |
Calculate the modulus-squared of a vector. .
Definition at line 76 of file VectorN.hh.
References Vector< N >::mod2(), and Vector< N >::size().
bool operator!= | ( | const Vector< N > & | a | ) | const [inline, inherited] |
Definition at line 104 of file VectorN.hh.
References Vector< N >::_vec.
00104 { 00105 return _vec != a._vec; 00106 }
double operator* | ( | const FourVector & | v | ) | const [inline, inherited] |
Contract two 4-vectors, with metric signature (+ - - -).
Definition at line 147 of file Vector4.hh.
References FourVector::contract().
00147 { 00148 return contract(v); 00149 }
FourVector& operator*= | ( | double | a | ) | [inline, inherited] |
Multiply by a scalar.
Definition at line 152 of file Vector4.hh.
References Vector< N >::_vec, Vector< 4 >::_vec, and FourVector::multiply.
FourVector& operator+= | ( | const FourVector & | v | ) | [inline, inherited] |
Definition at line 163 of file Vector4.hh.
References Vector< N >::_vec, Vector< 4 >::_vec, and FourVector::add.
FourVector operator- | ( | ) | const [inline, inherited] |
Invert the vector.
Reimplemented from Vector< 4 >.
Definition at line 173 of file Vector4.hh.
References Vector< 4 >::_vec, and Vector< N >::_vec.
00173 { 00174 FourVector result; 00175 result._vec = -_vec; 00176 return result; 00177 }
FourVector& operator-= | ( | const FourVector & | v | ) | [inline, inherited] |
Definition at line 168 of file Vector4.hh.
References Vector< N >::_vec, Vector< 4 >::_vec, and FourVector::add.
FourVector& operator/= | ( | double | a | ) | [inline, inherited] |
Divide by a scalar.
Definition at line 158 of file Vector4.hh.
References Vector< N >::_vec, Vector< 4 >::_vec, and FourVector::multiply.
bool operator< | ( | const Vector< N > & | a | ) | const [inline, inherited] |
Definition at line 108 of file VectorN.hh.
References Vector< N >::_vec.
00108 { 00109 return _vec < a._vec; 00110 }
bool operator<= | ( | const Vector< N > & | a | ) | const [inline, inherited] |
Definition at line 112 of file VectorN.hh.
References Vector< N >::_vec.
00112 { 00113 return _vec <= a._vec; 00114 }
bool operator== | ( | const Vector< N > & | a | ) | const [inline, inherited] |
Definition at line 100 of file VectorN.hh.
References Vector< N >::_vec.
00100 { 00101 return _vec == a._vec; 00102 }
bool operator> | ( | const Vector< N > & | a | ) | const [inline, inherited] |
Definition at line 116 of file VectorN.hh.
References Vector< N >::_vec.
00116 { 00117 return _vec > a._vec; 00118 }
bool operator>= | ( | const Vector< N > & | a | ) | const [inline, inherited] |
Definition at line 120 of file VectorN.hh.
References Vector< N >::_vec.
00120 { 00121 return _vec >= a._vec; 00122 }
double& operator[] | ( | const size_t | index | ) | [inline, inherited] |
Direct access to vector elements by index.
Definition at line 47 of file VectorN.hh.
const double& operator[] | ( | const size_t | index | ) | const [inline, inherited] |
Direct access to vector elements by index.
Definition at line 42 of file VectorN.hh.
Vector3 p | ( | ) | const [inline] |
Get 3-momentum part, .
Definition at line 339 of file Vector4.hh.
References FourVector::vector3().
Referenced by ATLAS_2010_CONF_2010_049::analyze().
00339 { return vector3(); }
double perp | ( | ) | const [inline, inherited] |
Synonym for polarRadius.
Definition at line 91 of file Vector4.hh.
References Vector3::perp(), and FourVector::vector3().
Referenced by MC_TTBAR::analyze(), LHCB_2010_S8758301::analyze(), ATLAS_2010_CONF_2010_081::analyze(), and Rivet::perp().
00091 { 00092 return vector3().perp(); 00093 }
double perp2 | ( | ) | const [inline, inherited] |
Synonym for polarRadius2.
Definition at line 77 of file Vector4.hh.
References Vector3::perp2(), and FourVector::vector3().
Referenced by Rivet::perp2().
00077 { 00078 return vector3().perp2(); 00079 }
double phi | ( | const PhiMapping | mapping = ZERO_2PI |
) | const [inline, inherited] |
Synonym for azimuthalAngle.
Definition at line 106 of file Vector4.hh.
References Vector3::phi(), and FourVector::vector3().
Referenced by STAR_2009_UE_HELEN::analyze(), STAR_2008_S7993412::analyze(), MC_ZZJETS::analyze(), MC_WWJETS::analyze(), MC_TTBAR::analyze(), MC_SUSY::analyze(), MC_PHOTONJETUE::analyze(), MC_PHOTONJETS::analyze(), MC_LEADINGJETS::analyze(), MC_GENERIC::analyze(), MC_DIPHOTON::analyze(), D0_2010_S8570965::analyze(), D0_2009_S8349509::analyze(), D0_2008_S7719523::analyze(), D0_1996_S3324664::analyze(), CDF_2009_S8436959::analyze(), CDF_2008_LEADINGJETS::analyze(), CDF_2005_S6080774::analyze(), CDF_2001_S4751469::analyze(), CDF_1994_S2952106::analyze(), CDF_1993_S2742446::analyze(), ATLAS_2010_CONF_2010_081::analyze(), and Rivet::phi().
00106 { 00107 return vector3().phi(mapping); 00108 }
double polarAngle | ( | ) | const [inline, inherited] |
Angle subtended by the 3-vector and the z-axis.
Definition at line 111 of file Vector4.hh.
References Vector3::polarAngle(), and FourVector::vector3().
Referenced by H1_1994_S2919893::beamAngle(), FourMomentum::Et(), Rivet::polarAngle(), and DISKinematics::project().
00111 { 00112 return vector3().polarAngle(); 00113 }
double polarRadius | ( | ) | const [inline, inherited] |
Definition at line 86 of file Vector4.hh.
References Vector3::polarRadius(), and FourVector::vector3().
Referenced by Rivet::polarRadius().
00086 { 00087 return vector3().polarRadius(); 00088 }
double polarRadius2 | ( | ) | const [inline, inherited] |
Definition at line 72 of file Vector4.hh.
References Vector3::polarRadius2(), and FourVector::vector3().
Referenced by Rivet::polarRadius2().
00072 { 00073 return vector3().polarRadius2(); 00074 }
double pseudorapidity | ( | ) | const [inline, inherited] |
Definition at line 120 of file Vector4.hh.
References Vector3::pseudorapidity(), and FourVector::vector3().
Referenced by CDF_2004_S5839831::_calcTransCones(), UA5_1986_S1583476::analyze(), UA5_1982_S875503::analyze(), SFM_1984_S1178091::analyze(), MC_TTBAR::analyze(), H1_2000_S4129130::analyze(), H1_1994_S2919893::analyze(), D0_2008_S7863608::analyze(), D0_2008_S7719523::analyze(), D0_2008_S6879055::analyze(), D0_2006_S6438750::analyze(), CDF_2008_S7540469::analyze(), CDF_2004_S5839831::analyze(), CDF_1990_S2089246::analyze(), ALICE_2010_S8625980::analyze(), Rivet::cmpMomByAscAbsPseudorapidity(), Rivet::cmpMomByAscPseudorapidity(), Rivet::cmpMomByDescAbsPseudorapidity(), Rivet::cmpMomByDescPseudorapidity(), Rivet::cmpParticleByAscAbsPseudorapidity(), Rivet::cmpParticleByAscPseudorapidity(), Rivet::cmpParticleByDescAbsPseudorapidity(), Rivet::cmpParticleByDescPseudorapidity(), TriggerUA5::project(), TriggerCDFRun2::project(), TriggerCDFRun0Run1::project(), SVertex::project(), ClusteredPhotons::project(), and Rivet::pseudorapidity().
00120 { 00121 return vector3().pseudorapidity(); 00122 }
double pT | ( | ) | const [inline] |
Calculate transverse momentum .
Definition at line 388 of file Vector4.hh.
References FourMomentum::pT2().
Referenced by SVertex::_applyVtxTrackCuts(), CDF_2004_S5839831::_calcTransCones(), FinalState::accept(), UA1_1990_S2044935::analyze(), STAR_2009_UE_HELEN::analyze(), STAR_2008_S7993412::analyze(), STAR_2008_S7869363::analyze(), STAR_2006_S6870392::analyze(), STAR_2006_S6860818::analyze(), STAR_2006_S6500200::analyze(), MC_ZZJETS::analyze(), MC_ZJETS::analyze(), MC_WWJETS::analyze(), MC_WJETS::analyze(), MC_SUSY::analyze(), MC_PHOTONJETUE::analyze(), MC_PHOTONJETS::analyze(), MC_LEADINGJETS::analyze(), MC_HJETS::analyze(), MC_GENERIC::analyze(), MC_DIPHOTON::analyze(), MC_DIJET::analyze(), D0_2010_S8671338::analyze(), D0_2010_S8570965::analyze(), D0_2009_S8349509::analyze(), D0_2008_S7863608::analyze(), D0_2008_S7719523::analyze(), D0_2007_S7075677::analyze(), D0_2006_S6438750::analyze(), D0_2001_S4674421::analyze(), D0_1998_S3711838::analyze(), CDF_2009_S8233977::analyze(), CDF_2008_S8095620::analyze(), CDF_2008_S7782535::analyze(), CDF_2008_S7541902::analyze(), CDF_2008_NOTE_9351::analyze(), CDF_2008_LEADINGJETS::analyze(), CDF_2005_S6080774::analyze(), CDF_2004_S5839831::analyze(), CDF_2002_S4796047::analyze(), CDF_2001_S4751469::analyze(), CDF_2000_S4155203::analyze(), CDF_1994_S2952106::analyze(), CDF_1991_S2313472::analyze(), CDF_1988_S1865951::analyze(), ATLAS_2010_S8591806::analyze(), ATLAS_2010_CONF_2010_081::analyze(), ATLAS_2010_CONF_2010_031::analyze(), ALICE_2010_S8706239::analyze(), JetShape::calc(), Rivet::cmpMomByAscPt(), Rivet::cmpMomByPt(), Rivet::cmpParticleByAscPt(), Rivet::cmpParticleByPt(), JetAlg::jets(), STARRandomFilter::operator()(), VetoedFinalState::project(), ClosestJetShape::project(), and Rivet::pT().
00388 { 00389 return sqrt(pT2()); 00390 }
double pT2 | ( | ) | const [inline] |
Calculate squared transverse momentum .
Definition at line 383 of file Vector4.hh.
References Vector3::polarRadius2(), and FourVector::vector3().
Referenced by byPTAscending::operator()(), FourMomentum::pT(), and Rivet::pT2().
00383 { 00384 return vector3().polarRadius2(); 00385 }
double px | ( | ) | const [inline] |
Get x-component of momentum .
Definition at line 342 of file Vector4.hh.
References FourVector::x().
Referenced by OPAL_1993_S2692198::analyze(), CDF_2008_S7541902::analyze(), FourMomentum::boostVector(), FastJets::calc(), Rivet::get2dDecayLength(), Rivet::get3dDecayLength(), and DISKinematics::project().
00342 { return x(); }
double py | ( | ) | const [inline] |
Get y-component of momentum .
Definition at line 345 of file Vector4.hh.
References FourVector::y().
Referenced by OPAL_1993_S2692198::analyze(), CDF_2008_S7541902::analyze(), FourMomentum::boostVector(), FastJets::calc(), Rivet::get2dDecayLength(), and Rivet::get3dDecayLength().
00345 { return y(); }
double pz | ( | ) | const [inline] |
Get z-component of momentum .
Definition at line 348 of file Vector4.hh.
References FourVector::z().
Referenced by SFM_1984_S1178091::analyze(), OPAL_1993_S2692198::analyze(), FourMomentum::boostVector(), FastJets::calc(), Rivet::get3dDecayLength(), DISLepton::project(), FourMomentum::rapidity(), and Rivet::sqrtS().
00348 { return z(); }
double rapidity | ( | ) | const [inline] |
Calculate rapidity.
Definition at line 378 of file Vector4.hh.
References FourMomentum::E(), and FourMomentum::pz().
Referenced by CDF_1996_S3349578::_fiveJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), STAR_2008_S7869363::analyze(), STAR_2006_S6860818::analyze(), STAR_2006_S6500200::analyze(), MC_PHOTONJETS::analyze(), MC_GENERIC::analyze(), LHCB_2010_S8758301::analyze(), D0_2010_S8566488::analyze(), D0_2009_S8349509::analyze(), D0_2009_S8320160::analyze(), D0_2008_S7863608::analyze(), D0_2008_S7719523::analyze(), CDF_2008_S8095620::analyze(), CDF_2008_S8093652::analyze(), CDF_2008_S7541902::analyze(), CDF_2006_S6653332::analyze(), CDF_1997_S3541940::analyze(), ATLAS_2010_S8817804::analyze(), JetShape::calc(), Rivet::cmpMomByAscAbsRapidity(), Rivet::cmpMomByAscRapidity(), Rivet::cmpMomByDescAbsRapidity(), Rivet::cmpMomByDescRapidity(), Rivet::cmpParticleByAscAbsRapidity(), Rivet::cmpParticleByAscRapidity(), Rivet::cmpParticleByDescAbsRapidity(), Rivet::cmpParticleByDescRapidity(), Rivet::deltaR(), JetAlg::jets(), CentralEtHCM::project(), and Rivet::rapidity().
double rho | ( | ) | const [inline, inherited] |
Synonym for polarRadius.
Definition at line 96 of file Vector4.hh.
References Vector3::rho(), and FourVector::vector3().
Referenced by Rivet::rho().
00096 { 00097 return vector3().rho(); 00098 }
double rho2 | ( | ) | const [inline, inherited] |
Synonym for polarRadius2.
Definition at line 82 of file Vector4.hh.
References Vector3::rho2(), and FourVector::vector3().
Referenced by Rivet::rho2().
00082 { 00083 return vector3().rho2(); 00084 }
Vector<N>& set | ( | const size_t | index, | |
const double | value | |||
) | [inline, inherited] |
Set indexed value.
Definition at line 52 of file VectorN.hh.
References Vector< N >::_vec.
00052 { 00053 if (index >= N) { 00054 throw std::runtime_error("Tried to access an invalid vector index."); 00055 } else { 00056 _vec[index] = value; 00057 } 00058 return *this; 00059 }
FourMomentum& setE | ( | double | E | ) | [inline] |
Set energy (time component of momentum).
Definition at line 351 of file Vector4.hh.
References FourVector::setT().
Referenced by FourMomentum::FourMomentum().
FourMomentum& setPx | ( | double | px | ) | [inline] |
Set x-component of momentum .
Definition at line 354 of file Vector4.hh.
References FourVector::setX().
Referenced by FourMomentum::FourMomentum().
FourMomentum& setPy | ( | double | py | ) | [inline] |
Set y-component of momentum .
Definition at line 357 of file Vector4.hh.
References FourVector::setY().
Referenced by FourMomentum::FourMomentum().
FourMomentum& setPz | ( | double | pz | ) | [inline] |
Set z-component of momentum .
Definition at line 360 of file Vector4.hh.
References FourVector::setZ().
Referenced by FourMomentum::FourMomentum().
FourVector& setT | ( | const double | t | ) | [inline, inherited] |
Definition at line 54 of file Vector4.hh.
Referenced by FourVector::FourVector(), and FourMomentum::setE().
00054 { set(0, t); return *this; }
FourVector& setX | ( | const double | x | ) | [inline, inherited] |
Definition at line 55 of file Vector4.hh.
Referenced by FourVector::FourVector(), and FourMomentum::setPx().
00055 { set(1, x); return *this; }
FourVector& setY | ( | const double | y | ) | [inline, inherited] |
Definition at line 56 of file Vector4.hh.
Referenced by FourVector::FourVector(), and FourMomentum::setPy().
00056 { set(2, y); return *this; }
FourVector& setZ | ( | const double | z | ) | [inline, inherited] |
Definition at line 57 of file Vector4.hh.
Referenced by FourVector::FourVector(), and FourMomentum::setPz().
00057 { set(3, z); return *this; }
size_t size | ( | ) | const [inline, inherited] |
double t | ( | ) | const [inline, inherited] |
Definition at line 50 of file Vector4.hh.
Referenced by FourVector::contract(), FourMomentum::E(), FourVector::invariant(), and Rivet::toString().
double theta | ( | ) | const [inline, inherited] |
Synonym for polarAngle.
Definition at line 116 of file Vector4.hh.
References Vector3::theta(), and FourVector::vector3().
Referenced by D0_1996_S3214044::_fourJetAnalysis(), ALEPH_1996_S3196992::analyze(), and Rivet::theta().
00116 { 00117 return vector3().theta(); 00118 }
Vector3 vector3 | ( | ) | const [inline, inherited] |
Get the spatial part of the 4-vector as a 3-vector.
Definition at line 130 of file Vector4.hh.
Referenced by CDF_1996_S3349578::_fiveJetAnalysis(), D0_1996_S3214044::_fourJetAnalysis(), CDF_1996_S3349578::_fourJetAnalysis(), CDF_1997_S3541940::_psi(), CDF_1996_S3349578::_psi(), D0_1996_S3214044::_threeJetAnalysis(), CDF_1996_S3349578::_threeJetAnalysis(), OPAL_1998_S3780481::analyze(), OPAL_1993_S2692198::analyze(), MC_ZZJETS::analyze(), MC_WWJETS::analyze(), DELPHI_1995_S3137023::analyze(), CDF_1997_S3541940::analyze(), ALEPH_2004_S5765862::analyze(), ALEPH_1996_S3486095::analyze(), Rivet::angle(), FourVector::angle(), FourVector::azimuthalAngle(), Thrust::calc(), Sphericity::calc(), Rivet::deltaR(), FourVector::eta(), FourMomentum::p(), FourVector::perp(), FourVector::perp2(), FourVector::phi(), FourVector::polarAngle(), FourVector::polarRadius(), FourVector::polarRadius2(), Hemispheres::project(), FoxWolframMoments::project(), DISKinematics::project(), FourVector::pseudorapidity(), FourMomentum::pT2(), FourVector::rho(), FourVector::rho2(), and FourVector::theta().
double x | ( | ) | const [inline, inherited] |
Definition at line 51 of file Vector4.hh.
Referenced by FourVector::contract(), FourVector::invariant(), ALEPH_1996_S3196992::particleInJet(), FourMomentum::px(), and Rivet::toString().
double y | ( | ) | const [inline, inherited] |
Definition at line 52 of file Vector4.hh.
Referenced by FourVector::contract(), FourVector::invariant(), ALEPH_1996_S3196992::particleInJet(), FourMomentum::py(), and Rivet::toString().
double z | ( | ) | const [inline, inherited] |
Definition at line 53 of file Vector4.hh.
Referenced by H1_1994_S2919893::analyze(), FourVector::contract(), FourVector::invariant(), ALEPH_1996_S3196992::particleInJet(), FourMomentum::pz(), and Rivet::toString().
Eigen::Vector<double,N> _vec [protected, inherited] |
Vector.
Definition at line 135 of file VectorN.hh.
Referenced by FourVector::operator*=(), FourVector::operator+=(), FourVector::operator-(), FourVector::operator-=(), and FourVector::operator/=().