Library

Public API

SuperVFM.GetSchwarzTempCoeffsMethod
GetSchwarzTempCoeffs(Temp::Real)

Obtain the Schwarz α and α' parameters from observational data using the temperature. For values that lie within the observation data, compute a cubic spline interpolation to extract parameters for the desired temperature.

source
SuperVFM.check_timestepMethod
check_timestep(SP::SimulationParams)

Checks if the current timestep in SP is small enough to resolve the smallest Kelvin waves. Returns true if the timestep $\Delta t < \Delta t_{max}$.

The maximum timestep is given by

\[ \Delta t_{max} = \frac{(\delta/2)^2}{\kappa\log(\delta/(2\pi a_0))}\]

source
SuperVFM.computeEnergySpectrumMethod
computeEnergySpectrum(VortexData, kMap, ::Anisotropic)

Full anisotropic computation of the superfluid kinetic energy spectrum $E(k)$. The method uses the vortex data read from file and uses kMap to loop through integer k-shells and the assosciated k vectors to each shell and saves the final array as a 1D k vector.

Warning

The full anistropic form is very slow, and should only be used sparingly when the isotropic version does not suffice. Only use this for the time steps necessary! Consider running on GPU for vastly improved performance.

The full anisotropic formula for the 3D energy spectrum $E(\mathbf{k})$is given by:

\[E(\mathbf{k}) = \frac{\rho_s \kappa^2}{16\pi^3}\frac{1}{k^2}\int_{\mathcal{L}_1}\int_{\mathcal{L}_2}d\xi_1 d\xi_2 \mathbf{s}'(\xi_1)\cdot\mathbf{s}'(\xi_2) e^{i\mathbf{k}\cdot\left[\mathbf{s}(\xi_2) - \mathbf{s}(\xi_1) \right]}.\]

The 1D spectrum is computed by integration over concentric shells

\[E(k) = \sum\limits_{k-1<\left\lVert \mathbf{k}\right\rVert}\leq k E(\mathbf{k}),\]

such that the total superfluid kinetic energy $E$ is

\[E = \int_{\hat{\Omega}}E(\mathbf{k})d^3\mathbf{k} = \int_{0}^{\infty}E(k)dk.\]

source
SuperVFM.computeEnergySpectrumMethod
computeEnergySpectrum(VortexData, kMap, ::Isotropic)

Computes the superfluid kinetic energy spectrum $E(k)$ using the isotropic approximation formula:

\[E(k) = \frac{\rho_s \kappa^2}{4\pi^2}\int_{\mathcal{L}_1}\int_{\mathcal{L}_2}d\xi_1 d\xi_2 \mathbf{s}'(\xi_1)\cdot\mathbf{s}'(\xi_2)\frac{\sin\left(k|\mathbf{s}(\xi_2) -\mathbf{s}(\xi_1)| \right)}{k|\mathbf{s}(\xi_2) -\mathbf{s}(\xi_1)| \right)},\]

such that the total superfluid kinetic energy $E$ is

\[E = \int_{0}^{\infty}E(k)dk.\]

source
SuperVFM.generateMappingMethod
generateMapping(N::Int)

Generate a vector of combination k-maps. A combination map with id=n contains all of the 3D k-vectors $(k_x,k_y,k_z)$ which reside within two consective integer spherical shells of size n.

The mapping partititions the vector as a set $\Omega = \lbrace \Omega_k : 1\leq k \leq M = N/2\rbrace$ where $N$ is the resolution and

\[\Omega_k = \lbrace \mathbf{k} = (k_x,k_y,k_z) : k-1 < ||\mathbf{k}||_{2} \leq k \rbrace\]

source
SuperVFM.get_densityMethod
get_density(T::AbstractFloat)

Returns the normal fluid density $ρ_n$ and superfluid density $ρ_s$ for a given temperature $T$ in arbitrary units of temperature.

Requires $0.0\leq T \leq T_{\lambda} = 2.178$

