API

Index

Docstrings

UltraDark.GridsType
struct Grids <: UltraDark.AbstractGrids
Grids(length, resol::Int)
Grids(length_tuple, resol_tuple::Tuple{Int, Int, Int})

struct containing grids used in a simulation

Fields

  • x::Array{Float64, 3}: Array of x positions

  • y::Array{Float64, 3}: Array of y positions

  • z::Array{Float64, 3}: Array of z positions

  • kx::Array{Float64, 3}: Array of x Fourier modes

  • ky::Array{Float64, 3}: Array of y Fourier modes

  • kz::Array{Float64, 3}: Array of z Fourier modes

  • k::Array{Float64, 3}: Fourier space postition array

  • rkx::Array{Float64, 3}: Array of x Fourier modes for use with rfft

  • rky::Array{Float64, 3}: Array of y Fourier modes for use with rfft

  • rkz::Array{Float64, 3}: Array of z Fourier modes for use with rfft

  • rk::Array{Float64, 3}: Fourier space postition array for use with rfft

  • ψx::Array{ComplexF64, 3}: ψ field

  • ψk::Array{ComplexF64, 3}: ψ field in Fourier space

  • ρx::Array{Float64, 3}: density field ρ

  • ρk::Array{ComplexF64, 3}: density field ρ in Fourier space

  • Φx::Array{Float64, 3}: gravitational potential field Φ

  • Φk::Array{ComplexF64, 3}: gravitational potential field Φ in fourier space

  • fft_plan::Any: FFT plan for complex-to-complex transforms

  • rfft_plan::Any: FFT plan for real-to-complex transforms

  • k_vanish_indices::Vector{CartesianIndex{3}}: Indices at which k==0.0

  • rk_vanish_indices::Any: Indices at which rk==0.0

Examples

Create an empty grid with length length and resolution resol

julia> using UltraDark

julia> length = 1;

julia> resol = 16;

julia> Grids(length, resol);

julia> size(ans.ψx)
(16, 16, 16)

Create an empty length[1]xlength[2]xlength[3] grid with resolution resol[1]xresol[2]xresol[3].

julia> using UltraDark

julia> Grids((1.0, 1.0, 0.5), (64, 64, 32));

julia> size(ans.ψx)
(64, 64, 32)
source
UltraDark.PencilGridsType
PencilGrids(length, resol)
PencilGrids(length_tuple, resol_tuple::Tuple{Int, Int, Int})

struct containing grids used in a simulation

Each grid is a PencilArray, allowing multiprocess FFTs. This comes with significant overhead so is only useful when running in a multi-node environment.

Fields

  • x::Array{Float64, 3}: Array of x positions

  • y::Array{Float64, 3}: Array of y positions

  • z::Array{Float64, 3}: Array of z positions

  • kx::Array{Float64, 3}: Array of x Fourier modes

  • ky::Array{Float64, 3}: Array of y Fourier modes

  • kz::Array{Float64, 3}: Array of z Fourier modes

  • k::Any: Fourier space postition array

  • rkx::Array{Float64, 3}: Array of x Fourier modes for use with rfft

  • rky::Array{Float64, 3}: Array of y Fourier modes for use with rfft

  • rkz::Array{Float64, 3}: Array of z Fourier modes for use with rfft

  • rk::Any: Fourier space postition array for use with rfft

  • ψx::Any: ψ field

  • ψk::Any: ψ field in Fourier space

  • ρx::Any: density field ρ

  • ρk::Any: density field ρ in Fourier space

  • Φx::Any: gravitational potential field Φ

  • Φk::Any: gravitational potential field Φ in fourier space

  • fft_plan::Any: FFT plan for complex-to-complex transforms

  • rfft_plan::Any: FFT plan for real-to-complex transforms

  • k_vanish_indices::Vector{CartesianIndex{3}}: Indices at which k==0.0

  • rk_vanish_indices::Any: Indices at which rk==0.0

  • MPI_COMM::Any: MPI communicator

Examples

Create an empty grid with length length and resolution resol. Uses PencilFFTs to create PencilArrays.

julia> using UltraDark

julia> len = 1;

julia> resol = 16;

julia> PencilGrids(len, resol);

Create an empty length[1]xlength[2]xlength[3] grid with resolution resol[1]xresol[2]xresol[3].

