rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.0
Rivet::FourMomentum Class Reference

Specialized version of the FourVector with momentum/energy functionality. More...

#include <Vector4.hh>

Inheritance diagram for Rivet::FourMomentum:
Rivet::FourVector Rivet::Vector< 4 >

Public Types

using EVector = RivetEigen::Matrix< double, N, 1 >
 Vector.
 

Public Member Functions

template<typename V4TYPE , typename std::enable_if< HasXYZT< V4TYPE >::value, int >::type DUMMY = 0>
 FourMomentum (const V4TYPE &other)
 
 FourMomentum (const Vector< 4 > &other)
 
 FourMomentum (const double E, const double px, const double py, const double pz)
 
 operator fastjet::PseudoJet () const
 Cast operator to FastJet PseudoJet.
 
double t () const
 
double t2 () const
 
FourVectorsetT (const double t)
 
double x () const
 
double x2 () const
 
FourVectorsetX (const double x)
 
double y () const
 
double y2 () const
 
FourVectorsetY (const double y)
 
double z () const
 
double z2 () const
 
FourVectorsetZ (const double z)
 
double invariant () const
 
bool isNull () const
 
double angle (const FourVector &v) const
 Angle between this vector and another.
 
double angle (const Vector3 &v3) const
 Angle between this vector and another (3-vector)
 
double polarRadius2 () const
 Mod-square of the projection of the 3-vector on to the \( x-y \) plane This is a more efficient function than polarRadius, as it avoids the square root. Use it if you only need the squared value, or e.g. an ordering by magnitude.
 
double perp2 () const
 Synonym for polarRadius2.
 
double rho2 () const
 Synonym for polarRadius2.
 
double polarRadius () const
 Magnitude of projection of 3-vector on to the \( x-y \) plane.
 
double perp () const
 Synonym for polarRadius.
 
double rho () const
 Synonym for polarRadius.
 
Vector3 polarVec () const
 Projection of 3-vector on to the \( x-y \) plane.
 
Vector3 perpVec () const
 Synonym for polarVec.
 
Vector3 rhoVec () const
 Synonym for polarVec.
 
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
 Pseudorapidity (defined purely by the 3-vector components)
 
double eta () const
 Synonym for pseudorapidity.
 
double abspseudorapidity () const
 Get the \( |\eta| \) directly.
 
double abseta () const
 Get the \( |\eta| \) directly (alias).
 
Vector3 vector3 () const
 Get the spatial part of the 4-vector as a 3-vector.
 
 operator Vector3 () const
 Implicit cast to 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 (+ - - -).
 
FourVectoroperator+= (const FourVector &v)
 Add to this 4-vector.
 
FourVectoroperator-= (const FourVector &v)
 Subtract from this 4-vector. NB time as well as space components are subtracted.
 
const doubleget (const size_t index) const
 
doubleget (const size_t index)
 
const doubleoperator[] (const size_t index) const
 Direct access to vector elements by index.
 
doubleoperator[] (const size_t index)
 Direct access to vector elements by index.
 
Vector< N > & set (const size_t index, const double value)
 Set indexed value.
 
constexpr 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. \( \sum_{i=1}^N x_i^2 \).
 
double mod () const
 Calculate the modulus of a vector. \( \sqrt{\sum_{i=1}^N x_i^2} \).
 
bool operator== (const Vector< N > &a) const
 
bool operator!= (const Vector< N > &a) const
 
Coordinate setters
FourMomentumsetE (double E)
 Set energy \( E \) (time component of momentum).
 
FourMomentumsetPx (double px)
 Set x-component of momentum \( p_x \).
 
FourMomentumsetPy (double py)
 Set y-component of momentum \( p_y \).
 
FourMomentumsetPz (double pz)
 Set z-component of momentum \( p_z \).
 
FourMomentumsetPE (double px, double py, double pz, double E)
 Set the p coordinates and energy simultaneously.
 
FourMomentumsetXYZE (double px, double py, double pz, double E)
 Alias for setPE.
 
FourMomentumsetPM (double px, double py, double pz, double mass)
 Set the p coordinates and mass simultaneously.
 
FourMomentumsetXYZM (double px, double py, double pz, double mass)
 Alias for setPM.
 
FourMomentumsetEtaPhiME (double eta, double phi, double mass, double E)
 
FourMomentumsetEtaPhiMPt (double eta, double phi, double mass, double pt)
 
FourMomentumsetRapPhiME (double y, double phi, double mass, double E)
 
FourMomentumsetRapPhiMPt (double y, double phi, double mass, double pt)
 
FourMomentumsetThetaPhiME (double theta, double phi, double mass, double E)
 
FourMomentumsetThetaPhiMPt (double theta, double phi, double mass, double pt)
 
FourMomentumsetPtPhiME (double pt, double phi, double mass, double E)
 
Accessors
double E () const
 Get energy \( E \) (time component of momentum).
 
double E2 () const
 Get energy-squared \( E^2 \).
 
double px () const
 Get x-component of momentum \( p_x \).
 
double px2 () const
 Get x-squared \( p_x^2 \).
 