source
SuperVFM.get_densityMethod
get_density(T::Unitful.Temperature)

Returns the normal fluid density $ρ_n$ and superfluid density $ρ_s$ for a given temperature $T$ in arbitrary units of temperature.

Requires $0.0\leq T \leq T_{\lambda} = 2.178$

source
SuperVFM.get_dynamic_viscosityMethod
get_dynamic_viscosity(T::AbstractFloat)

Returns the dynamic viscosity $η$ at temperature $T$.

Requires $0.8\leq T < T_{\lambda} = 2.178$

source
SuperVFM.get_dynamic_viscosityMethod
get_dynamic_viscosity(T::Unitful.Temperature)

Returns the dynamic viscosity $η$ at temperature $T$.

Requires $0.8\leq T \leq T_{\lambda} = 2.178$

source
SuperVFM.print_characteristicsMethod
print_characteristics(io::IO,SimParams::SimulationParams)

Prints the characteristic time and length scales of the simulation in dimensional units

source
SuperVFM.print_density_dataFunction
print_density_data(io::IO=stdout)

Prints the normal fluid $ρ_n$ and superfluid $ρ_s$densities to the IO buffer. Defaults to stdout.

source
KernelAbstractions.CPUType
CPU(; static=false)

Instantiate a CPU (multi-threaded) backend.

Options:

  • static: Uses a static thread assignment, this can be beneficial for NUMA aware code. Defaults to false.
source
SuperVFM.OpenBoundaryType
OpenBoundary(dim::Int)

Initialises an open boundary in the chosen direction. Vortex loops that exceed the box size (typically $2π$) are not restricted. Set dim=1 for the $x$ direction, dim=2 for the $y$ direction and dim=3 fir the $z$ direction.

Example usage:

    boundary_x = OpenBoundary(1)
    boundary_y = OpenBoundary(2)
    boundary_z = OpenBoundary(3)
source
SuperVFM.PeriodicBoundaryType
PeriodicBoundary(dim::Int)

Initialises a periodic boundary in the chosen direction selected by dim. Vortex loops that exceed the box size (typically $2π$) are looped back periodically to the other side of the box.Set dim=1 for the $x$ direction, dim=2 for the $y$ direction and dim=3 fir the $z$ direction.

Example usage:

    boundary_x = PeriodicBoundary(1)
    boundary_y = PeriodicBoundary(2)
    boundary_z = PeriodicBoundary(3)
Warning

As of v1.0.2, only the periodic boundary condition is fully implemented and working correctly, open boundary conditions and solid walls will be implemented in a future release.

source

Internal API

Base.showMethod
Base.show(io::IO,SimParams::SimulationParams)

Lists the simulation parameters stored in SimParams in a stylistic way with a key.

source
KernelAbstractions.zerosMethod
KA.zeros(BE::Backend,::Type{SVector{S,T}},N::Int) where {S,T}

Initialise an array of static vectors of size S and type T

source
SuperVFM.SchwarzModel_kernel!Method
SchwarzModel_kernel!(u_mf, u_sup, f, f_infront, ghosti, ghostb, FM::SchwarzModel, normal_velocity)

Kernel launch for the Schwarz model.

source
SuperVFM.check_active_pcountMethod
check_active_pcount(f_infront)

Exits the run if there are not enough vortex points, if there are less than 5 points (needed to properly resolve the derivatives).

source
SuperVFM.compute_filament_velocity!Method
compute_filament_velocity!(u, u_loc, u_sup, FM::SchwarzModel, SP::SimulationParams{S,T}; kwargs...) where {S,T}

Compute vortex filament velocities using the Schwarz model.

\[ \frac{d\mathbf{s}}{dt} = \mathbf{v}_s + \alpha\left[\mathbf{s}'\times\left(\mathbf{v}_n - \mathbf{v}_s\right)\right] - \alpha'\left(\mathbf{s}'\times\left[\mathbf{s}'\times\left(\mathbf{v}_n - \mathbf{v}_s\right)\right] \right)\]

