Lorentz Vector
Interface registry
QEDbase.register_LorentzVectorLike
— Functionregister_LorentzVectorLike(T)
Function to register a custom type as a LorentzVectorLike.
Ensure the passed custom type has implemented at least the function getT, getX, getY, getZ
and enables getter functions of the lorentz vector library for the given type. If additionally the functions setT!, setX!, setY!, setZ!
are implemened for the passed custom type, also the setter functions of the Lorentz vector interface are enabled.
Built-in Lorentz vector types
QEDbase.AbstractLorentzVector
— Typeabstract type AbstractLorentzVector{T} <: StaticArraysCore.FieldVector{4, T}
Abstract type to model generic Lorentz vectors, i.e. elements of a minkowski-like space, where the component space is arbitray.
QEDbase.AbstractFourMomentum
— TypeAbstractFourMomentum
Abstract base type for four-momentas, representing one energy and three spacial components.
Also see: QEDcore.SFourMomentum
, QEDcore.MFourMomentum
Accessor functions
QEDbase.minkowski_dot
— Functionminkowski_dot(v1,v2)
Return the Minkowski dot product of two LorentzVectorLike
.
If (t1,x1,y1,z1)
and (t2,x2,y2,z2)
are two LorentzVectorLike
, this is equivalent to
t1*t2 - (x1*x2 + y1*y2 + z1*z2)
We use the mostly minus metric.
QEDbase.mdot
— FunctionFunction alias for minkowski_dot
.
QEDbase.getT
— FunctiongetT(lv)
Return the 0-component of a given LorentzVectorLike
.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to t
.
QEDbase.getX
— FunctiongetX(lv)
Return the 1-component of a given LorentzVectorLike
.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to x
.
QEDbase.getY
— FunctiongetY(lv)
Return the 2-component of a given LorentzVectorLike
.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to y
.
QEDbase.getZ
— FunctiongetZ(lv)
Return the 3-component of a given LorentzVectorLike
.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to z
.
QEDbase.getMagnitude2
— FunctiongetMagnitude2(lv)
Return the square of the magnitude of a given LorentzVectorLike
, i.e. the sum of the squared spatial components.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to x^2+ y^2 + z^2
.
This function differs from a similar function for the TLorentzVector
used in the famous ROOT
library.
QEDbase.getMag2
— FunctionFunctiom alias for getMagnitude2
.
QEDbase.getMagnitude
— FunctiongetMagnitude(lv)
Return the magnitude of a given LorentzVectorLike
, i.e. the euklidian norm spatial components.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to sqrt(x^2 + y^2 + z^2)
.
This function differs from a similar function for the TLorentzVector
used in the famous ROOT
library.
QEDbase.getMag
— FunctionFunctiom alias for getMagnitude
.
QEDbase.getInvariantMass2
— FunctiongetInvariantMass2(lv)
Return the squared invariant mass of a given LorentzVectorLike
, i.e. the minkowski dot with itself.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to t^2 - (x^2 + y^2 + z^2)
.
QEDbase.getMass2
— FunctionFunction alias for getInvariantMass2
QEDbase.getInvariantMass
— FunctiongetInvariantMass(lv)
Return the the invariant mass of a given LorentzVectorLike
, i.e. the square root of the minkowski dot with itself.
If (t,x,y,z)
is a LorentzVectorLike
, this is equivalent to sqrt(t^2 - (x^2 + y^2 + z^2))
.
If the squared invariant mass m2
is negative, -sqrt(-m2)
is returned.
QEDbase.getMass
— FunctionFunction alias for getInvariantMass
.
QEDbase.getE
— FunctiongetE(lv)
Return the energy component of a given LorentzVectorLike
, i.e. its 0-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to E
.
QEDbase.getEnergy
— FunctionFunction alias for getE
.
QEDbase.getPx
— FunctiongetPx(lv)
Return the $p_x$ component of a given LorentzVectorLike
, i.e. its 1-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to px
.
QEDbase.getPy
— FunctiongetPy(lv)
Return the $p_y$ component of a given LorentzVectorLike
, i.e. its 2-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to py
.
QEDbase.getPz
— FunctiongetPz(lv)
Return the $p_z$ component of a given LorentzVectorLike
, i.e. its 3-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to pz
.
QEDbase.getBeta
— FunctiongetBeta(lv)
Return magnitude of the beta vector for a given LorentzVectorLike
, i.e. the magnitude of the LorentzVectorLike
divided by its 0-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to sqrt(px^2 + py^2 + pz^2)/E
.
QEDbase.getGamma
— FunctiongetGamma(lv)
Return the relativistic gamma factor for a given LorentzVectorLike
, i.e. the inverse square root of one minus the beta vector squared.
If (E,px,py,pz)
is a LorentzVectorLike
with beta vector β
, this is equivalent to 1/sqrt(1- β^2)
.
QEDbase.getTransverseMomentum2
— FunctiongetTransverseMomentum2(lv)
Return the squared transverse momentum for a given LorentzVectorLike
, i.e. the sum of its squared 1- and 2-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to px^2 + py^2
.
The transverse components are defined w.r.t. to the 3-axis.
QEDbase.getPt2
— FunctionFunction alias for getTransverseMomentum2
.
QEDbase.getPerp2
— FunctionFunction alias for getTransverseMomentum2
.
QEDbase.getTransverseMomentum
— FunctiongetTransverseMomentum(lv)
Return the transverse momentum for a given LorentzVectorLike
, i.e. the magnitude of its transverse components.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to sqrt(px^2 + py^2)
.
The transverse components are defined w.r.t. to the 3-axis.
QEDbase.getPt
— FunctionFunction alias for getTransverseMomentum
.
QEDbase.getPerp
— FunctionFunction alias for getTransverseMomentum
.
QEDbase.getTransverseMass2
— FunctiongetTransverseMass2(lv)
Return the squared transverse mass for a given LorentzVectorLike
, i.e. the difference of its squared 0- and 3-component.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to E^2 - pz^2
.
The transverse components are defined w.r.t. to the 3-axis.
QEDbase.getMt2
— FunctionFunction alias for getTransverseMass2
QEDbase.getTransverseMass
— FunctiongetTransverseMass(lv)
Return the transverse momentum for a given LorentzVectorLike
, i.e. the square root of its squared transverse mass.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to sqrt(E^2 - pz^2)
.
The transverse components are defined w.r.t. to the 3-axis.
If the squared transverse mass mT2
is negative, -sqrt(-mT2)
is returned.
QEDbase.getMt
— FunctionFunction alias for getTransverseMass
QEDbase.getRapidity
— FunctiongetRapidity(lv)
Return the rapidity for a given LorentzVectorLike
.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to 0.5*log((E+pz)/(E-pz))
.
The transverse components are defined w.r.t. to the 3-axis.
QEDbase.getRho2
— FunctionFunction alias for getMagnitude2
QEDbase.getRho
— FunctionFunction alias for getMagnitude
QEDbase.getTheta
— FunctiongetTheta(lv)
Return the theta angle for a given LorentzVectorLike
, i.e. the polar angle of its spatial components in spherical coordinates.
If (E,px,py,pz)
is a LorentzVectorLike
with magnitude rho
, this is equivalent to arccos(pz/rho)
, which is also equivalent to arctan(sqrt(px^2+py^2)/pz)
.
The spherical coordinates are defined w.r.t. to the 3-axis.
QEDbase.getCosTheta
— FunctiongetCosTheta(lv)
Return the cosine of the theta angle for a given LorentzVectorLike
.
This is an equivalent but faster version of cos(getTheta(lv))
; see getTheta
.
QEDbase.getPhi
— FunctiongetPhi(lv)
Return the phi angle for a given LorentzVectorLike
, i.e. the azimuthal angle of its spatial components in spherical coordinates.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to atan(py,px)
.
The spherical coordinates are defined w.r.t. to the 3-axis.
QEDbase.getCosPhi
— FunctiongetCosPhi(lv)
Return the cosine of the phi angle for a given LorentzVectorLike
.
This is an equivalent but faster version of cos(getPhi(lv))
; see getPhi
.
QEDbase.getSinPhi
— FunctiongetSinPhi(lv)
Return the sine of the phi angle for a given LorentzVectorLike
.
This is an equivalent but faster version of sin(getPhi(lv))
; see getPhi
.
QEDbase.getPlus
— FunctiongetPlus(lv)
Return the plus component for a given LorentzVectorLike
in light-cone coordinates.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to (E+pz)/2
.
The light-cone coordinates are defined w.r.t. to the 3-axis.
The definition `p^+ := (E + p_z)/2
differs from the usual definition of light-cone coordinates in general relativity.
QEDbase.getMinus
— FunctiongetMinus(lv)
Return the minus component for a given LorentzVectorLike
in light-cone coordinates.
If (E,px,py,pz)
is a LorentzVectorLike
, this is equivalent to (E-pz)/2
.
The light-cone coordinates are defined w.r.t. to the 3-axis.
The definition `p^- := (E - p_z)/2
differs from the usual definition of light-cone coordinates in general relativity.
Setter functions
QEDbase.setE!
— FunctionsetE!(lv,value)
Sets the energy component of a given LorentzVectorLike
to a given value
.
The value
set with setE!
is then returned by getE
.
QEDbase.setEnergy!
— FunctionFunction alias for setE!
.
QEDbase.setPx!
— FunctionsetPx!(lv,value)
Sets the 1-component of a given LorentzVectorLike
to a given value
.
The value
set with setPx!
is then returned by getPx
.
QEDbase.setPy!
— FunctionsetPy!(lv,value)
Sets the 2-component of a given LorentzVectorLike
to a given value
.
The value
set with setPy!
is then returned by getPy
.
QEDbase.setPz!
— FunctionsetPz!(lv,value)
Sets the 3-component of a given LorentzVectorLike
to a given value
.
The value
set with setPz!
is then returned by getPz
.
QEDbase.setTheta!
— FunctionsetTheta!(lv,value)
Sets the theta angle of a LorentzVectorLike
to a given value
.
The value
set with setTheta!
is then returned by getTheta
. Since the theta angle is computed on the call of getTheta
, the setter setTheta!
changes several components of the given LorentzVectorLike
.
QEDbase.setCosTheta!
— FunctionsetCosTheta!(lv,value)
Sets the cosine of the theta angle of a LorentzVectorLike
to a given value
.
The value
set with setCosTheta!
is then returned by getCosTheta
. Since the cosine of the theta angle is computed on the call of getCosTheta
, the setter setCosTheta!
changes several components of the given LorentzVectorLike
.
QEDbase.setRho!
— FunctionsetRho!(lv,value)
Sets the magnitude of a LorentzVectorLike
to a given value
.
The value
set with setRho!
is then returned by getRho
. Since the magnitude is computed on the call of getRho
, the setter setRho!
changes several components of the given LorentzVectorLike
.
QEDbase.setPhi!
— FunctionsetPhi!(lv,value)
Sets the phi angle of a LorentzVectorLike
to a given value
.
The value
set with setPhi!
is then returned by getPhi
. Since the phi angle is computed on the call of getPhi
, the setter setPhi!
changes several components of the given LorentzVectorLike
.
QEDbase.setPlus!
— FunctionsetPlus!(lv,value)
Sets the plus component of a LorentzVectorLike
to a given value
.
The value
set with setPlus!
is then returned by getPlus
. Since the plus component is computed on the call of getPlus
, the setter setPlus!
changes several components of the given LorentzVectorLike
.
QEDbase.setMinus!
— FunctionsetMinus!(lv,value)
Sets the minus component of a LorentzVectorLike
to a given value
.
The value
set with setMinus!
is then returned by getMinus
. Since the minus component is computed on the call of getMinus
, the setter setMinus!
changes several components of the given LorentzVectorLike
.
QEDbase.setTransverseMomentum!
— FunctionsetTransverseMomentum!(lv,value)
Sets the transverse momentum of a LorentzVectorLike
to a given value
.
The value
set with setTransverseMomentum!
is then returned by getTransverseMomentum
. Since the transverse momentum is computed on the call of getTransverseMomentum
, the setter setTransverseMomentum!
changes several components of the given LorentzVectorLike
.
QEDbase.setPerp!
— FunctionFunction alias for setTransverseMomentum!
.
QEDbase.setPt!
— FunctionFunction alias for setTransverseMomentum!
.
QEDbase.setTransverseMass!
— FunctionsetTransverseMass!(lv,value)
Sets the transverse mass of a LorentzVectorLike
to a given value
.
The value
set with setTransverseMass!
is then returned by getTransverseMass
. Since the transverse mass is computed on the call of getTransverseMass
, the setter setTransverseMass!
changes several components of the given LorentzVectorLike
.
QEDbase.setMt!
— FunctionFunction alias for setTransverseMass!
.
QEDbase.setRapidity!
— FunctionsetRapidity!(lv,value)
Sets the rapidity of a LorentzVectorLike
to a given value
.
The value
set with setRapidity!
is then returned by getRapidity
. Since the rapidity is computed on the call of setRapidity
, the setter setRapidity!
changes several components of the given LorentzVectorLike
.
Utilities
QEDbase.isonshell
— Functionisonshell(mom, mass)
On-shell check of a given four-momentum mom
w.r.t. a given mass mass
.
For AbstactFourMomentum
, the element type is fixed to Float64
, limiting the precision of comparisons between elements. The current implementation has been tested within the boundaries for energy scales E
with 1e-9 <= E <= 1e5
. In those bounds, the mass error, which is correctly detected as off-shell, is 1e-4
times the mean value of the components, but has at most the value 0.01
, e.g. at the high energy end. The energy scales correspond to 0.5 meV
for the lower bound and 50 GeV
for the upper bound.
- if
AbstractFourMomentum
is updated to elementtypesT<:Real
, theAbstractLorentzVector
should be updated with theAbstractFourMomentum
.
QEDbase.assert_onshell
— Functionassert_onshell(mom, mass)
Assertion if a FourMomentum mom
is on-shell w.r.t a given mass mass
.
The precision of this functions is explained in isonshell
.
QEDbase.AbstractCoordinateTransformation
— TypeAbstractCoordinateTransformation
Abstract base type for coordinate transformations supposed to be acting on four-momenta. Every subtype of trafo::AbstractCoordinateTransformation
should implement the following interface functions:
QEDbase._transform(trafo,p)
: transformsp
Base.inv(trafo)
: returns the inverted transform
Example
Implementing the interface by defining the interface functions:
julia> using QEDbase
julia> using QEDcore
julia> struct TestTrafo{T} <: AbstractCoordinateTransformation
a::T
end
julia> QEDbase._transform(trafo::TestTrafo,p) = trafo.a*p
julia> Base.inv(trafo::TestTrafo) = TestTrafo(inv(trafo.a))
The TestTrafo
can then be used to transform four-momenta:
julia> trafo = TestTrafo(2.0)
TestTrafo{Float64}(2.0)
julia> p = SFourMomentum(4,3,2,1)
4-element SFourMomentum with indices SOneTo(4):
4.0
3.0
2.0
1.0
julia> trafo(p) # multiply every component with 2.0
4-element SFourMomentum with indices SOneTo(4):
8.0
6.0
4.0
2.0
julia> inv(trafo)(p) # divide every component by 2.0
4-element SFourMomentum with indices SOneTo(4):
2.0
1.5
1.0
0.5