double py () const
 Get y-component of momentum \( p_y \).
 
double py2 () const
 Get y-squared \( p_y^2 \).
 
double pz () const
 Get z-component of momentum \( p_z \).
 
double pz2 () const
 Get z-squared \( p_z^2 \).
 
double mass () const
 Get the mass \( m = \sqrt{E^2 - p^2} \) (the Lorentz self-invariant).
 
double mass2 () const
 Get the squared mass \( m^2 = E^2 - p^2 \) (the Lorentz self-invariant).
 
Vector3 p3 () const
 Get 3-momentum part, \( p \).
 
double p () const
 Get the modulus of the 3-momentum.
 
double p2 () const
 Get the modulus-squared of the 3-momentum.
 
double rapidity () const
 Calculate the rapidity.
 
double rap () const
 Alias for rapidity.
 
double absrapidity () const
 Absolute rapidity.
 
double absrap () const
 Absolute rapidity.
 
Vector3 pTvec () const
 Calculate the transverse momentum vector \( \vec{p}_T \).
 
Vector3 ptvec () const
 Synonym for pTvec.
 
double pT2 () const
 Calculate the squared transverse momentum \( p_T^2 \).
 
double pt2 () const
 Calculate the squared transverse momentum \( p_T^2 \).
 
double pT () const
 Calculate the transverse momentum \( p_T \).
 
double pt () const
 Calculate the transverse momentum \( p_T \).
 
double Et2 () const
 Calculate the transverse energy \( E_T^2 = E^2 \sin^2{\theta} \).
 
double Et () const
 Calculate the transverse energy \( E_T = E \sin{\theta} \).
 
Lorentz boost factors and vectors
double gamma () const
 
Vector3 gammaVec () const
 
double beta () const
 
Vector3 betaVec () const
 
Arithmetic operators
FourMomentumoperator*= (double a)
 Multiply by a scalar.
 
FourMomentumoperator/= (double a)
 Divide by a scalar.
 
FourMomentumoperator+= (const FourMomentum &v)
 Add to this 4-vector. NB time as well as space components are added.
 
FourMomentumoperator-= (const FourMomentum &v)
 Subtract from this 4-vector. NB time as well as space components are subtracted.
 
FourMomentum operator- () const
 Multiply all components (time and space) by -1.
 
FourMomentum reverse () const
 Multiply space components only by -1.
 

Static Public Member Functions

Factory functions
static FourMomentum mkXYZE (double px, double py, double pz, double E)
 Make a vector from (px,py,pz,E) coordinates.
 
static FourMomentum mkXYZM (double px, double py, double pz, double mass)
 Make a vector from (px,py,pz) coordinates and the mass.
 
static FourMomentum mkEtaPhiME (double eta, double phi, double mass, double E)
 Make a vector from (eta,phi,energy) coordinates and the mass.
 
static FourMomentum mkEtaPhiMPt (double eta, double phi, double mass, double pt)
 Make a vector from (eta,phi,pT) coordinates and the mass.
 
static FourMomentum mkRapPhiME (double y, double phi, double mass, double E)
 Make a vector from (y,phi,energy) coordinates and the mass.
 
static FourMomentum mkRapPhiMPt (double y, double phi, double mass, double pt)
 Make a vector from (y,phi,pT) coordinates and the mass.
 
static FourMomentum mkThetaPhiME (double theta, double phi, double mass, double E)
 Make a vector from (theta,phi,energy) coordinates and the mass.
 
static FourMomentum mkThetaPhiMPt (double theta, double phi, double mass, double pt)
 Make a vector from (theta,phi,pT) coordinates and the mass.
 
static FourMomentum mkPtPhiME (double pt, double phi, double mass, double E)
 Make a vector from (pT,phi,energy) coordinates and the mass.
 

Friends

FourMomentum multiply (const double a, const FourMomentum &v)
 
FourMomentum multiply (const FourMomentum &v, const double a)
 
FourMomentum add (const FourMomentum &a, const FourMomentum &b)
 
FourMomentum transform (const LorentzTransform &lt, const FourMomentum &v4)
 

Detailed Description

Specialized version of the FourVector with momentum/energy functionality.

Member Function Documentation

◆ beta()

double Rivet::FourMomentum::beta ( ) const
inline

Calculate the boost factor \( \beta \).

Note
\( \beta = pc/E \) so we rely on the c=1 convention

References E(), and p().

◆ betaVec()

Vector3 Rivet::FourMomentum::betaVec ( ) const
inline

Calculate the boost vector \( \vec{\beta} \).

Note
\( \beta = pc/E \) so we rely on the c=1 convention

References E(), and p3().

Referenced by Rivet::LorentzTransform::mkFrameTransform(), and Rivet::LorentzTransform::mkObjTransform().

◆ gamma()

double Rivet::FourMomentum::gamma ( ) const
inline

Calculate the boost factor \( \gamma \).

Note
\( \gamma = E/mc^2 \) so we rely on the c=1 convention

References E2(), Rivet::FourVector::eta(), and mass2().

Referenced by gammaVec().

◆ gammaVec()

