Library
Public API
SuperVFM.GetSchwarzTempCoeffs
— MethodGetSchwarzTempCoeffs(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.
SuperVFM.check_timestep
— Methodcheck_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))}\]
SuperVFM.computeEnergySpectrum
— MethodcomputeEnergySpectrum(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.
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.\]
SuperVFM.computeEnergySpectrum
— MethodcomputeEnergySpectrum(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.\]
SuperVFM.generateMapping
— MethodgenerateMapping(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\]
SuperVFM.get_density
— Methodget_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$
SuperVFM.get_density
— Methodget_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$
SuperVFM.get_dynamic_viscosity
— Methodget_dynamic_viscosity(T::AbstractFloat)
Returns the dynamic viscosity $η$ at temperature $T$.
Requires $0.8\leq T < T_{\lambda} = 2.178$
SuperVFM.get_dynamic_viscosity
— Methodget_dynamic_viscosity(T::Unitful.Temperature)
Returns the dynamic viscosity $η$ at temperature $T$.
Requires $0.8\leq T \leq T_{\lambda} = 2.178$
SuperVFM.print_characteristics
— Methodprint_characteristics(io::IO,SimParams::SimulationParams)
Prints the characteristic time and length scales of the simulation in dimensional units
SuperVFM.print_density_data
— Functionprint_density_data(io::IO=stdout)
Prints the normal fluid $ρ_n$ and superfluid $ρ_s$densities to the IO buffer. Defaults to stdout
.
SuperVFM.print_dynamic_visosity_data
— Functionprint_dynamic_visosity_data(io::IO=stdout)
Prints the dynamic viscosity η to the IO buffer. Defaults to stdout
.
KernelAbstractions.CPU
— TypeCPU(; 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.
SuperVFM.OpenBoundary
— TypeOpenBoundary(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)
SuperVFM.PeriodicBoundary
— TypePeriodicBoundary(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)
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.
Internal API
Base.show
— MethodBase.show(io::IO,SimParams::SimulationParams)
Lists the simulation parameters stored in SimParams
in a stylistic way with a key.
KernelAbstractions.zeros
— MethodKA.zeros(BE::Backend,::Type{SVector{S,T}},N::Int) where {S,T}
Initialise an array of static vectors of size S
and type T
SuperVFM.SchwarzModel_kernel!
— MethodSchwarzModel_kernel!(u_mf, u_sup, f, f_infront, ghosti, ghostb, FM::SchwarzModel, normal_velocity)
Kernel launch for the Schwarz model.
SuperVFM.all_periodic_enforce_kernel!
— Methodall_periodic_enforce_kernel!(f, f_infront)
Launch kernel to enforce periodic boundary conditions in all dimensions.
SuperVFM.check_active_pcount
— Methodcheck_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).
SuperVFM.check_output_folder_structure
— Methodcheck_output_folder_structure(io::IO)
Creates the output folder structure if it did not exist already. Prints result to buffer.
SuperVFM.compute_LIA_kernel!
— Methodcompute_LIA_kernel!(u_loc, f, ghosti, ghostb, κ, corea)
Kernel to compute the superfluid velocity.
SuperVFM.compute_filament_velocity!
— Methodcompute_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)\]
SuperVFM.compute_filament_velocity!
— Methodcompute_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\]
SuperVFM.compute_velocity!
— Methodcompute_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 \]
SuperVFM.create_info_file
— Methodcreate_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.
SuperVFM.enforce_boundary!
— Methodenforce_boundary!(f, boundary_x::OpenBoundary, boundary_y::OpenBoundary, boundary_z::OpenBoundary; kwargs...)
Enforces open boundary conditions across all 3 dimensions.
SuperVFM.enforce_boundary!
— Methodenforce_boundary!(f, boundary_x::PeriodicBoundary, boundary_y::PeriodicBoundary, boundary_z::PeriodicBoundary; kwargs...)
Enforces periodic boundary conditions across all 3 dimensions.
SuperVFM.generate_initial_files
— Methodgenerate_initial_files(SP::SimulationParams)
Generates the output folder structure and the initialisation files.
SuperVFM.get_curvature!
— Methodget_curvature!(curv; kwargs...)
Calculate the vortex line curvature $\zeta$.
\[ \zeta = |\mathbf{s}''| = \left|\frac{d\mathbf{s}}{d\xi}\right|\]
SuperVFM.get_curvature_kernel!
— Methodget_curvature_kernel!(curv, f, f_infront, ghosti, ghostb)
Kernel launch for computing vortex line curvature.
SuperVFM.get_deriv_1
— Methodget_deriv_1(f::AbstractArray, ghosti::AbstractArray, ghostb::AbstractArray)
Computes the first order derivative using a 2nd order finite difference adaptive scheme using ghost points.
SuperVFM.get_deriv_2
— Methodget_deriv_2(f::AbstractArray, ghosti::AbstractArray, ghostb::AbstractArray)
Computes the first order derivative using a 2nd order finite difference adaptive scheme using ghost points.
SuperVFM.get_Δξ
— Methodget_Δξ(f, ghosti, f_infront, pcount, SP::SimulationParams{S,T}) where {S,T}
Compute the seperation distance $\Delta\xi$ between each filament.
SuperVFM.get_Δξ_kernel!
— Methodget_Δξ_kernel!(Δξ, f, ghosti, f_infront)
Launch kernel to compute $\Delta\xi$.
SuperVFM.ghostp!
— Method ghostp!(ghosti, ghostb; kwargs...)
In place variant of ghostp
.
SuperVFM.ghostp
— Methodghostp(f, f_infront, f_behind, pcount, SP::SimulationParams{S,T}) where {S,T}
Compute ghost points.
SuperVFM.ghostp_kernel!
— Methodghostp_kernel!(ghosti, ghostb, f, fint, box_size)
Kernel for computation of ghost points.
SuperVFM.identify_loop!
— MethodIdentify 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.
SuperVFM.identify_loops
— MethodIdentify reference particles for different vortex loops.
SuperVFM.initVortex_kernel!
— MethodinitVortex_kernel!(f, f_infront, f_behind, pcount, initf::SimpleTrefoil)
Launch kernel to initialise a simple trefoil.
SuperVFM.initialiseVortex
— MethodinitialiseVortex(SP::SimulationParams{S,T}) where {S,T}
Uses SP.backend
to initialise the vortex structure according to the initial condition SP.initf
.
SuperVFM.premove!
— Methodpremove()
SuperVFM.save_vortex
— Methodsave_vortex(it;kwargs...)
Save the vortex and related information to file.
SuperVFM.static2field
— Methodstatic2field(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
SuperVFM.timestep!
— Methodtimestep!(f, u, u1, u2, f_infront, pcount, SP::SimulationParams)
Perform a single timestep using the second order Adams-Bashforth method.
SuperVFM.timestep_kernel!
— Methodtimestep_kernel!(f, u, u1, u2, f_infront, dt)
Kernel launch for timestep.