class Rivet::FourMomentum
Rivet::FourMomentum
Specialized version of the FourVector with momentum/energy functionality.
#include <Vector4.hh>
Inherits from Rivet::FourVector, Rivet::Vector< 4 >
Public Types
Name | |
---|---|
using Eigen::Matrix< double, N, 1 > | EVector Vector. |
Public Functions
Name | |
---|---|
FourMomentum & | setE(double E) Set energy ( E ) (time component of momentum). |
FourMomentum & | setPx(double px) Set x-component of momentum ( p_x ). |
FourMomentum & | setPy(double py) Set y-component of momentum ( p_y ). |
FourMomentum & | setPz(double pz) Set z-component of momentum ( p_z ). |
FourMomentum & | setPE(double px, double py, double pz, double E) Set the p coordinates and energy simultaneously. |
FourMomentum & | setXYZE(double px, double py, double pz, double E) Alias for setPE. |
FourMomentum & | setPM(double px, double py, double pz, double mass) Set the p coordinates and mass simultaneously. |
FourMomentum & | setXYZM(double px, double py, double pz, double mass) Alias for setPM. |
FourMomentum & | setEtaPhiME(double eta, double phi, double mass, double E) |
FourMomentum & | setEtaPhiMPt(double eta, double phi, double mass, double pt) |
FourMomentum & | setRapPhiME(double y, double phi, double mass, double E) |
FourMomentum & | setRapPhiMPt(double y, double phi, double mass, double pt) |
FourMomentum & | setThetaPhiME(double theta, double phi, double mass, double E) |
FourMomentum & | setThetaPhiMPt(double theta, double phi, double mass, double pt) |
FourMomentum & | setPtPhiME(double pt, double phi, double mass, double E) |
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} ). |
double | gamma() const |
Vector3 | gammaVec() const |
double | beta() const |
Vector3 | betaVec() const |
FourMomentum & | operator*=(double a) Multiply by a scalar. |
FourMomentum & | operator/=(double a) Divide by a scalar. |
FourMomentum & | operator+=(const FourMomentum & v) Add to this 4-vector. NB time as well as space components are added. |
FourMomentum & | operator-=(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. |
FourMomentum | mkXYZE(double px, double py, double pz, double E) Make a vector from (px,py,pz,E) coordinates. |
FourMomentum | mkXYZM(double px, double py, double pz, double mass) Make a vector from (px,py,pz) coordinates and the mass. |
FourMomentum | mkEtaPhiME(double eta, double phi, double mass, double E) Make a vector from (eta,phi,energy) coordinates and the mass. |
FourMomentum | mkEtaPhiMPt(double eta, double phi, double mass, double pt) Make a vector from (eta,phi,pT) coordinates and the mass. |
FourMomentum | mkRapPhiME(double y, double phi, double mass, double E) Make a vector from (y,phi,energy) coordinates and the mass. |
FourMomentum | mkRapPhiMPt(double y, double phi, double mass, double pt) Make a vector from (y,phi,pT) coordinates and the mass. |
FourMomentum | mkThetaPhiME(double theta, double phi, double mass, double E) Make a vector from (theta,phi,energy) coordinates and the mass. |
FourMomentum | mkThetaPhiMPt(double theta, double phi, double mass, double pt) Make a vector from (theta,phi,pT) coordinates and the mass. |
FourMomentum | mkPtPhiME(double pt, double phi, double mass, double E) Make a vector from (pT,phi,energy) coordinates and the mass. |
FourMomentum() | |
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) | |
~FourMomentum() | |
double | t() const |
double | t2() const |
FourVector & | setT(const double t) |
double | x() const |
double | x2() const |
FourVector & | setX(const double x) |
double | y() const |
double | y2() const |
FourVector & | setY(const double y) |
double | z() const |
double | z2() const |
FourVector & | setZ(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 ( |
double | abseta() const Get the ( |
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 (+ - - -). |
FourVector & | operator+=(const FourVector & v) Add to this 4-vector. |
FourVector & | operator-=(const FourVector & v) Subtract from this 4-vector. NB time as well as space components are subtracted. |
const double & | get(const size_t index) const |
double & | get(const size_t index) |
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. |
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 |
Friends
Name | |
---|---|
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) |
Additional inherited members
Public Functions inherited from Rivet::FourVector
Name | |
---|---|
FourVector() | |
template <typename V4TYPE ,typename std::enable_if< HasXYZT< V4TYPE >::value, int >::type DUMMY =0> | FourVector(const V4TYPE & other) |
FourVector(const Vector< 4 > & other) | |
FourVector(const double t, const double x, const double y, const double z) | |
virtual | ~FourVector() |
Public Functions inherited from Rivet::Vector< 4 >
Name | |
---|---|
Vector() | |
Vector(const Vector< N > & other) |
Public Types Documentation
using EVector
using Rivet::Vector< N >::EVector = Eigen::Matrix<double,N,1>;
Vector.
Public Functions Documentation
function setE
inline FourMomentum & setE(
double E
)
Set energy ( E ) (time component of momentum).
function setPx
inline FourMomentum & setPx(
double px
)
Set x-component of momentum ( p_x ).
function setPy
inline FourMomentum & setPy(
double py
)
Set y-component of momentum ( p_y ).
function setPz
inline FourMomentum & setPz(
double pz
)
Set z-component of momentum ( p_z ).
function setPE
inline FourMomentum & setPE(
double px,
double py,
double pz,
double E
)
Set the p coordinates and energy simultaneously.
function setXYZE
inline FourMomentum & setXYZE(
double px,
double py,
double pz,
double E
)
Alias for setPE.
function setPM
inline FourMomentum & setPM(
double px,
double py,
double pz,
double mass
)
Set the p coordinates and mass simultaneously.
function setXYZM
inline FourMomentum & setXYZM(
double px,
double py,
double pz,
double mass
)
Alias for setPM.
function setEtaPhiME
inline FourMomentum & setEtaPhiME(
double eta,
double phi,
double mass,
double E
)
Set the vector state from (eta,phi,energy) coordinates and the mass
eta = -ln(tan(theta/2)) -> theta = 2 atan(exp(-eta))
function setEtaPhiMPt
inline FourMomentum & setEtaPhiMPt(
double eta,
double phi,
double mass,
double pt
)
Set the vector state from (eta,phi,pT) coordinates and the mass
eta = -ln(tan(theta/2)) -> theta = 2 atan(exp(-eta))
function setRapPhiME
inline FourMomentum & setRapPhiME(
double y,
double phi,
double mass,
double E
)
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)
function setRapPhiMPt
inline FourMomentum & setRapPhiMPt(
double y,
double phi,
double mass,
double pt
)
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]
function setThetaPhiME
inline FourMomentum & setThetaPhiME(
double theta,
double phi,
double mass,
double E
)
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)
function setThetaPhiMPt
inline FourMomentum & setThetaPhiMPt(
double theta,
double phi,
double mass,
double pt
)
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)
function setPtPhiME
inline FourMomentum & setPtPhiME(
double pt,
double phi,
double mass,
double E
)
Set the vector state from (pT,phi,energy) coordinates and the mass
pz = sqrt(E^2 - mass^2 - pt^2)
function E
inline double E() const
Get energy ( E ) (time component of momentum).
function E2
inline double E2() const
Get energy-squared ( E^2 ).
function px
inline double px() const
Get x-component of momentum ( p_x ).
function px2
inline double px2() const
Get x-squared ( p_x^2 ).
function py
inline double py() const
Get y-component of momentum ( p_y ).
function py2
inline double py2() const
Get y-squared ( p_y^2 ).
function pz
inline double pz() const
Get z-component of momentum ( p_z ).
function pz2
inline double pz2() const
Get z-squared ( p_z^2 ).
function mass
inline double mass() const
Get the mass ( m = \sqrt{E^2 - p^2} ) (the Lorentz self-invariant).
For spacelike momenta, the mass will be -sqrt(|mass2|).
function mass2
inline double mass2() const
Get the squared mass ( m^2 = E^2 - p^2 ) (the Lorentz self-invariant).
function p3
inline Vector3 p3() const
Get 3-momentum part, ( p ).
function p
inline double p() const
Get the modulus of the 3-momentum.
function p2
inline double p2() const
Get the modulus-squared of the 3-momentum.
function rapidity
inline double rapidity() const
Calculate the rapidity.
function rap
inline double rap() const
Alias for rapidity.
function absrapidity
inline double absrapidity() const
Absolute rapidity.
function absrap
inline double absrap() const
Absolute rapidity.
function pTvec
inline Vector3 pTvec() const
Calculate the transverse momentum vector ( \vec{p}_T ).
function ptvec
inline Vector3 ptvec() const
Synonym for pTvec.
function pT2
inline double pT2() const
Calculate the squared transverse momentum ( p_T^2 ).
function pt2
inline double pt2() const
Calculate the squared transverse momentum ( p_T^2 ).
function pT
inline double pT() const
Calculate the transverse momentum ( p_T ).
function pt
inline double pt() const
Calculate the transverse momentum ( p_T ).
function Et2
inline double Et2() const
Calculate the transverse energy ( E_T^2 = E^2 \sin^2{\theta} ).
function Et
inline double Et() const
Calculate the transverse energy ( E_T = E \sin{\theta} ).
function gamma
inline double gamma() const
Note: ( \gamma = E/mc^2 ) so we rely on the c=1 convention
Calculate the boost factor ( \gamma ).
function gammaVec
inline Vector3 gammaVec() const
Note: ( \gamma = E/mc^2 ) so we rely on the c=1 convention
Calculate the boost vector ( \vec{\gamma} ).
function beta
inline double beta() const
Note: ( \beta = pc/E ) so we rely on the c=1 convention
Calculate the boost factor ( \beta ).
function betaVec
inline Vector3 betaVec() const
Note: ( \beta = pc/E ) so we rely on the c=1 convention
Calculate the boost vector ( \vec{\beta} ).
function operator*=
inline FourMomentum & operator*=(
double a
)
Multiply by a scalar.
function operator/=
inline FourMomentum & operator/=(
double a
)
Divide by a scalar.
function operator+=
inline FourMomentum & operator+=(
const FourMomentum & v
)
Add to this 4-vector. NB time as well as space components are added.
function operator-=
inline FourMomentum & operator-=(
const FourMomentum & v
)
Subtract from this 4-vector. NB time as well as space components are subtracted.
function operator-
inline FourMomentum operator-() const
Multiply all components (time and space) by -1.
function reverse
inline FourMomentum reverse() const
Multiply space components only by -1.
function mkXYZE
static inline FourMomentum mkXYZE(
double px,
double py,
double pz,
double E
)
Make a vector from (px,py,pz,E) coordinates.
function mkXYZM
static inline FourMomentum mkXYZM(
double px,
double py,
double pz,
double mass
)
Make a vector from (px,py,pz) coordinates and the mass.
function mkEtaPhiME
static inline FourMomentum mkEtaPhiME(
double eta,
double phi,
double mass,
double E
)
Make a vector from (eta,phi,energy) coordinates and the mass.
function mkEtaPhiMPt
static inline FourMomentum mkEtaPhiMPt(
double eta,
double phi,
double mass,
double pt
)
Make a vector from (eta,phi,pT) coordinates and the mass.
function mkRapPhiME
static inline FourMomentum mkRapPhiME(
double y,
double phi,
double mass,
double E
)
Make a vector from (y,phi,energy) coordinates and the mass.
function mkRapPhiMPt
static inline FourMomentum mkRapPhiMPt(
double y,
double phi,
double mass,
double pt
)
Make a vector from (y,phi,pT) coordinates and the mass.
function mkThetaPhiME
static inline FourMomentum mkThetaPhiME(
double theta,
double phi,
double mass,
double E
)
Make a vector from (theta,phi,energy) coordinates and the mass.
function mkThetaPhiMPt
static inline FourMomentum mkThetaPhiMPt(
double theta,
double phi,
double mass,
double pt
)
Make a vector from (theta,phi,pT) coordinates and the mass.
function mkPtPhiME
static inline FourMomentum mkPtPhiME(
double pt,
double phi,
double mass,
double E
)
Make a vector from (pT,phi,energy) coordinates and the mass.
function FourMomentum
inline FourMomentum()
function FourMomentum
template <typename V4TYPE ,
typename std::enable_if< HasXYZT< V4TYPE >::value, int >::type DUMMY =0>
inline FourMomentum(
const V4TYPE & other
)
function FourMomentum
inline FourMomentum(
const Vector< 4 > & other
)
function FourMomentum
inline FourMomentum(
const double E,
const double px,
const double py,
const double pz
)
function ~FourMomentum
inline ~FourMomentum()
function t
inline double t() const
function t2
inline double t2() const
function setT
inline FourVector & setT(
const double t
)
function x
inline double x() const
function x2
inline double x2() const
function setX
inline FourVector & setX(
const double x
)
function y
inline double y() const
function y2
inline double y2() const
function setY
inline FourVector & setY(
const double y
)
function z
inline double z() const
function z2
inline double z2() const
function setZ
inline FourVector & setZ(
const double z
)
function invariant
inline double invariant() const
function isNull
inline bool isNull() const
function angle
inline double angle(
const FourVector & v
) const
Angle between this vector and another.
function angle
inline double angle(
const Vector3 & v3
) const
Angle between this vector and another (3-vector)
function polarRadius2
inline 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.
function perp2
inline double perp2() const
Synonym for polarRadius2.
function rho2
inline double rho2() const
Synonym for polarRadius2.
function polarRadius
inline double polarRadius() const
Magnitude of projection of 3-vector on to the ( x-y ) plane.
function perp
inline double perp() const
Synonym for polarRadius.
function rho
inline double rho() const
Synonym for polarRadius.
function polarVec
inline Vector3 polarVec() const
Projection of 3-vector on to the ( x-y ) plane.
function perpVec
inline Vector3 perpVec() const
Synonym for polarVec.
function rhoVec
inline Vector3 rhoVec() const
Synonym for polarVec.
function azimuthalAngle
inline double azimuthalAngle(
const PhiMapping mapping =ZERO_2PI
) const
Angle subtended by the 3-vector’s projection in x-y and the x-axis.
function phi
inline double phi(
const PhiMapping mapping =ZERO_2PI
) const
Synonym for azimuthalAngle.
function polarAngle
inline double polarAngle() const
Angle subtended by the 3-vector and the z-axis.
function theta
inline double theta() const
Synonym for polarAngle.
function pseudorapidity
inline double pseudorapidity() const
Pseudorapidity (defined purely by the 3-vector components)
function eta
inline double eta() const
Synonym for pseudorapidity.
function abspseudorapidity
inline double abspseudorapidity() const
Get the ( |\eta| ) directly.
function abseta
inline double abseta() const
Get the ( |\eta| ) directly (alias).
function vector3
inline Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
function operator Vector3
inline operator Vector3() const
Implicit cast to a 3-vector.
function contract
inline double contract(
const FourVector & v
) const
Contract two 4-vectors, with metric signature (+ - - -).
function dot
inline double dot(
const FourVector & v
) const
Contract two 4-vectors, with metric signature (+ - - -).
function operator*
inline double operator*(
const FourVector & v
) const
Contract two 4-vectors, with metric signature (+ - - -).
function operator+=
inline FourVector & operator+=(
const FourVector & v
)
Add to this 4-vector.
function operator-=
inline FourVector & operator-=(
const FourVector & v
)
Subtract from this 4-vector. NB time as well as space components are subtracted.
function get
inline const double & get(
const size_t index
) const
function get
inline double & get(
const size_t index
)
function operator[]
inline const double & operator[](
const size_t index
) const
Direct access to vector elements by index.
function operator[]
inline double & operator[](
const size_t index
)
Direct access to vector elements by index.
function set
inline Vector< N > & set(
const size_t index,
const double value
)
Set indexed value.
function size
inline constexpr size_t size() const
Vector dimensionality.
function isZero
inline bool isZero(
double tolerance =1E-5
) const
Check for nullness, allowing for numerical precision.
function mod2
inline double mod2() const
Calculate the modulus-squared of a vector. ( \sum_{i=1}^N x_i^2 ).
function mod
inline double mod() const
Calculate the modulus of a vector. ( \sqrt{\sum_{i=1}^N x_i^2} ).
function operator==
inline bool operator==(
const Vector< N > & a
) const
function operator!=
inline bool operator!=(
const Vector< N > & a
) const
Friends
friend multiply
friend FourMomentum multiply(
const double a,
const FourMomentum & v
);
friend multiply
friend FourMomentum multiply(
const FourMomentum & v,
const double a
);
friend add
friend FourMomentum add(
const FourMomentum & a,
const FourMomentum & b
);
friend transform
friend FourMomentum transform(
const LorentzTransform & lt,
const FourMomentum & v4
);
Updated on 2022-08-07 at 20:17:17 +0100