julia> using UltraDark

julia> PencilGrids((1.0, 1.0, 0.5), (64, 64, 32));
source
UltraDark.E_gravity_densityMethod
E_gravity_density(psi::Number, Phi::Real) -> Real

Gravitational energy density of field psi in gravitational potential Phi

source
UltraDark.E_totalMethod
E_total(grids; constants)

Total energy of the scalar field: the sum of the kinetic and quantum energies.

source
UltraDark.actual_time_stepMethod
actual_time_step(max_time_step, time_interval, time_step_options)

Actual size and number of time steps that should be taken if the maximum is max_time_step. No more than time_step_options.update_period steps should be taken, and they should fit in time_interval.

Examples

julia> using UltraDark: actual_time_step, TimeStepOptions

julia> actual_time_step(0.11, 1, TimeStepOptions())
(0.1, 10)
source
UltraDark.add_external_potential!Method
add_external_potential!(t, grids, constants)

Add a gravitational potential to the grid. By default this does nothing, but can be overridden in multiple dispatch.

source
UltraDark.angular_momentumMethod
angular_momentum(grids)
angular_momentum(grids, ψx, ρx)

Calculate total angular momentum

Returns

L: AbstractArray with length 3

source
UltraDark.angular_momentum_densityMethod
angular_momentum_density(grids)
angular_momentum_density(grids, ψx, ρx)

Calculate angular momentum density at each grid point

Returns

L: AbstractArray with dimensions 3 x resolx x resoly x resol_z

source
UltraDark.auxiliary_step!Method
auxiliary_step!(Δt, grids, t, constants)
auxiliary_step!(Δt, grids, t, constants, s; a = 1.0)

Do an auxiliary inner step. By default this does nothing, but can be overridden in multiple dispatch.

source
UltraDark.azimuthal_angleMethod
azimuthal_angle(x, y, z)
azimuthal_angle(grids)
azimuthal_angle(grids, r0)

Calculate the azimuthal angle in spherical or cylindrical coordinates

This is \phi in conventional physics notation.

source
UltraDark.dVMethod
dV(grids) -> Any

Calculate the volume of each grid cell

Examples

julia> using UltraDark

julia> box_length = 1.0;

julia> resol = 16;

julia> g = Grids(box_length, resol);

julia> dV(g) * resol^3 == box_length^3
true
source
UltraDark.evolve_to!Method
evolve_to!(
    t_start,
    t_end,
    grids,
    output_config,
    sim_config;
    constants,
    external_states
) -> Any

Evolve grids forward from t_start to t_end

source
UltraDark.inner_step!Method
inner_step!(Δt, grids, constants; a=1.0)
inner_step!(Δt, grids, constants, s; a=1.0)

Perform the "inner" time step in the symmetrized split-step Fourier method.

This step applies the diffusion terms and updates the gravitational potential.

source
UltraDark.k_vecMethod
k_vec(
    lengths,
    resols
) -> Tuple{AbstractFFTs.Frequencies, AbstractFFTs.Frequencies, AbstractFFTs.Frequencies}

Calculate the Fourier frequencies of a box with side lengths lengths and resolutions resols

Examples

julia> using UltraDark: k_vec

julia> kvec = k_vec((2π, 2π, 2π), (4, 4, 4));

julia> kvec[1]
4-element AbstractFFTs.Frequencies{Float64}:
  0.0
  1.0
 -2.0
 -1.0
source
UltraDark.massMethod
mass(grids)
mass(grids, rho)

Calculate total mass of a density field

Examples

julia> using UltraDark

julia> g = Grids(1.0, 16);

julia> g.ρx .= 0.0;

julia> g.ρx[1, 1, 1] = 1.0;

julia> UltraDark.mass(g) == 1.0 * (1.0 / 16)^3
true
source
UltraDark.max_normed_phase_diffMethod
normed_max_phase_grad(grids)

Compute maximum phase gradient of a grid

Normalised to ignore large gradients in regions with low density. These tend to be anomalous.

source
UltraDark.max_time_stepMethod
max_time_step(grids, a, external_states)

Consolidate the maximum time step implied by grids and each member of external_states.

source
UltraDark.max_time_step_gridsMethod
max_time_step_grids(grids, a)
max_time_step_grids(grids::PencilGrids, a)

