Phase Space Points
Types and Aliases
QEDcore.ParticleStateful
— TypeParticleStateful <: AbstractParticle
Representation of a particle with a state. It has four fields:
dir::ParticleDirection
: The direction of the particle,Incoming()
orOutgoing()
.species::AbstractParticleType
: The species of the particle,Electron()
,Positron()
etc.mom::AbstractFourMomentum
: The momentum of the particle.
Overloads for is_fermion
, is_boson
, is_particle
, is_anti_particle
, is_incoming
, is_outgoing
, mass
, and charge
are provided, delegating the call to the correct field and thus implementing the AbstractParticle
interface.
julia> using QEDcore
julia> ParticleStateful(Incoming(), Electron(), SFourMomentum(1, 0, 0, 0))
ParticleStateful: incoming electron
momentum: [1.0, 0.0, 0.0, 0.0]
julia> ParticleStateful(Outgoing(), Photon(), SFourMomentum(1, 0, 0, 0))
ParticleStateful: outgoing photon
momentum: [1.0, 0.0, 0.0, 0.0]
QEDcore.PhaseSpacePoint
— TypePhaseSpacePoint
Representation of a point in the phase space of a process. Contains the process (AbstractProcessDefinition
), the model (AbstractModelDefinition
), the phase space definition (AbstractPhasespaceDefinition
), and stateful incoming and outgoing particles (ParticleStateful
).
The legality of the combination of the given process and the incoming and outgoing particles is checked on construction. If the numbers of particles mismatch, the types of particles mismatch (note that order is important), or incoming particles have an Outgoing
direction, an error is thrown.
julia> using QEDcore; using QEDprocesses
julia> PhaseSpacePoint(
Compton(),
PerturbativeQED(),
PhasespaceDefinition(SphericalCoordinateSystem(), ElectronRestFrame()),
(
ParticleStateful(Incoming(), Electron(), SFourMomentum(1, 0, 0, 0)),
ParticleStateful(Incoming(), Photon(), SFourMomentum(1, 0, 0, 0))
),
(
ParticleStateful(Outgoing(), Electron(), SFourMomentum(1, 0, 0, 0)),
ParticleStateful(Outgoing(), Photon(), SFourMomentum(1, 0, 0, 0))
)
)
PhaseSpacePoint:
process: one-photon Compton scattering
model: perturbative QED
phasespace definition: spherical coordinates in electron rest frame
incoming particles:
-> incoming electron: [1.0, 0.0, 0.0, 0.0]
-> incoming photon: [1.0, 0.0, 0.0, 0.0]
outgoing particles:
-> outgoing electron: [1.0, 0.0, 0.0, 0.0]
-> outgoing photon: [1.0, 0.0, 0.0, 0.0]
PhaseSpacePoint
s can be constructed with only one of their in- or out-channel set. For this, see the special constructors InPhaseSpacePoint
and OutPhaseSpacePoint
. The InPhaseSpacePoint
and OutPhaseSpacePoint
type definitions can be used to dispatch on such PhaseSpacePoint
s. Note that a full PhaseSpacePoint
containing both its in- and out-channel matches both, .i.e. psp isa InPhaseSpacePoint
and psp isa OutPhaseSpacePoint
both evaluate to true if psp contains both channels. A completely empty PhaseSpacePoint
is not allowed.
QEDcore.InPhaseSpacePoint
— TypeInPhaseSpacePoint
A partial type specialization on PhaseSpacePoint
which can be used for dispatch in functions requiring only the in channel of the phase space to exist, for example implementations of _incident_flux
. No restrictions are imposed on the out-channel, which may or may not exist.
See also: OutPhaseSpacePoint
QEDcore.OutPhaseSpacePoint
— TypeOutPhaseSpacePoint
A partial type specialization on PhaseSpacePoint
which can be used for dispatch in functions requiring only the out channel of the phase space to exist. No restrictions are imposed on the in-channel, which may or may not exist.
See also: InPhaseSpacePoint
Accessors
Base.getindex
— FunctionBase.getindex(psp::PhaseSpacePoint, dir::Incoming, n::Int)
Overload for the array indexing operator []
. Returns the nth incoming particle in this phase space point.
Base.getindex(psp::PhaseSpacePoint, dir::Outgoing, n::Int)
Overload for the array indexing operator []
. Returns the nth outgoing particle in this phase space point.
QEDcore._momentum_type
— Function_momentum_type(psp::PhaseSpacePoint)
_momentum_type(type::Type{PhaseSpacePoint})
Returns the element type of the PhaseSpacePoint
object or type, e.g. SFourMomentum
.
julia> using QEDcore; using QEDprocesses
julia> psp = PhaseSpacePoint(Compton(), PerturbativeQED(), PhasespaceDefinition(SphericalCoordinateSystem(), ElectronRestFrame()), Tuple(rand(SFourMomentum) for _ in 1:2), Tuple(rand(SFourMomentum) for _ in 1:2));
julia> QEDcore._momentum_type(psp)
SFourMomentum
julia> QEDcore._momentum_type(typeof(psp))
SFourMomentum