Skip to content

API Reference

Types

Magnetostatics.AbstractMagneticField Type
julia
AbstractMagneticField

Abstract supertype for all magnetic field representations.

source
Magnetostatics.AbstractCurrentSource Type
julia
AbstractCurrentSource

Abstract supertype for all current sources.

source
Magnetostatics.AbstractSolver Type
julia
AbstractSolver

Abstract supertype for all magnetic field solvers.

source
Magnetostatics.Wire Type
julia
Wire{T} <: AbstractCurrentSource

A unified structure for a current carrying wire. It can represent a straight wire or a loop depending on the geometry of points.

Fields

  • points::Vector{SVector{3, T}}: The points defining the wire geometry.

  • current::T: The current flowing through the wire in Amperes.

source
Magnetostatics.CurrentLoop Type
julia
CurrentLoop{T}

A circular current loop.

Fields

  • radius::T: Radius of the loop [m].

  • current::T: Current in the loop [A].

  • center::SVector{3, T}: Position of the loop center [m].

  • normal::SVector{3, T}: Unit normal vector of the loop.

source
Magnetostatics.HarrisSheet Type
julia
HarrisSheet{T} <: AbstractMagneticField

Magnetic field of a Harris current sheet: B(z) = B0 * tanh(z/L) * x_hat.

Fields

  • B0::T: Asymptotic magnetic field strength.

  • L::T: Half-width of the current sheet.

source
Magnetostatics.Dipole Type
julia
Dipole{T} <: AbstractMagneticField

Magnetic field of a dipole moment M at the origin.

Fields

  • M::SVector{3, T}: Magnetic dipole moment.
source
Magnetostatics.InfiniteWire Type
julia
InfiniteWire{T} <: AbstractMagneticField

Magnetic field of an infinite straight wire.

Fields

  • I::T: Current in the wire [A].

  • r::T: Radius of the wire [m].

  • center::SVector{3, T}: A point on the wire [m].

  • direction::SVector{3, T}: Unit vector pointing in the direction of the current.

source
Magnetostatics.CurrentLoopAnalytic Type
julia
CurrentLoopAnalytic{T} <: AbstractMagneticField

Analytical magnetic field of a circular current loop.

source

Boundary Conditions

Magnetostatics.AbstractBoundary Type
julia
AbstractBoundary

Abstract supertype for all magnetic field boundary conditions.

source
Magnetostatics.OpenBoundary Type
julia
OpenBoundary <: AbstractBoundary

Free-space boundary condition: the field decays naturally to zero at infinity. This is the default for point-source solvers such as BiotSavart.

source
Magnetostatics.PeriodicBoundary Type
julia
PeriodicBoundary <: AbstractBoundary

Periodic boundary condition: sources are replicated at integer multiples of the given period in each Cartesian direction before field summation.

Fields

  • period::SVector{3,Float64}: box side lengths (Lx, Ly, Lz) [m].
source
Magnetostatics.ConductingWall Type
julia
ConductingWall <: AbstractBoundary

Planar perfectly-conducting wall boundary condition, implemented via the method of images. The wall is an infinite plane whose normal points along one of the Cartesian axes.

Fields

  • axis::Int: axis perpendicular to the wall (1 = x, 2 = y, 3 = z).

  • position::Float64: coordinate of the wall along that axis [m].

source
Magnetostatics.ConductingSphere Type
julia
ConductingSphere <: AbstractBoundary

Perfectly-conducting spherical boundary condition solved numerically via a boundary-element method (BEM). The sphere surface is discretized into patches; a linear system is solved for the induced surface current K that cancels the normal component of the incident field at every patch centre.

Fields

  • center::SVector{3,Float64}: centre of the sphere [m].

  • radius::Float64: radius of the sphere [m].

  • n_theta::Int: number of patches in the polar (latitude) direction.

  • n_phi::Int: number of patches in the azimuthal (longitude) direction.

source
Magnetostatics.PrescribedBoundary Type
julia
PrescribedBoundary{F} <: AbstractBoundary

Dirichlet-type background boundary condition: a prescribed magnetic field field_func(r) is added to the source-computed field everywhere in space.

Fields

  • field_func::F: callable with signature (r::SVector{3}) -> SVector{3}.
source

Solvers

Magnetostatics.BiotSavart Type
julia
BiotSavart{BC <: AbstractBoundary} <: AbstractSolver

Solver that uses the Biot-Savart law to compute the magnetic field. The boundary condition BC controls how sources outside the physical domain are treated (e.g., image sources for a conducting wall).

The default constructor BiotSavart() uses OpenBoundary.

source
Magnetostatics.FFTSolver Type
julia
FFTSolver{BC <: AbstractBoundary} <: AbstractSolver

Solver that uses Fast Fourier Transforms to compute the magnetic field B from a current distribution J.

The default constructor FFTSolver() uses PeriodicBoundary, which is inherent to the FFT method. Use FFTSolver(OpenBoundary()) to enable zero-padding for non-periodic (free-space) fields.

source
Magnetostatics.VectorPotential Type
julia
VectorPotential{BC <: AbstractBoundary} <: AbstractSolver

