Internal Functions
BoltzmannCollisionIntegral.DoesConserve
— MethodDoesConserve(SMatrix3,SMatrix4,TMatrix1,TMatrix2,Parameters)
Function prints the ratio of the sum of the S and T matrices and their differences, for all interaction paths, as to check number and energy conservation for a particular interaction. Arguments are as outputted by the fload_All
function.
BoltzmannCollisionIntegral.InvarientFlux
— MethodInvarientFlux(s,mu12,mu22)
returns the value of the invarient flux with 's' mandelstram variable and masses 'mass1' and 'mass2'
BoltzmannCollisionIntegral.InvarientFlux2
— MethodInvarientFlux2(s,mass12,mass22)
returns the value of the squared invarient flux with 's' mandelstram variable and masses 'mass1' and 'mass2'
BoltzmannCollisionIntegral.InvarientFlux2Small
— MethodInvarientFluxSmall(sSmol,mass12,mass22)
returns the value of the squared invarient flux with smalled 's' mandelstram variable (sSmol = s - (m1+m2)^2)
BoltzmannCollisionIntegral.InvarientFluxSmall
— MethodInvarientFluxSmall(sSmol,mu12,mu22)
returns the value of the invarient flux with smalled 's' mandelstram variable (sSmol = s - (m1+m2)^2)
BoltzmannCollisionIntegral.Momentum3Value!
— MethodMomentum3Value!(p3v,p3pv,p1v,p2v)
Takes set of random initial particle states 'p1v' and 'p2v' and random output states angles 'p3v[2:3]' and modifies outputs 'p3v' and 'p3pv' values with calculated output momentum and corrects angles if momentum is negative. Function also returns a two bools 'p3physical' and 'p3pphysical' indicating if p3 and p3p are physical momentum states given the inputs. Function also returns a Int 'NumStates' indicating the number of valid output states found.
Requrires normalised masses (mu1,mu2,mu3,mu4) to be defined in advance in Init.jl as const.
Examples
julia> mu1 = 1836.1528e0
julia> mu2 = 1836.1528e0
julia> mu3 = 1836.1528e0
julia> mu4 = 1836.1528e0
julia> p1v = [1e0, 0.5e0, 1.8e0]
julia> p2v = [2e0, 0.2e0, 0.7e0]
julia> p3v = [0e0, 0.3e0, 0.7e0]
julia> p3pv = zeros(Float64,3)
julia> p3pv .= p3v
julia> Momentum3Value!(p3v,p3pv,p1v,p2v,mu1,mu2,mu3,mu4)
(true,true,2)
julia> p3v
3-element Vector{Float64}:
2.04505
0.3
0.7
julia> p3pv
3-element Vector{Float64}
0.691423
-0.3
1.7
BoltzmannCollisionIntegral.PhaseSpaceFactors1!
— MethodPhaseSpaceFactors1!(SMatrix3,SMatrix4,TMatrix,u3_bounds,u4_bounds,p1_bounds,u1_bounds,p2_bounds,u2_bounds,Indistinguishable_12)
Applies phase space volume element factors for 'SMatrix' and 'TMatrix' terms in order to correctly apply 'STSymmetry' corrections.
BoltzmannCollisionIntegral.PhaseSpaceFactors2!
— MethodPhaseSpaceFactors2!(SMatrix3,SMatrix4,TMatrix,p3_bounds,u3_bounds,p4_bounds,u4_bounds,p1_bounds,u1_bounds,p2_bounds,u2_bounds)
To follow 'PhaseSpaceFactors1' and 'STSymmetry'. Corrects phase space factors on 'SMatrix' and 'TMatrix' for use in kinetic codes. Assumes f(x,p,μ)= constant
BoltzmannCollisionIntegral.PhaseSpaceFactorsSync1!
— MethodPhaseSpaceFactorsSync1!(SMatrix,p1val,t1val,p2val,t2val)
Applies phase space volume element factors for 'SMatrix' terms in order to correctly apply 'SyncSymmetry' corrections.
BoltzmannCollisionIntegral.PhaseSpaceFactorsSync2!
— MethodPhaseSpaceFactorsSync2!(SMatrix,p1val,t1val)
To follow 'PhaseSpaceFactorsSync1' and 'SyncSymmetry'. Correct phase space factors on 'SMatrix' for use in kinetic codes.
BoltzmannCollisionIntegral.RPointLogMomentum!
— MethodRPointLogMomentum!(pu,pl,pv,num)
Edits the first element of pv
with a random real-space momentum value between $10^{pl}$ and $10^{pu}$. This sample is chosen by first randomly picking a momentum bin in the range 1:num
and then uniformly sampling a momentum point in real-space (rather than log10 space) between pl and pu which are the momentum values at start and end of that bin. Sampling is done such there will be a constant number of points per momentum-space volume. As the momentum space between $10^{pl}$ and $10^{pu}$ it is a spherical shell hence the correct sampling is $p = (U*(10^{pu})^3+(1-U)*(10^{pl})^3)^{1/3}$ with uniform $U ∈ [0~~1]$.
Assumes $f(x,p,μ)=f(x,\vec{p})*(2πp^2)=const$ in bin, therefore momentum space volume element is $\mathrm{d}p$ and as such uniform sampling corresponds to $U*10^{u}+(1-U)*10^{l}$ where $U$ is a uniform random number between 0 and 1.
If instead $f(x,\vec{p})=const$ in bin, momentum space volume element is $p^2 \mathrm{d}p$ and uniform sampling corresponds to $(10^pu)*\sqrt[3]{U+(1-U)*10^{3pl-3pu}}$ where $U$ is a uniform random number between 0 and 1.
BoltzmannCollisionIntegral.RPointSphereCosTheta!
— MethodRPointSphereTheta!()
Assigns the second (cos(theta)) element of 'a' with a randomly, uniformly sampled values of spherical angles cos(theta).
BoltzmannCollisionIntegral.RPointSphereCosThetaPhi!
— MethodRPointSphereThetaPhi!()
Assigns the second (cos(theta)) and third (phi) elements of 'a' with a randomly, uniformly sampled values of spherical angles cos(theta) and phi (phi normalised by pi).
BoltzmannCollisionIntegral.STMonteCarloAxi_MultiThread!
— MethodSTMonteCarloAxi_MultiThread!(SAtotal3,SAtotal4,TAtotal,SAtally3,SAtally4,TAtally,p3Max,u3MinMax,p4Max,u4MinMax,sigma,dsigmadt,Parameters,numTiterPerThread,numSiterPerThread)
Arguments
SAtotal3::Array{Float64,6}
: Array of stored integration totals for S matrix for 12->34 interactionSAtotal4::Array{Float64,6}
: Array of stored integration totals for S matrix for 12->43 interactionTAtotal::Array{Float64,4}
: Array of stored integration totals for T matrixSAtally3::Array{UInt32,5}
: Array of stored integration tallies for S matrix for 12->34 interactionSAtally4::Array{UInt32,5}
: Array of stored integration tallies for S matrix for 12->43 interactionTAtally::Array{UInt32,4}
: Array of stored integration tallies for T matrixp3Max::Array{Float64,5}
: Array of maximum momentum values for species 3u3MinMax::Array{Float64,6}
: Array of minimum and maximum theta values for species 3p4Max::Array{Float64,5}
: Array of maximum momentum values for species 4u4MinMax::Array{Float64,6}
: Array of minimum and maximum theta values for species 4sigma::Function
: Cross section function for the interactiondsigmadt::Function
: Differential cross section function for the interactionParameters::Tuple{Float64,Float64,Float64,Float64,Float64,Float64,Int64,Float64,Float64,Int64,Float64,Float64,Int64,Float64,Float64,Int64,Int64,Int64,Int64,Int64}
: Tuple of parameters for the interactionnumTiterPerThread::Int64
: Number of T iterations per threadnumSiterPerThread::Int64
: Number of S iterations per thread
Output:
- Argument arrays SAtotal,TAtotal,SAtally,TAtally are mutated to include the results of the Monte Carlo Integration.
Calculation In Brief
- Set up worker threads
- Random Sample points in each of these domains
- RandomPointSphere for theta and phi (for species 1,2,3,4)
- RandomPointMomentum for p ( species 1,2 only)
- Take random points (u3,h3,p1,p2,u1,u2,h1,h2) and calculate valid p3 point/points
- Find position in local S and T arrays and allocated tallies and totals accordingly.
- Take random points (u4,h3,p1,p2,u1,u2,h1,h2) and calculate valid p4 point/points
- Find position in local S and T arrays and allocated tallies and totals accordingly.
- Update global S and T arrays with locks to prevent data races
BoltzmannCollisionIntegral.STMonteCarloAxi_Serial!
— MethodSTMonteCarloAxi_Serial!(SAtotal3,SAtotal4,TAtotal,SAtally3,SAtally4,TAtally,p3Max,p4Max,t3MinMax,t4MinMax,sigma,dsigmadt,Parameters,numTiter,numSiter)
Arguments
SAtotal3::Array{Float64,6}
: Array of stored integration totals for S matrix for 12->34 interactionSAtotal4::Array{Float64,6}
: Array of stored integration totals for S matrix for 12->43 interactionTAtotal::Array{Float64,4}
: Array of stored integration totals for T matrixSAtally3::Array{UInt32,5}
: Array of stored integration tallies for S matrix for 12->34 interactionSAtally4::Array{UInt32,5}
: Array of stored integration tallies for S matrix for 12->43 interactionTAtally::Array{UInt32,4}
: Array of stored integration tallies for T matrixp3Max::Array{Float64,5}
: Array of maximum momentum values for species 3t3MinMax::Array{Float64,6}
: Array of minimum and maximum theta values for species 3p4Max::Array{Float64,5}
: Array of maximum momentum values for species 4t4MinMax::Array{Float64,6}
: Array of minimum and maximum theta values for species 4sigma::Function
: Cross section function for the interactiondsigmadt::Function
: Differential cross section function for the interactionParameters::Tuple{Float64,Float64,Float64,Float64,Float64,Float64,Int64,Float64,Float64,Int64,Float64,Float64,Int64,Float64,Float64,Int64,Int64,Int64,Int64,Int64}
: Tuple of parameters for the interactionnumTiter::Int64
: Number of T iterationsnumSiter::Int64
: Number of S iterations
Output:
- Argument arrays SAtotal,TAtotal,SAtally,TAtally are mutated to include the results of the Monte Carlo Integration.
Calculation In Breif
- Random Sample points in each of these domains
- RandomPointSphere for theta and phi (for species 1,2,3,4)
- RandomPointMomentum for p ( species 1,2 only)
- Take random points (u3,h3,p1,pu2_grid2,u1,u2,h1,h2) and calculate valid p3 point/points
- Find position in local S and T arrays and allocated tallies and totals accordingly.
- Take random points (u4,h3,p1,p2,u1,u2,h1,h2) and calculate valid p4 point/points
- Find position in local S and T arrays and allocated tallies and totals accordingly.
BoltzmannCollisionIntegral.STSymmetry!
— MethodSTSymmetry!(SMatrix3,SMatrix4,TMatrix,u3_bounds,mu1,mu2)
To follow 'PhaseSpaceFactors1'. Physical nature of binary interaction has certain symmetries. 'STSymmetry' uses these symmetries to improve MC sampling of 'SMatrix' and 'TMatrix'.
BoltzmannCollisionIntegral.SValue3
— MethodSValue3(p3v,p1v,p2v,dsigmadt,mu1,mu2,mu3)
Returns Sval
from MC integration based on initial momentum states p1v
and p2v
and final state p3v
and differential cross section dsigmadt
based on particle selection 12->34. Assumes f(x,p,μ)=constant over bin
BoltzmannCollisionIntegral.SValue4
— MethodSValue4(p3v,p1v,p2v,dsigmadt,mu1,mu2,mu3,mu4)
Returns Sval
from MC integration based on initial momentum states p1v
and p2v
and final state p4v
and differential cross section dsigmadt
based on particle selection 12->34. Assumes f(x,p,μ)=constant over bin
BoltzmannCollisionIntegral.SpectraEvaluateMultiThread
— MethodSpectraEvaluateMultiThread(userInputMultiThread)
Function to run the Monte Carlo integration of the S and T arrays in a multi-threaded environment. The function will run the Monte Carlo integration in parallel across the number of threads specified in the global variable nThreads. The function will then calculate the S and T matrices and save the results to a file.
BoltzmannCollisionIntegral.SpectraEvaluateSerial
— MethodSpectraEvaluateSerial(userInputSerial)
Function to run the Monte Carlo integration of the S and T arrays in a serial environment. The function will run the Monte Carlo integration in serial and then calculate the S and T matrices and save the results to a file.
BoltzmannCollisionIntegral.SyncEvaluateMultiThread
— MethodSyncEvaluateMultiThread!(userSyncInputMultiThread)
Function to run the Monte Carlo integration of the S array in a serial enviroment.
BoltzmannCollisionIntegral.SyncEvaluateSerial
— MethodSyncEvaluateSerial!(userSyncInputSerial)
Function to run the Monte Carlo integration of the S array in a serial environment.
BoltzmannCollisionIntegral.SyncKernel
— MethodSyncKernel(p1v,p2v,m2,z2,B)
Returns the emission rate for a single photon $p1v$ state emitted by a charged particle in state $p2v$ with charge $z2$ relative to the fundamental charge and mass $m2$ relative to the mass of the electron, in a uniform magnetic field $B$.
BoltzmannCollisionIntegral.SyncMonteCarloAxi_MultiThread!
— MethodSyncMonteCarloAxi_MultiThread!(SAtotal,SAtally,pMax,tMinMax,Parameters,numSiter)
Arguments
SAtotal::Array{Float64,4}
: Array of stored integration totals for S matrix for 2+B->2+B+1 interactionSAtally::Array{UInt32,4}
: Array of stored integration tallies for S matrix for 2+B->2+B+1 interactionpMax::Array{Float64,3}
: Array of maximum momentum values for species 2tMinMax::Array{Float64,3}
: Array of minimum and maximum theta values for species 2Parameters::Tuple{Float64,Float64,Float64,Int64,Float64,Float64,Int64,Int64,Int64}
: Tuple of parameters for the interactionnumSiter::Int64
: Number of S iterations
Output:
- Argument arrays SAtotal,SAtally are mutated to include the results of the Monte Carlo Integration.
Calculation In Brief
- Random Sample points in each of these domains
- RandomPointSphere for theta (for species 1,2)
- RandomPointMomentum for p ( species 1,2)
- Take random points (p1,p2,t1,t2) and calculate Synchrotron emissivity
- Find position in S arrays and allocated tallies and totals accordingly.
BoltzmannCollisionIntegral.SyncMonteCarloAxi_Serial!
— MethodSyncMonteCarloAxi_Serial!(SAtotal,SAtally,pMax,tMinMax,Parameters,numSiter)
Arguments
SAtotal::Array{Float64,4}
: Array of stored integration totals for S matrix for 2+B->2+B+1 interactionSAtally::Array{UInt32,4}
: Array of stored integration tallies for S matrix for 2+B->2+B+1 interactionpMax::Array{Float64,3}
: Array of maximum momentum values for species 2tMinMax::Array{Float64,3}
: Array of minimum and maximum theta values for species 2Parameters::Tuple{Float64,Float64,Float64,Int64,Float64,Float64,Int64,Int64,Int64}
: Tuple of parameters for the interactionnumSiter::Int64
: Number of S iterations
Output:
- Argument arrays SAtotal,SAtally are mutated to include the results of the Monte Carlo Integration.
Calculation In Brief
- Random Sample points in each of these domains
- RandomPointSphere for theta (for species 1,2)
- RandomPointMomentum for p ( species 1,2)
- Take random points (p1,p2,t1,t2) and calculate Synchrotron emissivity
- Find position in S arrays and allocated tallies and totals accordingly.
BoltzmannCollisionIntegral.SyncSymmetry!
— MethodSyncSymmetry!(SMatrix)
To follow 'PhaseSpaceFactorsSync1'. Synchrotron emission has a symmetry with respect to cos(theta) -> -cos(theta) for both initial particle and photon momenta.
BoltzmannCollisionIntegral.TValue
— MethodTValue(p1v,p2v,sigma,mu1,mu2)
returns Tval
with its Tval from MC integration based on initial momentum states p1v
and p2v
and cross section sigma
based on particle selection. If initial state fails sCheck
, i.e. cannot generate a physical output state, Tval is set to 0e0. Assumes f(x,p,μ)=constant over bin
BoltzmannCollisionIntegral.bounds
— Methodbounds(up_bound,low_bound,num,spacing)
Returns a num+1
long Vector{Float}
of grid bounds. These grid bounds can spaced either by - linear spacing: spacing = "u"
- log10 spacing: spacing = "l"
- binary (1/2^n) spacing: spacing = "b"
BoltzmannCollisionIntegral.bounds_p
— Methodbounds_p(pl,pu,nump)
Returns a nump+1
long Vector{Float}
of p-space grid bounds NOT in Log10 space.
Examples
julia> bounds_p(-5e0,4e0,9)
10-element Vector{Float64}:
1.0e-5
1.0e-4
1.0e-3
0.01
0.1
1.0
10.0
100.0
1000.0
10000.0
BoltzmannCollisionIntegral.bounds_t
— Methodbounds_t(numt)
Returns a numt+1
long Vector{Float}
of theta-space grid bounds in terms of cos(theta). Upper and lower bounds [tl tu] are defined as CONST in Init.jl as [-1 1], type returned is that of tl, tu.
Examples
julia> bounds_t(8)
9-element Vector{Float64}:
-1.0
-0.75
-0.5
-0.25
0.0
0.25
0.5
0.75
1.0
BoltzmannCollisionIntegral.deltaEVector
— MethoddeltaEVector(pr,mu)
Inputs a num+1
long Vector{Float}
of p grid boundaries and the particle mu
value (normalised mass) and returns a num
long Vector{Float}
of average energy values per grid cell.
Examples
julia> deltaEVector([1.0e0, 10.0e0, 100.0e0, 1000.0e0], 1.0e0)
3-element Vector{Float64}:
50.600693
4951.15
495001.16
BoltzmannCollisionIntegral.deltaEkinVector
— MethoddeltaEkinVector(pr,mu)
Inputs a num+1
long Vector{Float}
of p grid boundaries and the particle mu
value (normalised mass) and returns a num
long Vector{Float}
of average kinetic energy values per grid cell.
Examples
julia> deltaEkinVector([1.0e0, 10.0e0, 100.0e0, 1000.0e0], 1.0e0)
3-element Vector{Float64}:
46.10069600605712
4906.1506753523645
494551.15128635924
BoltzmannCollisionIntegral.deltaVector
— MethoddeltaVector(valr)
Inputs a num+1
long Vector{Float}
quantity values (domain bounds) and returns a num
long Vector{Float}
of differences (domain widths).
Examples
julia> deltaVector([1.0e0, 10.0e0, 100.0e0, 1000.0e0])
3-element Vector{Float64}:
9.0
90.0
900.0
BoltzmannCollisionIntegral.dsigmadt_ElePhoElePho
— Methoddsigmadt_ElePhoElePho(sSmol,sBig,tSmol,tBig,uSmol,uBig)
returns the differential cross section for electron-photon scattering (Compton) scattering. Berestetskii 1982 (86.6). Masses and momenta are normalised by the rest mass of the electron $m_{\text{Ele}}$ and the cross section is normalised by $σ_T$.
\[\frac{d\sigma_{e\gamma\rightarrow e\gamma}}{dt}(s,t)=\frac{3}{(s-1)^2}\left[\left(\frac{1}{s-1}+\frac{1}{u-1}\right)^2+\left(\frac{1}{s-1}+\frac{1}{u-1}\right)-\frac{1}{4}\left(\frac{s-1}{u-1}+\frac{u-1}{s-1}\right)\right]\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2 = 1 ∴ s = sSmol + 1$tSmol::Float64
: $t - tBig$tBig::Float64
: $(m_3-m_1)^2 = 0 ∴ t = tSmol$uSmol::Float64
: $u - uBig$uBig::Float64
: $(m_2-m_3)^2 = 1 ∴ u = uSmol + 1$
BoltzmannCollisionIntegral.dsigmadt_ElePosPhoPho
— Methoddsigmadt_ElePosPhoPho(sSmol,sBig,tSmol,tBig,uSmol,uBig)
returns the differential cross section for electron positron annihilation to two photons. Berestetskii 1982 (88.4). Masses and momenta are normalised by the rest mass of the electron $m_{\text{Ele}}$ and the cross section is normalised by $σ_T$.
\[\frac{dσ_{e^+e^-\rightarrow\gamma\gamma}}{dt} = -\frac{3}{s(s-4)}\left(\left(\frac{1}{t-1}+\frac{1}{u-1}\right)^2+\left(\frac{1}{t-1}+\frac{1}{u-1}\right)-\frac{1}{4}\left(\frac{t-1}{u-1}+\frac{u-1}{t-1}\right)\right)\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2 = 4 ∴ s = sSmol + 4$tSmol::Float64
: $t - tBig$tBig::Float64
: $(m_3-m_1)^2 = 1 ∴ t = tSmol + 1$uSmol::Float64
: $u - uBig$uBig::Float64
: $(m2-m3)^2 = 1 ∴ u = uSmol + 1$
BoltzmannCollisionIntegral.dsigmadt_PhoPhoElePos
— Methoddsigmadt_PhoPhoElePos(sSmol,sBig,tSmol,tBig,uSmol,uBig)
returns the differential cross section for photon-photon annihilation to electron-positron pair. (Inverse proceess of electron positron annihilation to two photons). Masses and momenta are normalised by the rest mass of the electron $m_{\text{Ele}}$ and the cross section is normalised by $σ_T$.
\[\frac{dσ_{\gamma\gamma\rightarrow e^+e^-}}{dt} = -\frac{3}{s^2}\left(\left(\frac{1}{t-1}+\frac{1}{u-1}\right)^2+\left(\frac{1}{t-1}+\frac{1}{u-1}\right)-\frac{1}{4}\left(\frac{t-1}{u-1}+\frac{u-1}{t-1}\right)\right)\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2 = 0 ∴ s = sSmol$tSmol::Float64
: $t - tBig$tBig::Float64
: $(m_3-m_1)^2 = 1 ∴ t = tSmol + 1$uSmol::Float64
: $u - uBig$uBig::Float64
: $(m_2-m_3)^2 = 1 ∴ u = uSmol + 1$
BoltzmannCollisionIntegral.dsigmadt_SphSphSphSph
— Methoddsigmadt_SphSphSphSph(sSmol,sBig,tSmol,tBig,uSmol,uBig)
returns the differential cross section for the binary interaction of hard spheres with normalised masses $m_1,m_2,m_3,m_4=m_{\text{Sph}}$. Normalised by $πR_{Sph}^2$.
\[\frac{dσ}{dt} = \frac{1}{s-4m_{\text{Sph}}^2}\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2=4m_{\text{Sph}}^2$tSmol::Float64
: $t - tBig$tBig::Float64
: $(m_3-m_1)^2=0$uSmol::Float64
: $u - uBig$uBig::Float64
: $(m_2-m_3)^2=0$
BoltzmannCollisionIntegral.fload_All
— Methodfload_All(fileLocation,fileName)
Loads all the data stored in fileName
stored at fileLocation
.
Example
(Parameters,SAtot3,SAtot4,TAtot,SAtal3,SAtal4,TAtal,SMatrix3,SMatrix4,TMatrix1,TMatrix2,p3Max,p4Max,u3MinMax,u4MinMax,SConv3,SConv4,TConv) = fload_All(fileLocation,fileName);
Returns a tuple of the data stored in the file. The fields are as follows:
Parameters
: A tuple of the parameters used in the evaluation.Stot3
: A 6D matrix totalling all the emission spectrum values sampled for 12->34 interaction.Stot4
: A 6D matrix totalling all the emission spectrum values sampled for 12->43 interaction.Ttot
: A 4D matrix totalling all the absorption spectrum values sampled.Stal3
: A 5D matrix of tallies of the number of emission spectrum values sampled for 12->34 interaction.Stal4
: A 5D matrix of tallies of the number of emission spectrum values sampled for 12->43 interaction.Ttal
: A 4D matrix of tallies of the number of absorption spectrum values sampled.SMatrix3
: A 6D matrix of the emission spectrum for 12->34 interaction.SMatrix4
: A 6D matrix of the emission spectrum for 12->43 interaction.TMatrix1
: A 4D matrix of the absorption spectrum for 12->34 interaction.TMatrix2
: A 4D matrix of the absorption spectrum for 21->34 interaction i.e. by permutation of TMatrix1 and correct application of phase space factors if species 1 != species 2.p3Max
: The maximum value of the momentum space variable p3 sampled for each bin. (Useful for correcting numerical diffusion)u3MinMax
: The minimum and maximum values of the momentum space variable t3 sampled for each bin. (Useful for correcting numerical diffusion)p4Max
: The maximum value of the momentum space variable p4 sampled for each bin. (Useful for correcting numerical diffusion)u4MinMax
: The minimum and maximum values of the momentum space variable t4 sampled for each bin. (Useful for correcting numerical diffusion)SConv3
: A 6D matrix of the convergence of the emission spectrum compared to the previous run with givenParameters
for 12->34 interaction.SConv4
: A 6D matrix of the convergence of the emission spectrum compared to the previous run with givenParameters
for 12->43 interaction.TConv
: A 4D matrix of the convergence of the absorption spectrum compared to the previous run with givenParameters
.
BoltzmannCollisionIntegral.fload_All_Sync
— Methodfload_All_Sync(fileLocation,fileName)
Loads all the data stored in fileName
stored at fileLocation
.
Example
(Run_Parameters,SAtot,SAtal,SMatrix,#=pMax,tMinMax,=#SConv) = fload_All_Sync(fileLocation,fileName);
Returns a tuple of the data stored in the file. The fields are as follows:
Run_Parameters
: A tuple of the parameters used in the evaluation.Stot
: A 4D matrix totalling all the synchrotron emission spectrum valuesStal
: A 4D matrix of tallies of the number of synchrotron emission spectrum values sampledSMatrix
: A 4D matrix of the synchrotron emission spectrum.pMax
: The maximum value of the momentum space variable p1 (photon mommentum) sampled for each bin. (Useful for correcting numerical diffusion)tMinMax
: The minimum and maximum values of the momentum space variable t1 sampled for each bin. (Useful for correcting numerical diffusion)SConv
: A 4D matrix of the convergence of the synchrotron emission spectrum compared to the previous run with givenRun_Parameters
.
BoltzmannCollisionIntegral.fload_Matrix
— Methodfload_Matrix(fileLocation,fileName)
Loads just the S and T Matrices stored in fileName
stored at fileLocation
.
Example
Matrices = fload_All(fileLocation,fileName);
Returns a tuple of the data stored in the file. The fields are as follows:
Parameters
: A tuple of the parameters used in the evaluation.SMatrix3
: A 6D matrix of the emission spectrum for 12->34 interaction.SMatrix4
: A 6D matrix of the emission spectrum for 12->43 interaction.TMatrix1
: A 4D matrix of the absorption spectrum for 12->34 interaction.TMatrix2
: A 4D matrix of the absorption spectrum for 21->34 interaction.
If initial or final particles are identical then only one of the SMatrices or TMatrices will be returned for that state.
BoltzmannCollisionIntegral.fload_Matrix_ISO
— Methodfload_Matrix_ISO(fileLocation,fileName)
Loads just the S and T Matrices stored in fileName
stored at fileLocation
first converting them to an isotropic form by summing over angles. (The dimensions of the matrices stay the same i.e. 6D->6D with three dimensions having a size of 1)
Example
Matrices = fload_All_ISO(fileLocation,fileName);
Returns a tuple of the data stored in the file. The fields are as follows:
Parameters
: A tuple of the parameters used in the evaluation.SMatrix3
: A 6D matrix of the emission spectrum for 12->34 interaction.SMatrix4
: A 6D matrix of the emission spectrum for 12->43 interaction.TMatrix1
: A 4D matrix of the absorption spectrum for 12->34 interaction.TMatrix2
: A 4D matrix of the absorption spectrum for 21->34 interaction.
If initial or final particles are identical then only one of the SMatrices or TMatrices will be returned for that state.
BoltzmannCollisionIntegral.fload_Matrix_Sync
— Methodfload_Matrix_Sync(fileLocation,fileName)
Loads just the S and T Matrices stored in fileName
stored at fileLocation
.
Example
Matrices = fload_Matrix_Sync(fileLocation,fileName);
Returns a tuple of the data stored in the file. The fields are as follows:
SMatrix
: A 4D matrix of the emission spectrum for Synchrotron.
BoltzmannCollisionIntegral.fload_Matrix_SyncISO
— Methodfload_Matrix_SyncISO(fileLocation,fileName)
Loads just the S and T Matrices stored in fileName
stored at fileLocation
first converting them to an isotropic form by summing over angles. (The dimensions of the matrices stay the same i.e. 6D->6D with three dimensions having a size of 1)
Example
Matrices = fload_Matrix_SyncISO(fileLocation,fileName);
Returns a tuple of the data stored in the file. The fields are as follows:
SMatrix
: A 4D matrix of the emission spectrum for Synchrotron.
BoltzmannCollisionIntegral.location
— Methodlocation(low_bound,up_bound,num,val,spacing)
Returns the index of the bin in which 'val' is contained based the grid bounds of that variable with 'num' bins.
Implemented grid spacing types are: - uniform spacing: spacing = "u"
- log10 spacing: spacing = "l"
- up and low bounds should be supplied as log10 values - binary (1/2^n) spacing: spacing = "b"
- binary spacing is used for u=cos(theta) grids and therefore bounds should always be [-1 1] and num must be odd!
Examples
julia> location(0e0,10e0,9,2e0,"u")
2
BoltzmannCollisionIntegral.location_p
— Methodlocation_p(u,l,num,val)
Returns the index of the momentum bin in which 'val' is contatined based on the 'num' of bins and their 'u' upper and 'l' lower bound including overflow and underflow possibilities. Overflow are assigned to num+1 while underflow are assigned to lowest bin i.e. 1.
Examples
julia> location_p(10e0,1e0,9,2e0)
2
julia> location_p(10e0,1e0,9,11e0) # overflow
10
julia> location_p(10e0,1e0,9,0.5e0) # underflow
1
BoltzmannCollisionIntegral.location_t
— Methodlocation_t(numt,val)
Returns the index of the bin in which the costheta 'val' is contatined based on the 'numt' of bins. Bounds [tl tu] are defined as CONST in Init.jl
Examples
julia> location_t(8,0.5e0)
6
BoltzmannCollisionIntegral.meanVector
— MethodmeanVector(valr)
Inputs a num+1
long Vector{Float}
of domain bounds and returns a num
long Vector{Float}
of mean value in domain range.
Examples
julia> meanVector([1.0e0, 10.0e0, 100.0e0, 1000.0e0])
3-element Vector{Float64}:
5.5
55.0
550.0
BoltzmannCollisionIntegral.p4Vector!
— Methodp4Vector!(p4v,p3v,p1v,p2v)
Returns the p4 vector (in standard form [p,cos(theta),phi/pi]) given the p1, p2 and p3 vectors using conservation of momentum.
BoltzmannCollisionIntegral.sCheck
— MethodsCheck(sSmol,sBig,mu3,mu4)
Returns 'true' if 's' mandelstram generated from inital system state can generate a physical output state.
BoltzmannCollisionIntegral.sigma_ElePhoElePho
— Methodsigma_ElePhoElePho(sSmol,sBig)
returns the total cross section for electron-photon (Compton) scattering. Berestetskii 1982 (86.16). Masses and momenta are normalised by the rest mass of the electron $m_{\text{Ele}}$ and the cross section is normalised by $σ_T$.
\[\sigma_{e\gamma\rightarrow e\gamma}(s)=\frac{3}{4(s-1)}\left[(1-\frac{4}{\left(s-1\right)}-\frac{8m_e^4}{\left(s-1\right)^2})\log\left(s\right)+\frac{1}{2}+\frac{8}{s-1}-\frac{1}{2s^2}\right]\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2 = 1 ∴ s = sSmol + 1$
BoltzmannCollisionIntegral.sigma_ElePosPhoPho
— Methodsigma_ElePosPhoPho(sSmol,sBig)
returns the total cross section for electron positron annihilation to two photons. Berestetskii 1982 (88.6). Masses and momenta are normalised by the rest mass of the electron $m_{\text{Ele}}$ and the cross section is normalised by $σ_T$.
\[σ_{e^+e^-\rightarrow\gamma\gamma} = \frac{3}{4s^2(s-4)}\left((s^2+4s-8)\log\left(\frac{\sqrt{s}+\sqrt{s-4}}{\sqrt{s}-\sqrt{s-4}}\right)-(s+4)\sqrt{s(s-4)}\right)\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2 = 4 ∴ s = sSmol + 4$
BoltzmannCollisionIntegral.sigma_PhoPhoElePos
— Methodsigma_PhoPhoElePos(sSmol,sBig)
returns the total cross section for photon-photon annihilation to electron-positron pair. Masses and momenta are normalised by the rest mass of the electron $m_{\text{Ele}}$ and the cross section is normalised by $σ_T$.
\[σ_{\gamma\gamma\rightarrow e^+e^-} = \frac{3}{2s^3}\left((s^2+4s-8)\log\left(\frac{\sqrt(s)+\sqrt(s-4)}{\sqrt(s)-\sqrt(s-4)}\right)-(s+4)\sqrt{s(s-4)}\right)\]
Arguments
sSmol::Float64
: $s - sBig$sBig::Float64
: $(m_1+m_2)^2 = 0 ∴ s = sSmol$
BoltzmannCollisionIntegral.sigma_SphSphSphSph
— Methodsigma_SphSphSphSph(sSmol,sBig)
returns the total cross section for the binary interaction of hard spheres with normalised masses (wrt electron mass) $m_1,m_2,m_3,m_4=m_\text{Sph}$. Normalised by $πR_{Sph}^2$.
\[σ = \frac{1}{2}\]
Arguments
sSmol::Float64
: s - sBigsBig::Float64
: (m1+m2)^2
BoltzmannCollisionIntegral.vectorLocation
— MethodvectorLocation(pu,pl,nump,numt,vector)
Returns a tuple of bin location for (log10momentum,cos(theta)) based on an input 'vector' and bounds 'u,l' of their domains and the 'num' of uniformly spaced bins. costheta bounds [tl tu] are defined as CONST in Init.jl
Examples
```julia-repl julia> vectorLocation(4e0,-5e0,9,8,[1e0,0.5e0,1.5e0]) (5,6)