QEDprocesses
Documentation for QEDprocesses.
QEDprocesses.Compton
QEDprocesses.ScatteringProcess
QEDbase._averaging_norm
QEDbase.in_phase_space_dimension
QEDprocesses._base_component_type
QEDprocesses.isphysical
QEDprocesses.Compton
— TypeCompton(
in_spin [= AllSpin()]
in_pol [= AllPol()]
out_spin [= AllSpin()]
out_pol [= AllPol()]
)
QEDprocesses.ScatteringProcess
— TypeScatteringProcess <: AbstractProcessDefinition
Generic implementation for scattering processes of arbitrary particles. Currently, only calculations in combination with PerturbativeQED
are supported. However, this is supposed to describe scattering processes with any number of incoming and outgoing particles, and any combination of spins or polarizations for the particles.
The isphysical
function can be used to check whether the process is possible in perturbative QED.
The computation of cross sections and probabilities is currently unimplemented.
Constructors
ScatteringProcess(
in_particles::Tuple{AbstractParticleType},
out_particles::Tuple{AbstractParticleType},
[in_sp::Tuple{AbstractSpinOrPolarization},
out_sp::Tuple{AbstractSpinOrPolarization}]
)
Constructor for a ScatteringProcess with the given incoming and outgoing particles and their respective spins and pols. The constructor asserts that the particles are compatible with their respective spins and polarizations. If the assertion fails, an InvalidInputError
is thrown.
The in_sp
and out_sp
parameters can be omitted in which case all spins and polarizations will be set to AllSpin
and AllPol
for every fermion and boson, respectively.
QEDbase._averaging_norm
— Method_averaging_norm(proc::Compton)
We average over the initial spins and pols, and sum over final.
QEDbase.in_phase_space_dimension
— Methodin_phase_space_dimension(proc::AbstractProcessDefinition, ::PerturbativeQED)
Return the number of degrees of freedom to determine the incoming phase space for processes in PerturbativeQED.
The current implementation only supports the case where two of the incoming particles collide head-on.
QEDprocesses._base_component_type
— Method_base_component_type(array_of_lv::AbstractArray{LV}) where {LV<:AbstractLorentzVector}
Return the type of the components of given Lorentz vectors, which are by themself elements of an AbstractArray
.
Examples
julia> using QEDbase
julia> using QEDprocesses
julia> v = Vector{SFourMomentum}(undef,10)
julia> QEDprocesses._base_component_type(v)
Float64
QEDprocesses.isphysical
— Methodisphysical(proc::AbstractProcessDefinition, model::PerturbativeQED)
A utility function that returns whether a given AbstractProcessDefinition
conserves the number and charge of fermions and has at least 2 participating particles.