Calculate an upper bound on the time step from grid properties

This time step depends on the gravitational potential and the resolution.

source
UltraDark.outer_step!Method
outer_step!(Δt, grids, constants; a=1.0)
outer_step!(Δt, grids, constants, s; a=1.0)

Perform the "outer" time step in the symmetrized split-step Fourier method.

This step only updates the phase of ψ applying accelerations due to gravity, the amplitude is not changed.

source
UltraDark.phase_diffMethod
phase_diff(field, dir)

Compute point-to-point difference of phase on a grid along a direction

Returns an array of size (size(field)[1], size(field)[2], size(field)[3]) containing differences in direction dir.

source
UltraDark.phase_diffMethod
phase_diff(field::PencilArray, dir)

Compute point-to-point difference of phase on a grid along a direction

Returns an array of size (size(field)[1], size(field)[2], size(field)[3]) containing differences in direction dir.

source
UltraDark.polar_angleMethod
polar_angle(x, y, z)
polar_angle(grids)
polar_angle(grids, r0)

Calculate the polar angle in spherical coordinates

This is \theta in conventional physics notation.

source
UltraDark.radius_cylindricalMethod
radius_cylindrical(x, y, z)
radius_cylindrical(grids)
radius_cylindrical(grids, r0)

Calculate the radial coordinate in cylindrical coordinates

Examples

julia> using UltraDark

julia> import UltraDark: radius_cylindrical, azimuthal_angle

julia> box_length = 1.0;

julia> resol = 16;

julia> g = Grids(box_length, resol);

julia> all(radius_cylindrical(g) .* cos.(azimuthal_angle(g)) .≈ g.x)
true

julia> all(
           radius_cylindrical(g, (0.0, 0.0, 1.0)) .* cos.(azimuthal_angle(g, (0.0, 0.0, 1.0))) .≈
           g.x,
       )
true
source
UltraDark.radius_sphericalMethod
radius_spherical(x, y, z)
radius_spherical(grids)
radius_spherical(grids, r0)

Calculate the radial coordinate in a spherical coordinate system

Examples

julia> using UltraDark

julia> import UltraDark: radius_spherical, polar_angle, azimuthal_angle

julia> box_length = 1.0;

julia> resol = 16;

julia> g = Grids(box_length, resol);

julia> all(radius_spherical(g) .* sin.(polar_angle(g)) .* cos.(azimuthal_angle(g)) .≈ g.x)
true

julia> all(
           radius_spherical(g, (1.0, 0.0, 0.0)) .* sin.(polar_angle(g, (1.0, 0.0, 0.0))) .*
           cos.(azimuthal_angle(g, (1.0, 0.0, 0.0))) .+ 1.0 .≈ g.x,
       )
true
source
UltraDark.rk_vecMethod
rk_vec(
    lengths,
    resols
) -> Tuple{AbstractFFTs.Frequencies, AbstractFFTs.Frequencies, AbstractFFTs.Frequencies}
source
UltraDark.simulate!Method
simulate!(grids, output_config::OutputConfig; constants = nothing, external_states = (),)
simulate!(grids, sim_config, output_config::OutputConfig; constants = nothing, external_states = (),)

simulate! is the main entry point into a simulation

simulate! evolves grids and external_states forward, writing output as specified by output_config.

Arguments

grids::AbstractGrids Contains coordinate systems and the scalar field ψ.

sim_config::UltraDark.Config gives further control over the way the simulation is done.

output_config::OutputConfig controls what output is written and when.

constants::Any can be used to overload other functions in the evolution, for example to introduce self-interactions.

external_states::Tuple{Any} can be used in combination with overloading to add other dynamical objects to a simulation.

source
UltraDark.take_steps!Method

Take n steps with time step Δt

Examples

julia> using UltraDark: take_steps!, Grids, OutputConfig, Config

julia> take_steps!(
           Grids(1.0, 16),
           0,
           0.5,
           10,
           OutputConfig(mktempdir(), []),
           Config.constant_scale_factor,
           nothing,
           (),
       )
5.0
source
UltraDark.update_gravitational_potential!Method
update_gravitational_potential!(grids; a = 1.0)
update_gravitational_potential!(grids, constants; a = 1.0)

Update density grids.ρx and the gravitational potential grids.Φx based on grids.ψx

source