source
SuperVFM.compute_filament_velocity!Method
compute_filament_velocity!(u, u_loc, u_sup, ::ZeroTemperature, SP::SimulationParams{S,T}; kwargs...) where {S,T}

Compute vortex filament velocities in the zero temperature limit. Filaments are evolved by

\[ \frac{d\mathbf{s}}{dt} = \mathbf{v}_s\]

source
SuperVFM.compute_velocity!Method
compute_velocity!(u_loc, u_sup, ::LIA, SP::SimulationParams{S,T},; kwargs...) where {S,T}

Computes the superfluid velocity using the Local Induction Approximation (LIA)

\[\mathbf{v}_i = \beta \mathbf{s}'_i \times \mathbf{s}''_i \]

source
SuperVFM.create_info_fileMethod
create_info_file(::SimulationParams{S,T}) where {S,T}

Create file to print the simulation precision of floating point and integers. To be used by vortex reading methods for correct parsing.

source
SuperVFM.enforce_boundary!Method
enforce_boundary!(f, boundary_x::OpenBoundary, boundary_y::OpenBoundary, boundary_z::OpenBoundary; kwargs...)

Enforces open boundary conditions across all 3 dimensions.

source
SuperVFM.enforce_boundary!Method
enforce_boundary!(f, boundary_x::PeriodicBoundary, boundary_y::PeriodicBoundary, boundary_z::PeriodicBoundary; kwargs...)

Enforces periodic boundary conditions across all 3 dimensions.

source
SuperVFM.get_curvature!Method
get_curvature!(curv; kwargs...)

Calculate the vortex line curvature $\zeta$.

\[ \zeta = |\mathbf{s}''| = \left|\frac{d\mathbf{s}}{d\xi}\right|\]

source
SuperVFM.get_deriv_1Method
get_deriv_1(f::AbstractArray, ghosti::AbstractArray, ghostb::AbstractArray)

Computes the first order derivative using a 2nd order finite difference adaptive scheme using ghost points.

source
SuperVFM.get_deriv_2Method
get_deriv_2(f::AbstractArray, ghosti::AbstractArray, ghostb::AbstractArray)

Computes the first order derivative using a 2nd order finite difference adaptive scheme using ghost points.

source
SuperVFM.get_ΔξMethod
get_Δξ(f, ghosti, f_infront, pcount, SP::SimulationParams{S,T}) where {S,T}

Compute the seperation distance $\Delta\xi$ between each filament.

source
SuperVFM.ghostpMethod
ghostp(f, f_infront, f_behind, pcount, SP::SimulationParams{S,T}) where {S,T}

Compute ghost points.

source
SuperVFM.identify_loop!Method

Identify a single vortex loop.

The input n_start indicates from which particle index to start looking for non-visited particles.

Returns the index of a reference particle in the loop and marks all particles in the loop as visited. Also returns the number of particles in the loop.

If no loop is found (because all particles have already been visited), returns zero.

source
SuperVFM.initialiseVortexMethod
initialiseVortex(SP::SimulationParams{S,T}) where {S,T}

Uses SP.backend to initialise the vortex structure according to the initial condition SP.initf.

source
SuperVFM.static2fieldMethod
static2field(field::AbstractArray{SVector{S,T}}, N::Int) where {S<:Int,T<:AbstractFloat}

Convert array of static vectors of size S and type T of length N, to a matrix of with dimensions (S,N). Return type is an child type of AbstractMatrix.

julia> N = 128;

julia> f = rand(SVector{3,Float32},N);

julia> typeof(static2field(f,N)) <: AbstractMatrix
true
source
SuperVFM.timestep!Method
timestep!(f, u, u1, u2, f_infront, pcount, SP::SimulationParams)

Perform a single timestep using the second order Adams-Bashforth method.

source