Solver that computes the magnetic vector potential A. The default constructor VectorPotential() uses OpenBoundary.

source
Magnetostatics.PoissonSolver Type
julia
PoissonSolver{ALG} <: AbstractSolver

Grid-based solver for the magnetostatic vector potential equation ∇²A = −μ₀ J on a uniform Cartesian grid with homogeneous Dirichlet boundary conditions (A = 0 on all domain faces).

Each Cartesian component is solved independently with the same 7-point finite-difference Laplacian matrix via LinearSolve.jl. Pass anyLinearSolve algorithm as alg; the default is KrylovJL_CG().

Fields

  • xs, ys, zs: grid node positions as AbstractRange objects.

  • alg: LinearSolve algorithm (e.g. KrylovJL_CG(), KLUFactorization(), UMFPACKFactorization()).

Example

julia
xs = range(-1.0, 1.0; length = 32)
ys = range(-1.0, 1.0; length = 32)
zs = range(-0.5, 0.5; length = 16)
solver = PoissonSolver(xs, ys, zs) # default: CG
solver_direct = PoissonSolver(xs, ys, zs; alg = KLUFactorization())
A = solve(solver, J_grid)  # J_grid :: (3,nx,ny,nz)
source
Magnetostatics.solve Function
julia
solve(solver::BiotSavart, source::Wire, r) where T

Compute the magnetic field at position r generated by a Wire source using the Biot-Savart law. It assumes the wire is composed of straight segments between the points.

The magnetic field for a finite wire segment is calculated using the algebraic form:

B=μ0I4πdl×a|dl×a|2(dla|a|dlb|b|)

where:

a=rrstart

is the vector from the start of the segment to the observation point.

  • b=rrend

    is the vector from the end of the segment to the observation point.

  • dl=rendrstart

    is the segment vector.

source
julia
solve(solver::BiotSavart, mesh::SurfaceCurrentMesh, r)

Compute the magnetic field at r from a SurfaceCurrentMesh using the far-field Biot-Savart approximation for each patch:

julia
dB(r) = (μ₀/) * Aᵢ * (Kᵢ × Rᵢ) / |Rᵢ|³
source
julia
solve(solver::FFTSolver, J::AbstractArray{T, 4}, dx::Real) where T

Compute the magnetic field B from a discrete current distribution J using FFT.

Arguments

  • J: 4D array of size (3, Nx, Ny, Nz) representing the current density components (Jx, Jy, Jz).

  • dx: Grid spacing (assumed isotropic for now).

Returns

  • B: 4D array of size (3, Nx, Ny, Nz) representing the magnetic field components (Bx, By, Bz).
source
julia
solve(solver::FFTSolver{OpenBoundary}, J, dx)

Open-boundary (free-space) variant: zero-pads J to twice the size in each spatial dimension before the FFT to suppress the periodic wrap-around artefact inherent to the standard FFT kernel. The result is cropped back to the original domain size.

source
julia
solve(solver::VectorPotential, source::Wire, r) where T

Compute the vector potential A at position r generated by a Wire source. Integral of (mu0 I / 4pi) * dl / |r - r'|.

source
julia
solve(solver::VectorPotential, source::Dipole, r) where T

Compute vector potential A for a magnetic dipole. A = (mu0 / 4pi) * (m x r) / r^3

source
julia
solve(solver::VectorPotential, source::CurrentLoop, r::AbstractVector)

Compute vector potential A for a circular current loop.

source
julia
solve(solver::PoissonSolver, J::AbstractArray{T,4}) where T

Solve ∇²A = −μ₀ J on the grid defined by solver. Grid spacings are derived automatically from solver.xs, solver.ys, solver.zs, supporting anisotropic grids (dx ≠ dy ≠ dz).

Arguments

  • J: 4-D array of size (3, nx, ny, nz) (same convention as FFTSolver).

Returns

  • A: 4-D array of size (3, nx, ny, nz) with the vector potential.
source

Surface Current (BEM)

Magnetostatics.SurfaceCurrentMesh Type
julia
SurfaceCurrentMesh{T} <: AbstractCurrentSource

Discretized surface-current distribution on a closed surface. Each patch i carries a surface current density K [A/m] expressed in the local tangent frame (tangents1[i], tangents2[i]).

Created by compute_surface_current for a ConductingSphere boundary.

source
Magnetostatics.compute_surface_current Function
julia
compute_surface_current(bc::ConductingSphere, source::Wire) -> SurfaceCurrentMesh

Solve the BEM linear system for the induced surface current K on the perfectly-conducting sphere defined by bc that cancels the normal component of the incident field from source.

The incident field is computed with the free-space BiotSavart kernel.

source

Utilities

Magnetostatics.discretize_loop Function
julia
discretize_loop(radius, n_segments, current; center=[0,0,0], normal=[0,0,1])

Discretize a circular current loop into a Wire object.

source
julia
discretize_loop(loop::CurrentLoop, n_segments)

Discretize a CurrentLoop object into a Wire object.

source
Magnetostatics.getB_loop Function
julia
getB_loop(r, loop::CurrentLoop)