Vector3 Rivet::FourMomentum::gammaVec ( ) const
inline

Calculate the boost vector \( \vec{\gamma} \).

Note
\( \gamma = E/mc^2 \) so we rely on the c=1 convention

References gamma(), p3(), and Rivet::Vector3::unit().

◆ mass()

double Rivet::FourMomentum::mass ( ) const
inline

◆ operator fastjet::PseudoJet()

Rivet::FourVector::operator fastjet::PseudoJet ( ) const
inherited

Cast operator to FastJet PseudoJet.

Needed, since otherwise the PseudoJet template constructor assumes the indices [0-3] mean px,py,pz,E... but Rivet uses E,px,py,pz ordering.

◆ setEtaPhiME()

FourMomentum & Rivet::FourMomentum::setEtaPhiME ( double  eta,
double  phi,
double  mass,
double  E 
)
inline

Set the vector state from (eta,phi,energy) coordinates and the mass

eta = -ln(tan(theta/2)) -> theta = 2 atan(exp(-eta))

References E(), Rivet::FourVector::eta(), mass(), Rivet::FourVector::phi(), setThetaPhiME(), and Rivet::FourVector::theta().

Referenced by mkEtaPhiME().

◆ setEtaPhiMPt()

FourMomentum & Rivet::FourMomentum::setEtaPhiMPt ( double  eta,
double  phi,
double  mass,
double  pt 
)
inline

Set the vector state from (eta,phi,pT) coordinates and the mass

eta = -ln(tan(theta/2)) -> theta = 2 atan(exp(-eta))

References E(), Rivet::FourVector::eta(), mass(), p(), Rivet::FourVector::phi(), pt(), setThetaPhiME(), Rivet::sqr(), and Rivet::FourVector::theta().

Referenced by mkEtaPhiMPt().

◆ setPtPhiME()

FourMomentum & Rivet::FourMomentum::setPtPhiME ( double  pt,
double  phi,
double  mass,
double  E 
)
inline

Set the vector state from (pT,phi,energy) coordinates and the mass

pz = sqrt(E^2 - mass^2 - pt^2)

References E(), Rivet::FourVector::eta(), mass(), Rivet::FourVector::phi(), pt(), px(), py(), pz(), setPE(), and Rivet::sqr().

Referenced by mkPtPhiME().

◆ setRapPhiME()

FourMomentum & Rivet::FourMomentum::setRapPhiME ( double  y,
double  phi,
double  mass,
double  E 
)
inline

Set the vector state from (y,phi,energy) coordinates and the mass

y = 0.5 * ln((E+pz)/(E-pz)) -> (E^2 - pz^2) exp(2y) = (E+pz)^2 & (E^2 - pz^2) exp(-2y) = (E-pz)^2 -> E = sqrt(pt^2 + m^2) cosh(y) -> pz = sqrt(pt^2 + m^2) sinh(y) -> sqrt(pt^2 + m^2) = E / cosh(y)

References E(), Rivet::FourVector::eta(), mass(), Rivet::FourVector::phi(), pt(), px(), py(), pz(), setPE(), and Rivet::sqr().

Referenced by mkRapPhiME(), and setRapPhiMPt().

◆ setRapPhiMPt()

FourMomentum & Rivet::FourMomentum::setRapPhiMPt ( double  y,
double  phi,
double  mass,
double  pt 
)
inline

Set the vector state from (y,phi,pT) coordinates and the mass

y = 0.5 * ln((E+pz)/(E-pz)) -> E = sqrt(pt^2 + m^2) cosh(y) [see above]

References E(), Rivet::FourVector::eta(), mass(), Rivet::FourVector::phi(), pt(), setRapPhiME(), and Rivet::sqr().

Referenced by mkRapPhiMPt().

◆ setThetaPhiME()

FourMomentum & Rivet::FourMomentum::setThetaPhiME ( double  theta,
double  phi,
double  mass,
double  E 
)
inline

Set the vector state from (theta,phi,energy) coordinates and the mass

p = sqrt(E^2 - mass^2) pz = p cos(theta) pt = p sin(theta)

References E(), Rivet::FourVector::eta(), mass(), p(), Rivet::FourVector::phi(), pt(), px(), py(), pz(), setPE(), Rivet::sqr(), and Rivet::FourVector::theta().

Referenced by mkThetaPhiME(), setEtaPhiME(), and setEtaPhiMPt().

◆ setThetaPhiMPt()

FourMomentum & Rivet::FourMomentum::setThetaPhiMPt ( double  theta,
double  phi,
double  mass,
double  pt 
)
inline

Set the vector state from (theta,phi,pT) coordinates and the mass

p = pt / sin(theta) pz = p cos(theta) E = sqrt(p^2 + mass^2)

References E(), Rivet::FourVector::eta(), mass(), p(), Rivet::FourVector::phi(), pt(), px(), py(), pz(), setPE(), Rivet::sqr(), and Rivet::FourVector::theta().

Referenced by mkThetaPhiMPt().


The documentation for this class was generated from the following file:
  • /Users/chrisg/software/rivet/include/Rivet/Math/Vector4.hh