API Reference
Simulation
SimulationAnalysis.SingleComponentSimulation
— TypeSingleComponentSimulation
A mutable struct that holds all data for a single-component particle simulation.
Fields
N::Int64
: Number of particles.Ndims::Int64
: Number of spatial dimensions.Nt::Int64
: Number of time steps.dt::Float64
: Time step size.r_array::Array{Float64, 3}
: Particle positions. The dimensions of the array are(N, Ndims, N_time_steps)
.v_array::Array{Float64, 3}
: Particle velocities. The dimensions of the array are(N, Ndims, N_time_steps)
.F_array::Array{Float64, 3}
: Forces on particles. The dimensions of the array are(N, Ndims, N_time_steps)
.D_array::Array{Float64, 1}
: Diameter for each particle.t_array::Array{Float64, 1}
: Vector of time points corresponding to each time step.box_sizes::Vector{Float64}
: A vector holding the simulation box dimensions, e.g.,[Lx, Ly, Lz]
.dt_array::Array{Float64, 1}
: An array of time step differences (Δt
) used for calculating time correlation functions.t1_t2_pair_array::Vector{Array{Int64, 2}}
: An array of(t1, t2)
index pairs, used for efficient calculation of time correlations.filepath::String
: Path to the file from which the simulation data was loaded.
SimulationAnalysis.MultiComponentSimulation
— TypeMultiComponentSimulation
A mutable struct that holds all data for a multi-component particle simulation.
Fields
N::Int64
: Total number of particles across all species.Ndims::Int64
: Number of spatial dimensions.N_species::Int64
: Number of different particle species.Nt::Int64
: Number of time steps.dt::Float64
: Time step size.N_particles_per_species::Vector{Int}
: A vector containing the number of particles for each species.r_array::Vector{Array{Float64, 3}}
: A vector of particle position arrays. Each elementr_array[i]
is a(N_particles_per_species[i], Ndims, N_time_steps)
array for speciesi
.v_array::Vector{Array{Float64, 3}}
: A vector of particle velocity arrays for each species.F_array::Vector{Array{Float64, 3}}
: A vector of particle force arrays for each species.t_array::Vector{Float64}
: Vector of time points corresponding to each time step.box_sizes::Vector{Float64}
: A vector holding the simulation box dimensions, e.g.,[Lx, Ly, Lz]
.dt_array::Vector{Float64}
: An array of time step differences (Δt
) used for calculating time correlation functions.t1_t2_pair_array::Vector{Array{Int64, 2}}
: An array of(t1, t2)
index pairs, used for efficient calculation of time correlations.filepath::String
: Path to the file from which the simulation data was loaded.
SimulationAnalysis.SelfPropelledVoronoiSimulation
— TypeSelfPropelledVoronoiSimulation
A mutable struct that holds all data for a self-propelled Voronoi particle simulation.
Fields
N::Int64
: Number of particles.Ndims::Int64
: Number of spatial dimensions.Nt::Int64
: Number of time steps.dt::Float64
: Time step size.r_array::Array{Float64, 3}
: Particle positions.u_array::Array{Float64, 2}
: Particle orientations (radians).F_array::Array{Float64, 3}
: Forces on particles.perimeter_array::Array{Float64, 2}
: Perimeters of Voronoi cells.area_array::Array{Float64, 2}
: Areas of Voronoi cells.Epot_array::Array{Float64, 1}
: Potential energies.t_array::Vector{Float64}
: Vector of time points.box_sizes::Vector{Float64}
: Dimensions of the simulation box.dt_array::Vector{Float64}
: Time step sizes.t1_t2_pair_array::Vector{Array{Int64, 2}}
: Array of(t1, t2)
index pairs.filepath::String
: Path to the file from which the simulation data was loaded.
SimulationAnalysis.MCSPVSimulation
— TypeMCSPVSimulation
A mutable struct that holds all data for a multi-component self-propelled Voronoi particle simulation.
Fields
N::Int64
: Number of particles.Ndims::Int64
: Number of spatial dimensions.N_species::Int64
: Number of different particle species.Nt::Int64
: Number of time steps.dt::Float64
: Time step size.N_particles_per_species::Vector{Int}
: A vector containing the number of particles for each species.r_array::Vector{Array{Float64, 3}}
: A vector of particle position arrays. Each elementr_array[i]
is a(N_particles_per_species[i], Ndims, N_time_steps)
array for speciesi
.u_array::Vector{Array{Float64, 2}}
: A vector of particle orientation arrays for each species.F_array::Vector{Array{Float64, 3}}
: A vector of particle force arrays for each species.perimeter_array::Vector{Array{Float64, 2}}
: A vector of Voronoi cell perimeter arrays for each species.area_array::Vector{Array{Float64, 2}}
: A vector of Voronoi cell area arrays for each species.Epot_array::Vector{Array{Float64, 1}}
: A vector of potential energy arrays for each species.t_array::Vector{Float64}
: Vector of time points corresponding to each time step.box_sizes::Vector{Float64}
: A vector holding the simulation box dimensions, e.g.,[Lx, Ly, Lz]
.dt_array::Vector{Float64}
: An array of time step differences (Δt
) used for calculating time correlation functions.t1_t2_pair_array::Vector{Array{Int64, 2}}
: An array of(t1, t2)
index pairs, used for efficient calculation of time correlations.filepath::String
: Path to the file from which the simulation data was loaded.
K-Space
SimulationAnalysis.KSpace
— TypeKSpace{dims, S, OA}
A struct to hold the reciprocal space (k-space) vectors for a simulation.
This struct is usually created via construct_k_space
.
Fields
s::S
: TheSimulation
object for which the k-space is defined.Nk::Int64
: The total number of k-vectors.k_lengths::Vector{Float64}
: A vector of lengthNk
containing the magnitude of each k-vector.k_array::Matrix{Float64}
: A(dims, Nk)
matrix where each column is a k-vector(kx, ky, kz, ...)
.kfactor::Int64
: An integer scaling factor used to determine the resolution of the k-space grid. A higher factor means a courser grid.cartesian_to_linear::OA
: AnOffsetArray
that maps the Cartesian grid indices(ix, iy, iz, ...)
of a k-vector to its linear index ink_lengths
andk_array
.
SimulationAnalysis.construct_k_space
— Functionconstruct_k_space(s::Simulation, bounds; kfactor=1, negative=false, rectangular=false)
Constructs the k-space for a given simulation.
This function generates a set of k-vectors based on the simulation box and specified parameters. The k-vectors are chosen from a grid with spacing 2π/L * kfactor
, where L
is the box size in each dimension.
Arguments
s::Simulation
: The simulation object (SingleComponentSimulation
orMultiComponentSimulation
).bounds
: A tuple or vector(kmin, kmax)
specifying the range of magnitudes for the k-vectors. k-vectors with magnitudek
such thatkmin < k < kmax
are included.
Keyword Arguments
kfactor::Int=1
: An integer scaling factor for the k-vector grid resolution.negative::Bool=false
: Iftrue
, includes k-vectors with negative components for the firstd-1
dimensions.rectangular::Bool=false
: Iftrue
, a rectangular k-space grid is generated, and all k-vectors within the grid up tokmax
are included, ignoringkmin
. Iffalse
, k-vectors are selected from a spherical shell.
Returns
- A
KSpace
object containing the generated k-vectors and related information.
See also: KSpace
Correlation Function
SimulationAnalysis.find_correlation_function
— Functionfind_correlation_function(A, B, s)
Computes the correlation function between two lists of observables. The observables are assumed to be complex-valued. The correlation function is a complex number. The correlation function is averaged over the first dimension of A and B.
Arguments
A::Array{Complex{Float64}, 2}
: The first list of observables.B::Array{Complex{Float64}, 2}
: The second list of observables.s::Simulation
: The simulation.
Returns
C::Array{Complex{Float64}, 1}
: The correlation function.
SimulationAnalysis.find_correlation_matrix
— Functionfind_correlation_matrix(observables_list1, observables_list2, kspace, s)
Computes the correlation matrix between two lists of observables. The observables are assumed to be complex-valued. The correlation matrix is a matrix of complex numbers.
Arguments
observables_list1::Array{Array{Complex{Float64}, 2}, 1}
: The first list of observables.observables_list2::Array{Array{Complex{Float64}, 2}, 1}
: The second list of observables.s::Simulation
: The simulation.
Returns
matrix::Array{Array{Complex{Float64}, 1}, 2}
: The correlation matrix.
Density Modes
SimulationAnalysis.find_density_modes
— Functionfind_density_modes(s::SingleComponentSimulation, kspace::KSpace; verbose=true)
Calculates the density modes ρ(k,t) = Σ_j exp(i * k ⋅ r_j(t))
for a single-component simulation.
The calculation is performed for all k-vectors in kspace
and for all time steps in the simulation. This function is computationally intensive. The verbose
option can be used to monitor progress.
Arguments
s::SingleComponentSimulation
: The simulation data.kspace::KSpace
: The k-space vectors for which to calculate the density modes.verbose::Bool=true
: Iftrue
, prints progress and performance information to the console.
Returns
SingleComponentDensityModes
: A struct containing the real and imaginary parts of the density modes.
See also: SingleComponentDensityModes
, find_density_modes(::Union{MultiComponentSimulation,MCSPVSimulation}, ::KSpace)
find_density_modes(s::Union{MultiComponentSimulation,MCSPVSimulation}, kspace::KSpace; verbose=true)
Calculates the density modes ρ_s(k,t) = Σ_j exp(i * k ⋅ r_j(t))
for each species s
in a multi-component simulation.
The calculation is performed for all k-vectors in kspace
and for all time steps in the simulation. This function is computationally intensive. The verbose
option can be used to monitor progress.
Arguments
s::Union{MultiComponentSimulation,MCSPVSimulation}
: The simulation data.kspace::KSpace
: The k-space vectors for which to calculate the density modes.verbose::Bool=true
: Iftrue
, prints progress and performance information to the console.
Returns
MultiComponentDensityModes
: A struct containing the real and imaginary parts of the density modes for each species.
See also: MultiComponentDensityModes
, find_density_modes(::SingleComponentSimulation, ::KSpace)
SimulationAnalysis.SingleComponentDensityModes
— TypeSingleComponentDensityModes
A struct to hold the density modes ρ(k,t)
for a single-component simulation.
The density modes are defined as ρ(k,t) = Σ_j exp(i * k ⋅ r_j(t))
, where the sum is over all particles j
. This struct stores the real and imaginary parts of ρ(k,t)
for various k-vectors and time steps.
Fields
Re::Matrix{Float64}
: The real part of the density modes,Σ_j cos(k ⋅ rj)
. Dimensions are(N_timesteps, Nk)
.Im::Matrix{Float64}
: The imaginary part of the density modes,Σ_j sin(k ⋅ rj)
. Dimensions are(N_timesteps, Nk)
.
SimulationAnalysis.MultiComponentDensityModes
— TypeMultiComponentDensityModes
A struct to hold the density modes ρ_s(k,t)
for each species s
in a multi-component simulation.
The density modes for each species are defined as ρ_s(k,t) = Σ_j exp(i * k ⋅ r_j(t))
, where the sum is over all particles j
of species s
. This struct stores the real and imaginary parts of ρ_s(k,t)
for various k-vectors and time steps.
Fields
Re::Vector{Matrix{Float64}}
: A vector whereRe[s]
is the real part of the density modes for speciess
. Each matrix has dimensions(N_timesteps, Nk)
.Im::Vector{Matrix{Float64}}
: A vector whereIm[s]
is the imaginary part of the density modes for speciess
. Each matrix has dimensions(N_timesteps, Nk)
.
Forces
SimulationAnalysis.Weysser
— TypeWeysser <: InteractionPotential
A struct for the Weysser interaction potential.
Fields
ϵ::Float64
: The interaction strength.δ::Float64
: The interaction range.
SimulationAnalysis.Berthier
— TypeBerthier <: InteractionPotential
A struct for the Berthier interaction potential.
Fields
c0::Float64
: The c0 coefficient.c2::Float64
: The c2 coefficient.c4::Float64
: The c4 coefficient.ζ::Float64
: The zeta coefficient.σ_ratio::Float64
: The ratio of the diameters.
SimulationAnalysis.WCA
— TypeWCA <: InteractionPotential
A struct for the WCA interaction potential.
Fields
ϵ::Float64
: The interaction strength.σ::Float64
: The interaction range.
SimulationAnalysis.KAWCA
— TypeKAWCA <: InteractionPotential
A struct for the Kob-Andersen WCA interaction potential.
Fields
ϵ11::Float64
: The interaction strength between particles of type 1.ϵ12::Float64
: The interaction strength between particles of type 1 and 2.ϵ22::Float64
: The interaction strength between particles of type 2.σ11::Float64
: The interaction range between particles of type 1.σ12::Float64
: The interaction range between particles of type 1 and 2.σ22::Float64
: The interaction range between particles of type 2.
SimulationAnalysis.calculate_forces!
— Functioncalculate_forces!(s, forcestype::Bool)
Calculate the forces for a simulation.
This is a placeholder function that throws an error if forcestype
is true.
Arguments
s
: The simulation object.forcestype::Bool
: Whether to calculate forces.
calculate_forces!(s, U; friction=false, cutoff=1.25)
Recalculates the total force on all particles according to the langevin equation F = -∇U - γv. This function also updates the total potential energy in the output datastructure.
calculate_forces!(s, U::KAWCA; cutoff=2.0^(1.0/6.0))
Recalculates the total force on all particles in a multi-component simulation s
according to the Kob-Andersen Weeks-Chandler-Andersen system U, and updates the force array in s
with the calculated values.
Arguments
s
: aMultiComponentSimulation
object that contains the position and force arrays for each particle.U::KAWCA
: the interaction type.cutoff
: the cutoff distance for the potential energy function. Default is2.0^(1.0/6.0)
.
Intermediate Scattering Function
SimulationAnalysis.find_intermediate_scattering_function
— Functionfind_intermediate_scattering_function(s::Simulation; kmin=7.0, kmax=7.4, kfactor=1)
Calculates the coherent intermediate scattering function F(k, t)
for a simulation, averaged over a shell in k-space.
This is a convenience function that first constructs the k-space and density modes, and then computes F(k, t)
. The function is defined as: F(k, t) = (1/N) * < Σ_{i,j} exp(i * k ⋅ (r_i(t) - r_j(0))) >
. The average is taken over time origins and over k-vectors with magnitudes k
such that kmin < k < kmax
.
Arguments
s::Simulation
: The simulation data.kmin::Float64=7.0
: The minimum magnitude of the k-vectors to be included in the average.kmax::Float64=7.4
: The maximum magnitude of the k-vectors to be included in the average.kfactor::Int=1
: The resolution factor for the k-space grid.
Returns
Fk
: The intermediate scattering function. For aSingleComponentSimulation
, this is aVector{Float64}
. For aMultiComponentSimulation
, this is aMatrix{Vector{Float64}}
.
find_intermediate_scattering_function(s::Simulation, kspace::KSpace, ρkt, k_sample_array::AbstractVector; k_binwidth=0.1)
Calculates the coherent intermediate scattering function F(k, t)
for a list of specified k
values.
For each k
in k_sample_array
, this function computes F(k, t)
by averaging over a k-shell of width k_binwidth
centered at k
.
Arguments
s::Simulation
: The simulation data.kspace::KSpace
: The pre-computed k-space.ρkt
: The pre-computed density modes (SingleComponentDensityModes
orMultiComponentDensityModes
).k_sample_array::AbstractVector
: A vector of k-magnitudes for which to computeF(k, t)
.k_binwidth::Float64=0.1
: The width of the k-shell to average over for each value ink_sample_array
.
Returns
F_array::Vector
: A vector whereF_array[i]
is theF(k,t)
corresponding tok_sample_array[i]
.
find_intermediate_scattering_function(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, kspace::KSpace, ρkt::SingleComponentDensityModes; kmin=0.0, kmax=10.0^10.0)
Calculates the coherent intermediate scattering function F(k, t)
for a single-component simulation.
This is the main implementation that computes F(k, t) = (1/N) * <ρ(k, t) ρ*(-k, 0)>
by correlating the pre-computed density modes ρkt
. The average is performed over time origins and k-vectors within the magnitude range [kmin, kmax]
.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.kspace::KSpace
: The pre-computed k-space.ρkt::SingleComponentDensityModes
: The pre-computed density modes.kmin::Float64=0.0
: The minimum magnitude of k-vectors to include in the average.kmax::Float64=10.0^10.0
: The maximum magnitude of k-vectors to include in the average.
Returns
Fk::Vector{Float64}
: A vector containingF(k, t)
for each time delayΔt
ins.dt_array
.
find_intermediate_scattering_function(s::Union{MultiComponentSimulation,MCSPVSimulation}, kspace::KSpace, ρkt::MultiComponentDensityModes; kmin=0.0, kmax=10.0^10.0)
Calculates the partial intermediate scattering functions F_αβ(k, t)
for a multi-component simulation.
This function computes F_αβ(k, t) = (1/N) * <ρ_α(k, t) ρ_β*(-k, 0)>
by correlating the pre-computed density modes ρkt
for species α
and β
. The average is performed over time origins and k-vectors within the magnitude range [kmin, kmax]
.
Arguments
s::Union{MultiComponentSimulation,MCSPVSimulation}
: The simulation data.kspace::KSpace
: The pre-computed k-space.ρkt::MultiComponentDensityModes
: The pre-computed density modes for all species.kmin::Float64=0.0
: The minimum magnitude of k-vectors to include in the average.kmax::Float64=10.0^10.0
: The maximum magnitude of k-vectors to include in the average.
Returns
Fk::Matrix{Vector{Float64}}
: A matrix of vectors, whereFk[α, β]
is the partial intermediate scattering functionF_αβ(k, t)
.
SimulationAnalysis.find_self_intermediate_scattering_function
— Functionfind_self_intermediate_scattering_function(s::Simulation; kmin=7.0, kmax=7.4, kfactor=1)
Calculates the self-intermediate scattering function Fs(k, t)
, averaged over a shell in k-space.
This is a convenience function that first constructs the k-space and then computes Fs(k, t)
. The function is defined as: Fs(k, t) = (1/N) * <Σ_j exp(i * k ⋅ (r_j(t) - r_j(0)))>
. The average is taken over time origins and over k-vectors with magnitudes k
such that kmin < k < kmax
. This function is only implemented for SingleComponentSimulation
.
Arguments
s::Simulation
: The simulation data.kmin::Float64=7.0
: The minimum magnitude of the k-vectors to be included in the average.kmax::Float64=7.4
: The maximum magnitude of the k-vectors to be included in the average.kfactor::Int=1
: The resolution factor for the k-space grid.
Returns
Fk::Vector{Float64}
: A vector containingFs(k, t)
for each time delayΔt
ins.dt_array
.Fks_per_particle::Matrix{Float64}
: A(Ndt, N)
matrix containing the self-intermediate scattering function for each particle.
find_self_intermediate_scattering_function(s::Simulation, kspace::KSpace, k_sample_array::AbstractVector; k_binwidth=0.1)
Calculates the self-intermediate scattering function Fs(k, t)
for a list of specified k
values.
For each k
in k_sample_array
, this function computes Fs(k, t)
by averaging over a k-shell of width k_binwidth
centered at k
. This function is only implemented for SingleComponentSimulation
.
Arguments
s::Simulation
: The simulation data.kspace::KSpace
: The pre-computed k-space.k_sample_array::AbstractVector
: A vector of k-magnitudes for which to computeFs(k, t)
.k_binwidth::Float64=0.1
: The width of the k-shell to average over for each value ink_sample_array
.
Returns
F_array::Vector
: A vector of tuples(Fk, Fks_per_particle)
for eachk
ink_sample_array
.
find_self_intermediate_scattering_function(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, kspace::KSpace; kmin=0.0, kmax=10.0^10.0)
Calculates the self-intermediate scattering function Fs(k, t)
for a single-component simulation.
This is the main implementation that computes Fs(k, t) = (1/N) * <Σ_j exp(i * k ⋅ (r_j(t) - r_j(0)))>
. The calculation is done without pre-calculating density modes to avoid large memory allocation. The average is performed over time origins and k-vectors within the magnitude range [kmin, kmax]
.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.kspace::KSpace
: The pre-computed k-space.kmin::Float64=0.0
: The minimum magnitude of k-vectors to include in the average.kmax::Float64=10.0^10.0
: The maximum magnitude of k-vectors to include in the average.
Returns
Fk::Vector{Float64}
: A vector containingFs(k, t)
for each time delayΔt
ins.dt_array
.Fks_per_particle::Matrix{Float64}
: A(Ndt, N)
matrix containing the self-intermediate scattering function for each particle.
Structure Factors
SimulationAnalysis.find_structure_factor
— Functionfind_structure_factor(s::Simulation; kmin=7.0, kmax=7.4, kfactor=1)
Calculates the static structure factor S(k)
for a simulation, averaged over a shell in k-space.
This is a convenience function that first constructs the k-space and density modes, and then computes S(k)
. The function is defined as: S(k) = (1/N) * <|Σ_j exp(i * k ⋅ r_j)|^2>
. The average is taken over time and over k-vectors with magnitudes k
such that kmin < k < kmax
.
Arguments
s::Simulation
: The simulation data.kmin::Float64=7.0
: The minimum magnitude of the k-vectors to be included in the average.kmax::Float64=7.4
: The maximum magnitude of the k-vectors to be included in the average.kfactor::Int=1
: The resolution factor for the k-space grid.
Returns
Sk
: The structure factor. For aSingleComponentSimulation
, this is aFloat64
. For aMultiComponentSimulation
, this is aMatrix{Float64}
.
find_structure_factor(s::Simulation, kspace::KSpace, ρkt::AbstractDensityModes, k_sample_array::AbstractVector; k_binwidth=0.1)
Calculates the static structure factor S(k)
for a list of specified k
values.
For each k
in k_sample_array
, this function computes S(k)
by averaging over a k-shell of width k_binwidth
centered at k
.
Arguments
s::Simulation
: The simulation data.kspace::KSpace
: The pre-computed k-space.ρkt::AbstractDensityModes
: The pre-computed density modes.k_sample_array::AbstractVector
: A vector of k-magnitudes for which to computeS(k)
.k_binwidth::Float64=0.1
: The width of the k-shell to average over for each value ink_sample_array
.
Returns
S_array::Vector
: A vector whereS_array[i]
is theS(k)
corresponding tok_sample_array[i]
.
find_structure_factor(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, kspace::KSpace, ρkt::SingleComponentDensityModes; kmin=0.0, kmax=10.0^10.0)
Calculates the static structure factor S(k)
for a single-component simulation.
This is the main implementation that computes S(k) = (1/N) * <|ρ(k)|^2>
from the pre-computed density modes ρkt
. The average is performed over time and k-vectors within the magnitude range [kmin, kmax]
.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.kspace::KSpace
: The pre-computed k-space.ρkt::SingleComponentDensityModes
: The pre-computed density modes.kmin::Float64=0.0
: The minimum magnitude of k-vectors to include in the average.kmax::Float64=10.0^10.0
: The maximum magnitude of k-vectors to include in the average.
Returns
Sk::Float64
: The value of the structure factorS(k)
.
find_structure_factor(s::Union{MultiComponentSimulation,MCSPVSimulation}, kspace::KSpace, ρkt::MultiComponentDensityModes; kmin=0.0, kmax=10.0^10.0)
Calculates the partial static structure factors S_αβ(k)
for a multi-component simulation.
This function computes S_αβ(k) = (1/N) * <ρ_α(k) ρ_β*(-k)>
from the pre-computed density modes ρkt
for species α
and β
. The average is performed over time and k-vectors within the magnitude range [kmin, kmax]
.
Arguments
s::Union{MultiComponentSimulation,MCSPVSimulation}
: The simulation data.kspace::KSpace
: The pre-computed k-space.ρkt::MultiComponentDensityModes
: The pre-computed density modes for all species.kmin::Float64=0.0
: The minimum magnitude of k-vectors to include in the average.kmax::Float64=10.0^10.0
: The maximum magnitude of k-vectors to include in the average.
Returns
Sk::Matrix{Float64}
: A matrix whereSk[α, β]
is the partial structure factorS_αβ(k)
.
Mean Squared Displacement
SimulationAnalysis.find_mean_squared_displacement
— Functionfind_mean_squared_displacement(simulation::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}; per_particle=false)
Calculates the mean squared displacement (MSD) for a given simulation.
The MSD is defined as <|r(t) - r(0)|^2>
, where the average is taken over all particles and time origins.
Arguments
simulation::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.per_particle::Bool=false
: Iftrue
, the function returns the MSD calculated for each particle individually, in addition to the total MSD.power::Int=2
: The power to which the displacement is raised. Default is 2 for MSD. Set to 4 for the mean quartic displacement of the displacement for example.
Returns
If per_particle
is false
(default):
msd::Vector{Float64}
: A vector containing the MSD for each time delayΔt
insimulation.dt_array
.
If per_particle
is true
:
msd::Vector{Float64}
: The total MSD vector.msd_per_particle::Matrix{Float64}
: A(N, N_dt)
matrix of MSDs for each particle.
find_mean_squared_displacement(simulation::Union{MultiComponentSimulation, MCSPVSimulation}; per_particle=false)
Calculates the mean squared displacement (MSD) for a given simulation.
The MSD is defined as <|r(t) - r(0)|^2>
, where the average is taken over all particles and time origins.
Arguments
simulation::Union{MultiComponentSimulation, MCSPVSimulation}
: The simulation data.power::Int=2
: The power to which the displacement is raised. Default is 2 for MSD. Set to 4 for the mean quartic displacement of the displacement for example.
Returns
msd::Vector{Float64}
: A vector containing the MSD for each time delayΔt
insimulation.dt_array
.
find_mean_squared_displacement(r, dt_array, t1_t2_pair_array)
Calculates the mean squared displacement (MSD) using pre-calculated time pairs.
This is the core implementation for calculating the MSD. It iterates through particles and time pairs to compute the displacement.
Arguments
r::Array{Float64, 3}
: A(Ndims, N, N_timesteps)
array of particle positions.dt_array::Vector{Int}
: An array of time step differences (Δt
).t1_t2_pair_array::Vector{Array{Int64, 2}}
: An array of(t1, t2)
index pairs for eachΔt
.power::Int=2
: The power to which the displacement is raised.
Returns
msd::Vector{Float64}
: A vector containing the MSD for each time delayΔt
indt_array
.msd_per_particle::Matrix{Float64}
: A(N, N_dt)
matrix of MSDs for each particle.
SimulationAnalysis.find_non_gaussian_parameter
— Functionfind_non_gaussian_parameter(simulation::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}; per_particle=false)
Calculates the non-Gaussian parameter α2 for a given simulation.
Arguments
simulation::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.per_particle::Bool=false
: Iftrue
, the function returns the non-Gaussian parameter calculated for each particle individually, in addition to the total non-Gaussian parameter.
Returns
α2::Float64
: The non-Gaussian parameter for the entire simulation.α2_per_particle::Matrix{Float64}
: A(N, N_dt)
matrix of non-Gaussian parameters for each particle, ifper_particle
istrue
.
Overlap Function
SimulationAnalysis.find_overlap_function
— Functionfind_overlap_function(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}; a=0.5)
Calculates the overlap function Q(t)
.
The overlap function is a measure of the similarity between the particle configurations at time 0
and time t
. It is defined as: Q(t) = (1/N) * Σ_i θ(a - |r_i(t) - r_i(0)|)
, where θ
is the Heaviside step function, a
is a cutoff distance, and the average is taken over all particles i
and time origins. A value of 1 means the configurations are identical (within the cutoff), and a value of 0 means they are completely different.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.a::Float64=0.5
: The cutoff distance for calculating the overlap. Typically this is a fraction of the particle diameter.
Returns
Fs::Vector{Float64}
: A vector containingQ(t)
for each time delayΔt
ins.dt_array
.Fs_pp::Matrix{Float64}
: A(Ndt, N)
matrix containing the overlap function for each particle.
Neighborlists
SimulationAnalysis.find_relative_distance_neighborlists
— Functionfind_relative_distance_neighborlists(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, rc::Float64; ζ::Float64 = 0.2)
Finds neighbor lists for all time steps based on a relative distance criterion.
This method is suitable for polydisperse systems where particle diameters vary. Two particles i
and j
are considered neighbors if the distance between them is less than a cutoff that depends on their respective diameters D_i
and D_j
: distance(i, j) < rc * (D_i + D_j)/2 * (1 - ζ * abs(D_i - D_j))
This function requires the D_array
(particle diameters) to be loaded in the Simulation
object.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data. Must contain particle diameters.rc::Float64
: The relative cutoff distance.ζ::Float64=0.2
: A non-additivity parameter to modulate the cutoff for particles of different sizes.
Returns
neighborlists::Vector{Vector{Vector{Int}}}
: A vector over time steps. For each time step, a vector over particles, whereneighborlists[t][i]
contains the list of neighbors of particlei
.
SimulationAnalysis.find_absolute_distance_neighborlists
— Functionfind_absolute_distance_neighborlists(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, rc::Float64)
Finds neighbor lists for all time steps based on a fixed, absolute distance cutoff.
Two particles i
and j
are considered neighbors if the distance between them is less than rc
. Periodic boundary conditions are taken into account.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.rc::Float64
: The absolute cutoff distance.
Returns
neighborlists::Vector{Vector{Vector{Int}}}
: A vector over time steps. For each time step, a vector over particles, whereneighborlists[t][i]
contains the list of neighbors of particlei
.
SimulationAnalysis.find_voronoi_neighborlists
— Functionfind_voronoi_neighborlists(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}; max_distance_from_boundary=3.0, verbose=true, indices=eachindex(s.t_array))
Computes neighbor lists based on Voronoi tessellation for each time step.
Two particles are considered neighbors if their Voronoi cells are adjacent (i.e., share a facet). This method uses the Quickhull.jl
package to perform the Delaunay triangulation, which is the dual of the Voronoi tessellation.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.max_distance_from_boundary::Float64=3.0
: How far from the central box to consider periodic images for the tessellation. This should be large enough to avoid boundary effects.verbose::Bool=true
: Iftrue
, prints progress information.indices=eachindex(s.t_array)
: The range of time step indices to analyze. Defaults to all time steps.
Returns
neighborlists::Vector{Vector{Vector{Int}}}
: A vector over time steps. For each time step, a vector over particles, whereneighborlists[t][i]
contains the list of neighbors of particlei
.
Bond Breaking Parameter
SimulationAnalysis.find_CB
— Functionfind_CB(s::Simulation, neighbourlists1::Vector{<:Vector}, neighbourlists2::Vector{<:Vector})
Calculates the bond-breaking correlation function C_B(t)
for all particles.
This function measures the fraction of neighbors that a particle at time t_0
still has at a later time t_0 + t
. It is averaged over all particles and time origins t_0
.
Arguments
s::Simulation
: The simulation data, used for time information.neighbourlists1::Vector{<:Vector}
: A vector of neighbor lists at the initial times (t_0
).neighbourlists2::Vector{<:Vector}
: A vector of neighbor lists at the final times (t_0 + t
). Often, this is the same asneighbourlists1
.
Returns
CB::Matrix{Float64}
: A(Ndt, N)
matrix whereCB[i, j]
is the bond-breaking correlation for particlej
at time delaydt_array[i]
.
Utils
SimulationAnalysis.find_relaxation_time
— Functionfind_relaxation_time(t::AbstractVector, F::AbstractVector; threshold::Float64=exp(-1))
Calculates the relaxation time τ
of a decaying function F(t)
.
The relaxation time is defined as the time at which F(t)
first drops below a certain threshold
. The function finds this time by performing linear interpolation between the two points straddling the threshold. The interpolation is linear in F
and log(t)
.
Arguments
t::AbstractVector
: A vector of time points.F::AbstractVector
: A vector of the function valuesF(t)
, corresponding to the time points int
. This function should be monotonically decreasing.threshold::Float64=exp(-1)
: The threshold value forF(t)
. The relaxation time ist
such thatF(t) = threshold
.
Returns
τ::Float64
: The calculated relaxation time. Returns0.0
if the function is already below the threshold at the first time point. ReturnsInf
if the function never drops below the threshold.
Radial Distribution Function
SimulationAnalysis.find_radial_distribution_function
— Functionfind_radial_distribution_function(s::Union{MultiComponentSimulation,MCSPVSimulation}, Nbins::Int, rmax::Float64; Ntasks=1, verbose=true)
Calculates the partial radial distribution functions g_αβ(r)
for a multi-component simulation.
The g_αβ(r)
describes the probability of finding a particle of species β
at a distance r
from a particle of species α
, relative to that of an ideal gas.
Arguments
s::Union{MultiComponentSimulation,MCSPVSimulation}
: The simulation data.Nbins::Int
: The number of bins to use for the histogram.rmax::Float64
: The maximum distance to consider. If negative, it's set to half the smallest box dimension.Ntasks::Int=1
: The number of parallel tasks to use for the calculation.verbose::Bool=true
: Iftrue
, prints progress to the console.
Returns
bin_centres::Vector{Float64}
: A vector of the centre points of the histogram bins.g_r::Matrix{Vector{Float64}}
: A matrix of vectors, whereg_r[α, β]
is the partialg_αβ(r)
.
find_radial_distribution_function(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, Nbins::Int, rmax::Float64; Ntasks::Int=1, verbose=true)
Calculates the radial distribution function g(r)
for a single-component simulation.
The g(r)
describes the probability of finding a particle at a distance r
from another particle, relative to that of an ideal gas.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.Nbins::Int
: The number of bins to use for the histogram.rmax::Float64
: The maximum distance to consider. If negative, it's set to half the smallest box dimension.Ntasks::Int=1
: The number of parallel tasks to use for the calculation.verbose::Bool=true
: Iftrue
, prints progress to the console.
Returns
bin_centres::Vector{Float64}
: A vector of the centre points of the histogram bins.g_r::Vector{Float64}
: A vector containing the values ofg(r)
.
find_radial_distribution_function(s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}, dr::Float64, rmax::Float64; Ntasks::Int=1, verbose=true)
Calculates the radial distribution function g(r)
for a single-component simulation.
This method defines the histogram bins by a fixed width dr
.
Arguments
s::Union{SingleComponentSimulation, SelfPropelledVoronoiSimulation}
: The simulation data.dr::Float64
: The width of the histogram bins.rmax::Float64
: The maximum distance to consider. If negative, it's set to half the smallest box dimension.Ntasks::Int=1
: The number of parallel tasks to use for the calculation.verbose::Bool=true
: Iftrue
, prints progress to the console.
Returns
bin_centres::Vector{Float64}
: A vector of the centre points of the histogram bins.g_r::Vector{Float64}
: A vector containing the values ofg(r)
.