Calculate the magnetic field B [T] at 3D point r from a CurrentLoop.

source
Magnetostatics.set_current_wire! Function
julia
set_current_wire!(J, x, y, z, point, direction, current, width)

Set the current density contribution of a straight wire with a Gaussian profile to J.

Arguments

  • J: 4D array of size (3, Nx, Ny, Nz) representing the current density components.

  • x, y, z: 1D arrays identifying the grid coordinates.

  • point: A point on the wire (Cartesian coordinates).

  • direction: Unit vector of the wire direction.

  • current: Total current in [A].

  • width: Gaussian width parameter.

source
Magnetostatics.set_current_wire Function
julia
set_current_wire(x, y, z, point, direction, current, width)

Create a new current density array J for a straight wire with a Gaussian profile. See set_current_wire!.

source
Magnetostatics.compute_curl! Function
julia
compute_curl!(B, A, xs, ys, zs)

Compute B = ∇ × A in-place using second-order central differences on a uniform Cartesian grid with potentially anisotropic spacing.

Interior points are updated; boundary layers are left unchanged (zero).

Arguments

  • B: output 4-D array of size (3, Nx, Ny, Nz).

  • A: input 4-D array of size (3, Nx, Ny, Nz).

  • xs, ys, zs: grid coordinate ranges (spacing extracted via step).

source
Magnetostatics.compute_curl Function
julia
compute_curl(A, xs, ys, zs)

Allocate and return B = ∇ × A. See compute_curl! for details.

source
Magnetostatics.dipole_fieldline Function
julia
dipole_fieldline(ϕ, L=2.5, nP=100)

Creates nP points on one field line of the magnetic field from a dipole. In a centered dipole magnetic field model, the path along a given L shell can be described as r = L*cos²λ, where r is the radial distance (in planetary radii) to a point on the line, λ is its co-latitude, and L is the L-shell of interest.

source
Magnetostatics.image_source Function
julia
image_source(source::Wire, bc::ConductingWall) -> Wire

Return the image Wire of source reflected across the planar conducting wall defined by bc. The image current is negated so that the tangential magnetic field vanishes at the wall surface.

source

Configurations

Magnetostatics.getB_mirror Function
julia
getB_mirror(x, y, z, distance, a, I1) :: SVector{3}
getB_mirror(r, distance, a, I1) :: SVector{3}

Get magnetic field at [x, y, z] from a magnetic mirror generated from two coils.

Arguments

  • r: location, vector of length 3 [m]

  • x,y,z: location [m]

  • distance: distance between solenoids in [m].

  • a: radius of each side coil in [m].

  • I1: current in the solenoid times number of windings in side coils.

source
Magnetostatics.getB_bottle Function
julia
getB_bottle(x, y, z, distance, a, b, I1, I2) :: SVector{3}
getB_bottle(r, distance, a, b, I1, I2) :: SVector{3}

Get magnetic field from a magnetic bottle. Reference: wiki

Arguments

  • r: location, vector of length 3 [m]

  • x,y,z: location [m]

  • distance::Float: distance between solenoids in [m].

  • a::Float: radius of each side coil in [m].

  • b::Float: radius of central coil in [m].

  • I1::Float: current in the solenoid times number of windings in side coils in [A].

  • I2::Float: current in the central solenoid times number of windings in the central loop in [A].

source
Magnetostatics.getB_tokamak_coil Function
julia
getB_tokamak_coil(x, y, z, a, b, ICoils, IPlasma) -> SVector{3}
getB_tokamak_coil(r, a, b, ICoils, IPlasma) -> SVector{3}

Get the magnetic field from a Tokamak topology consists of 16 coils. Original: Tokamak-Fusion-Reactor

Arguments

  • r: location, vector of length 3 [m]

  • x,y,z: location [m]

  • a: radius of each coil in [m].

  • b: radius of central region in [m].

  • ICoil: current in the coil times number of windings in [A].

  • IPlasma: current of the plasma in [A].

source
Magnetostatics.getB_tokamak_profile Function
julia
getB_tokamak_profile(x, y, z, q_profile, a, R₀, Bζ0) :: SVector{3}
getB_tokamak_profile(r, q_profile, a, R₀, Bζ0) :: SVector{3}

Reconstruct the magnetic field distribution from a safe factor(q) profile. Reference: Tokamak, 4th Edition, John Wesson.

Arguments

  • r: location, vector of length 3 [m]

  • x,y,z: location [m]

  • q_profile::Function: profile of q. The variable of this function must be the normalized radius.

  • a: minor radius [m].

  • R₀: major radius [m].

  • Bζ0: toroidal magnetic field on axis [T].

source
Magnetostatics.getB_zpinch Function
julia
getB_zpinch(x, y, z, I, a) -> SVector{3}
getB_zpinch(r, I, a) -> SVector{3}

Get magnetic field from a Z-pinch configuration. Reference: Z-pinch

Arguments

  • r: location, vector of length 3 [m]

  • x,y,z: location [m]

  • I::Float: current in the wire [A].

  • a::Float: radius of the wire [m].

source