Field Interpolation
A robust field interpolation is the prerequisite for pushing particles. This example demonstrates the construction of scalar/vector field interpolators for Cartesian/Spherical grids. If the field is analytic, you can directly pass the generated function to prepare.
using TestParticle
using Meshes
using StaticArrays
using Chairmarks
function setup_spherical_field(ns = 16)
r = logrange(0.1, 10.0, length = ns)
r_uniform = range(0.1, 10.0, length = ns)
θ = range(0, π, length = ns)
ϕ = range(0, 2π, length = ns)
B₀ = 1e-8 # [nT]
B = zeros(3, length(r), length(θ), length(ϕ)) # vector
A = zeros(length(r), length(θ), length(ϕ)) # scalar
for (iθ, θ_val) in enumerate(θ)
sinθ, cosθ = sincos(θ_val)
B[1, :, iθ, :] .= B₀ * cosθ
B[2, :, iθ, :] .= -B₀ * sinθ
A[:, iθ, :] .= B₀ * sinθ
end
B_field_nu = build_interpolator(StructuredGrid, B, r, θ, ϕ)
A_field_nu = build_interpolator(StructuredGrid, A, r, θ, ϕ)
B_field = build_interpolator(StructuredGrid, B, r_uniform, θ, ϕ)
A_field = build_interpolator(StructuredGrid, A, r_uniform, θ, ϕ)
return B_field_nu, A_field_nu, B_field, A_field
end
function setup_cartesian_field(ns = 16)
x = range(-10, 10, length = ns)
y = range(-10, 10, length = ns)
z = range(-10, 10, length = ns)
B = zeros(3, length(x), length(y), length(z)) # vector
B[3, :, :, :] .= 10e-9
A = zeros(length(x), length(y), length(z)) # scalar
A[:, :, :] .= 10e-9
B_field = build_interpolator(B, x, y, z)
A_field = build_interpolator(A, x, y, z)
return B_field, A_field
end
function setup_cartesian_nonuniform_field()
x = logrange(0.1, 10.0, length = 16)
y = range(-10, 10, length = 16)
z = range(-10, 10, length = 16)
B = zeros(3, length(x), length(y), length(z)) # vector
B[3, :, :, :] .= 10e-9
A = zeros(length(x), length(y), length(z)) # scalar
A[:, :, :] .= 10e-9
B_field = build_interpolator(RectilinearGrid, B, x, y, z)
A_field = build_interpolator(RectilinearGrid, A, x, y, z)
return B_field, A_field
end
function setup_time_dependent_field(ns = 16)
x = range(-10, 10, length = ns)
y = range(-10, 10, length = ns)
z = range(-10, 10, length = ns)
# Create two time snapshots
B0 = zeros(3, length(x), length(y), length(z))
B0[3, :, :, :] .= 1.0 # Bz = 1 at t=0
B1 = zeros(3, length(x), length(y), length(z))
B1[3, :, :, :] .= 2.0 # Bz = 2 at t=1
times = [0.0, 1.0]
function loader(i)
if i == 1
# For demonstration, we assume we load from disk here
return build_interpolator(CartesianGrid, B0, x, y, z)
elseif i == 2
return build_interpolator(CartesianGrid, B1, x, y, z)
else
error("Index out of bounds")
end
end
# B_field_t(x, t)
B_field_t = LazyTimeInterpolator(times, loader)
return B_field_t
end
function setup_mixed_precision_field(ns = 11, order = 1, bc = FillExtrap(NaN); coeffs = OnTheFly())
x = range(0.0f0, 10.0f0, length = ns)
y = range(0.0f0, 10.0f0, length = ns)
z = range(0.0f0, 10.0f0, length = ns)
B = fill(0.0f0, 3, ns, ns, ns)
B[3, :, :, :] .= 1.0f-8
itp = build_interpolator(B, x, y, z, order, bc; coeffs)
return itp
end
B_sph_nu, A_sph_nu, B_sph, A_sph = setup_spherical_field();
B_car, A_car = setup_cartesian_field();
B_car_nu, A_car_nu = setup_cartesian_nonuniform_field();
B_td = setup_time_dependent_field();
itp_f32 = setup_mixed_precision_field();
loc = SA[1.0, 1.0, 1.0];
loc_f32 = SA[1.0f0, 1.0f0, 1.0f0];
loc_f64 = SA[1.0, 1.0, 1.0];3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
1.0
1.0
1.0Gridded spherical interpolation
Input Location
For spherical data, the input location is still in Cartesian coordinates!
julia> @be B_sph_nu($loc)
Benchmark: 3149 samples with 235 evaluations
min 122.698 ns (2 allocs: 64 bytes)
median 124.447 ns (2 allocs: 64 bytes)
mean 126.188 ns (2 allocs: 64 bytes)
max 243.226 ns (2 allocs: 64 bytes)
julia> @be A_sph_nu($loc)
Benchmark: 3112 samples with 293 evaluations
min 97.280 ns (2 allocs: 48 bytes)
median 99.164 ns (2 allocs: 48 bytes)
mean 101.852 ns (2 allocs: 48 bytes)
max 275.126 ns (2 allocs: 48 bytes)Uniform spherical interpolation
julia> @be B_sph($loc)
Benchmark: 3168 samples with 239 evaluations
min 120.226 ns (2 allocs: 64 bytes)
median 121.611 ns (2 allocs: 64 bytes)
mean 123.239 ns (2 allocs: 64 bytes)
max 264.682 ns (2 allocs: 64 bytes)
julia> @be A_sph($loc)
Benchmark: 3185 samples with 312 evaluations
min 91.391 ns (2 allocs: 48 bytes)
median 92.929 ns (2 allocs: 48 bytes)
mean 94.090 ns (2 allocs: 48 bytes)
max 191.577 ns (2 allocs: 48 bytes)Uniform Cartesian interpolation
julia> @be B_car($loc)
Benchmark: 3223 samples with 559 evaluations
min 49.574 ns (2 allocs: 64 bytes)
median 50.830 ns (2 allocs: 64 bytes)
mean 51.611 ns (2 allocs: 64 bytes)
max 94.039 ns (2 allocs: 64 bytes)
julia> @be A_car($loc)
Benchmark: 3215 samples with 599 evaluations
min 46.898 ns (2 allocs: 48 bytes)
median 47.970 ns (2 allocs: 48 bytes)
mean 48.623 ns (2 allocs: 48 bytes)
max 96.526 ns (2 allocs: 48 bytes)Non-uniform Cartesian interpolation
julia> @be B_car_nu($loc)
Benchmark: 2628 samples with 553 evaluations
min 50.998 ns (2 allocs: 64 bytes)
median 53.664 ns (2 allocs: 64 bytes)
mean 64.919 ns (2 allocs: 64 bytes)
max 436.588 ns (2 allocs: 64 bytes)
julia> @be A_car_nu($loc)
Benchmark: 3208 samples with 587 evaluations
min 48.097 ns (2 allocs: 48 bytes)
median 49.359 ns (2 allocs: 48 bytes)
mean 50.037 ns (2 allocs: 48 bytes)
max 92.233 ns (2 allocs: 48 bytes)Based on the benchmarks, for the same grid size, gridded interpolation (StructuredGrid with non-uniform ranges, RectilinearGrid) is 2x slower than uniform mesh interpolation (StructuredGrid with uniform ranges, CartesianGrid).
Mixed precision interpolation
Numerical field data from files is often stored in Float32. TestParticle supports constructing interpolators from Float32 data and ranges, which can then be queried with both Float32 and Float64 location vectors.
julia> @be itp_f32($loc_f32)
Benchmark: 3142 samples with 591 evaluations
min 48.161 ns (2 allocs: 64 bytes)
median 49.382 ns (2 allocs: 64 bytes)
mean 50.354 ns (2 allocs: 64 bytes)
max 93.898 ns (2 allocs: 64 bytes)
julia> @be itp_f32($loc_f64)
Benchmark: 3199 samples with 539 evaluations
min 52.455 ns (2 allocs: 64 bytes)
median 53.681 ns (2 allocs: 64 bytes)
mean 54.370 ns (2 allocs: 64 bytes)
max 103.440 ns (2 allocs: 64 bytes)Memory usage analysis
Large numerical field arrays can consume significant amounts of memory. It's important that interpolators are memory-efficient during both construction and storage.
For linear interpolation (order=1), construction is near zero-allocation because it creates a wrapper around a reinterpreted view of your existing array.
For higher-order interpolation (order = 3), FastInterpolations.jl provides two modes: OnTheFly() and PreCompute(). By default, OnTheFly() is used, which calculates interpolation coefficients during each query. This approach is memory-efficient as it does not require additional storage beyond the input data. Alternatively, PreCompute() precalculates and stores coefficients in an additional array of the same size, enabling faster
We can measure this difference using @be.
julia> # Order 1: Minimal allocations (Uses a view)
@be setup_mixed_precision_field(11, 1)
Benchmark: 658 samples with 7 evaluations
min 3.036 μs (8 allocs: 31.750 KiB)
median 3.729 μs (8 allocs: 31.750 KiB)
mean 25.866 μs (8 allocs: 31.750 KiB, 0.30% gc time)
max 10.811 ms (8 allocs: 31.750 KiB, 98.41% gc time)
julia> # Order 3: Minimal allocations (Uses a view)
@be setup_mixed_precision_field(11, 3)
Benchmark: 1524 samples with 6 evaluations
min 3.348 μs (14 allocs: 32.047 KiB)
median 4.252 μs (14 allocs: 32.047 KiB)
mean 10.691 μs (14 allocs: 32.047 KiB, 0.19% gc time)
max 4.605 ms (14 allocs: 32.047 KiB, 99.08% gc time)Comparing the ratios relative to the original array size illustrates the overhead:
julia> B = fill(0.0f0, 3, 11, 11, 11);
julia> size_B = Base.summarysize(B) # Original 4D field array size
16036
julia> size_itp1 = Base.summarysize(itp_f32) # Total size for order=1 (Essentially the input array size)
16172
julia> itp_f32_q = setup_mixed_precision_field(11, 3) # Total size for order=3
(::TestParticle.FieldInterpolator{FastInterpolations.HeteroInterpolantND{Float32, StaticArraysCore.SVector{3, Float32}, 3, Tuple{FastInterpolations._CachedRange{Float32, Float32}, FastInterpolations._CachedRange{Float32, Float32}, FastInterpolations._CachedRange{Float32, Float32}}, Tuple{FastInterpolations.CardinalInterp{Float64, FastInterpolations.NoBC}, FastInterpolations.CardinalInterp{Float64, FastInterpolations.NoBC}, FastInterpolations.CardinalInterp{Float64, FastInterpolations.NoBC}}, Tuple{FillExtrap{StaticArraysCore.SVector{3, Float32}}, FillExtrap{StaticArraysCore.SVector{3, Float32}}, FillExtrap{StaticArraysCore.SVector{3, Float32}}}, Tuple{FastInterpolations.AutoSearch, FastInterpolations.AutoSearch, FastInterpolations.AutoSearch}, Array{StaticArraysCore.SVector{3, Float32}, 3}}}) (generic function with 2 methods)
julia> size_itp3 = Base.summarysize(itp_f32_q)
16196
julia> # Ratios relative to raw data
size_itp1 / size_B
1.008480917934647
julia> size_itp3 / size_B
1.0099775505113495As a rule of thumb, linear interpolation and cubic interpolation with OnTheFly() coefficients have nearly zero memory overhead (ratio ≈ 1.0), as they both operate directly on the input data. When supported, cubic interpolation with PreCompute() coefficients increases the memory footprint by
On-the-fly vs Precomputed coefficients
Cubic interpolation (order = 3) requires high-order coefficients. By default, TestParticle uses OnTheFly() coefficients, which are calculated at query time. This saves memory but increases evaluation time. For maximum performance, you can use PreCompute(), which stores the coefficients in an additional array.
Status in v0.4.8
Support for PreCompute() with local Hermite cubic interpolation is currently under development and not yet available in FastInterpolations.jl v0.4.8.
julia> # Benchmark evaluation time
itp_fly = setup_mixed_precision_field(11, 3; coeffs = OnTheFly());
julia> # itp_pre = setup_mixed_precision_field(11, 3; coeffs = PreCompute()); # Not yet supported in v0.4.8
@be itp_fly($loc_f64)
Benchmark: 3134 samples with 55 evaluations
min 527.345 ns (2 allocs: 64 bytes)
median 532.818 ns (2 allocs: 64 bytes)
mean 539.466 ns (2 allocs: 64 bytes)
max 1.160 μs (2 allocs: 64 bytes)
julia> # Compare total object size
Base.summarysize(itp_fly)
16196As shown, OnTheFly() preserves memory efficiency while providing higher-order accuracy. Once available, PreCompute() will offer a faster alternative for memory-abundant systems.
Time-dependent field interpolation
For time-dependent fields, we can use LazyTimeInterpolator. It takes a list of time points and a loader function that returns a spatial interpolator for a given time index. The interpolator will linearly interpolate between the two nearest time points.
julia> @be B_td($loc, 0.5)
Benchmark: 3130 samples with 177 evaluations
min 164.090 ns (2 allocs: 64 bytes)
median 166.469 ns (2 allocs: 64 bytes)
mean 168.694 ns (2 allocs: 64 bytes)
max 302.825 ns (2 allocs: 64 bytes)Related API
TestParticle.build_interpolator Function
build_interpolator(gridtype, A, grids..., order::Int=1, bc=FillExtrap(NaN))
build_interpolator(A, grids..., order::Int=1, bc=FillExtrap(NaN))Return a function for interpolating field array A on the given grids.
Arguments
gridtype:CartesianGrid,RectilinearGridorStructuredGrid. Usually determined by the number of grids.A: field array. For vector field, the first dimension should be 3 if it's not an SVector wrapper.order::Int=1: order of interpolation in [0,1,3].bc=FillExtrap(NaN): boundary condition type fromFastInterpolations.jl.FillExtrap(NaN): Fill with NaN (default).ClampExtrap(): Clamp (flat extrapolation).WrapExtrap(): Exclusive periodic wrapping ().
coeffs=OnTheFly(): coefficient strategy for cubic interpolation (order=3). Default isOnTheFly().
Notes
- The input array
Amay be modified in-place for memory optimization.
TestParticle.prepare Function
prepare(args...; kwargs...) -> (q2m, m, E, B, F)
prepare(E, B, F = ZeroField(); kwargs...)
prepare(grid::CartesianGrid, E, B, F = ZeroField(); kwargs...)
prepare(x, E, B, F = ZeroField(); dir = 1, kwargs...)
prepare(x, y, E, B, F = ZeroField(); kwargs...)
prepare(x, y, z, E, B, F = ZeroField(); kwargs...)
prepare(B; E = ZeroField(), F = ZeroField(), kwargs...)Return a tuple consists of particle charge-mass ratio for a prescribed species of charge q and mass m, mass m for a prescribed species, analytic/interpolated EM field functions, and external force F.
Prescribed species are Electron and Proton; other species can be manually specified with m and q keywords or species = Ion(m̄, q̄), where m̄ and q̄ are the mass and charge numbers respectively.
Direct range input for uniform grid in 1/2/3D is supported. The grid vectors must be sorted. For 1D grid, an additional keyword dir is used for specifying the spatial direction, 1 -> x, 2 -> y, 3 -> z. For 3D grid, the default grid type is CartesianGrid. To use StructuredGrid (spherical) grid, an additional keyword gridtype is needed. For StructuredGrid (spherical) grid, dimensions of field arrays should be (Br, Bθ, Bϕ).
Keywords
order::Int=1: order of interpolation in [0,1,3].bc=FillExtrap(NaN): boundary condition type fromFastInterpolations.jl.species=Proton: particle species.q=nothing: particle charge.m=nothing: particle mass.gridtype:CartesianGrid,RectilinearGrid,StructuredGrid.