Skip to main content

Crate complex_bessel

Crate complex_bessel 

Source
Expand description

Pure Rust implementation of complex Bessel functions based on Amos Algorithm 644 (ACM TOMS 644).

§Features

  • f32 & f64 — all functions accept Complex<f64> or Complex<f32>
  • Complete function set — J, Y, I, K, H(1), H(2), Ai, Bi
  • Consecutive orders_seq variants return ν, ν+1, …, ν+n−1 in one call
  • Exponential scaling_scaled variants prevent overflow/underflow
  • Negative orders — supports ν < 0 via DLMF reflection formulas (not in Amos)
  • no_std support — 3-tier: bare no_std, alloc, std (default)

§Quick start

use complex_bessel::*;
use num_complex::Complex;

let z = Complex::new(1.0, 2.0);

let j = besselj(0.5, z).unwrap();
let k = besselk(1.0, z).unwrap();
let h = hankel1(0.0, z).unwrap();

// Scaled versions prevent overflow/underflow
let k_scaled = besselk_scaled(1.0, z).unwrap();

// Airy functions
let ai = airy(z).unwrap();
let ai_prime = airyprime(z).unwrap();

§Functions

FunctionDescriptionScaled variant returns
besselj(ν, z)
besselj_scaled(ν, z)
besselj_seq(ν, z, n, scaling)
Jν(z), Bessel first kindexp(−|Im(z)|) · Jν(z)
bessely(ν, z)
bessely_scaled(ν, z)
bessely_seq(ν, z, n, scaling)
Yν(z), Bessel second kindexp(−|Im(z)|) · Yν(z)
besseli(ν, z)
besseli_scaled(ν, z)
besseli_seq(ν, z, n, scaling)
Iν(z), modified first kindexp(−|Re(z)|) · Iν(z)
besselk(ν, z)
besselk_scaled(ν, z)
besselk_seq(ν, z, n, scaling)
Kν(z), modified second kindexp(z) · Kν(z)
hankel1(ν, z)
hankel1_scaled(ν, z)
hankel1_seq(ν, z, n, scaling)
Hν(1)(z), Hankel first kindexp(−iz) · Hν(1)(z)
hankel2(ν, z)
hankel2_scaled(ν, z)
hankel2_seq(ν, z, n, scaling)
Hν(2)(z), Hankel second kindexp(iz) · Hν(2)(z)
airy(z)
airy_scaled(z)
airy_raw(z, scaling)
Ai(z), Airy first kindexp(ζ) · Ai(z)
airyprime(z)
airyprime_scaled(z)
airyprime_raw(z, scaling)
Ai′(z), derivative of Airy first kindexp(ζ) · Ai′(z)
biry(z)
biry_scaled(z)
biry_raw(z, scaling)
Bi(z), Airy second kindexp(−|Re(ζ)|) · Bi(z)
biryprime(z)
biryprime_scaled(z)
biryprime_raw(z, scaling)
Bi′(z), derivative of Airy second kindexp(−|Re(ζ)|) · Bi′(z)

where ζ = (2/3) z√z.

§Exponential scaling

Each _scaled variant multiplies the result by an exponential factor that cancels the asymptotic growth (or decay), keeping values in a representable floating-point range. The exact factor is listed in the “Scaled variant returns” column above. Use scaled functions whenever |z| or |Im(z)| is large enough that the unscaled result would overflow or underflow.

§Negative orders

All functions accept any real order, including negative values. DLMF reflection formulas are applied automatically:

  • J: J−ν(z) = cos(νπ) Jν(z) − sin(νπ) Yν(z) (DLMF 10.2.3)
  • Y: Y−ν(z) = sin(νπ) Jν(z) + cos(νπ) Yν(z) (DLMF 10.2.3)
  • I: I−ν(z) = Iν(z) + (2/π) sin(νπ) Kν(z) (DLMF 10.27.2)
  • K: K−ν(z) = Kν(z) (DLMF 10.27.3)
  • H(1): H(1)−ν(z) = exp(νπi) H(1)ν(z) (DLMF 10.4.6)
  • H(2): H(2)−ν(z) = exp(−νπi) H(2)ν(z) (DLMF 10.4.6)

For integer orders, simplified identities are used (e.g., J−n(z) = (−1)n Jn(z)).

§Function variants

§Consecutive orders

The _seq variants (besselj_seq, besselk_seq, …) compute values at consecutive orders ν, ν+1, …, ν+n−1 in a single call and return a BesselResult that includes an Accuracy field. The _raw Airy variants (airy_raw, biry_raw, …) similarly return an AiryResult with Accuracy. All other functions return only the computed value without Accuracy.

Accuracy:

StatusMeaning
NormalFull machine precision
ReducedMore than half of significant digits may be lost

Reduced is extremely rare in practice. SciPy’s Bessel wrappers also silently discard the equivalent Amos IERR=3 flag by default.

To check accuracy status, use a _seq function (Bessel) or a _raw function (Airy):

use complex_bessel::*;
use num_complex::Complex;

let z = Complex::new(1.0, 2.0);
let result = besselk_seq(0.0, z, 1, Scaling::Unscaled).unwrap();
assert!(matches!(result.status, Accuracy::Normal));

let result = airy_raw(z, Scaling::Unscaled).unwrap();
assert!(matches!(result.status, Accuracy::Normal));

§Error handling

All functions return Result<_, Error>. The four error variants are:

VariantCause
InvalidInputz = 0 for K/Y/H, n < 1
Overflow|z| or |ν| too large (or |z| too small) for finite result
TotalPrecisionLoss|z| or |ν| too large for meaningful computation
ConvergenceFailureInternal algorithm did not converge

§no_std support

Cargo featuresAvailable API
default-features = false24 single-value functions
features = ["alloc"]full API (+ 6 _seq variants)
features = ["std"] (default)full API, with platform-native math

The 24 single-value functions include 12 Bessel (J/Y/I/K/H(1)/H(2) × unscaled/scaled) and 12 Airy (Ai/Ai’/Bi/Bi’ × unscaled/scaled/raw).

# Bare no_std — no allocator needed:
complex-bessel = { version = "0.1", default-features = false }

# no_std + alloc — full API:
complex-bessel = { version = "0.1", default-features = false, features = ["alloc"] }

Re-exports§

pub use machine::BesselFloat;
pub use types::BesselResult;
pub use types::Accuracy;
pub use types::AiryResult;
pub use types::Error;
pub use types::Scaling;

Modules§

machine
Machine constants and the BesselFloat trait.
types
Core types for Bessel and Airy function computation.

Functions§

airy
Airy function Ai(z).
airy_raw
Airy function Ai(z) with accuracy status.
airy_scaled
Scaled Airy function, exp(ζ) · Ai(z), where ζ = (2/3) z√z.
airyprime
Derivative of the Airy function, Ai’(z).
airyprime_raw
Derivative of the Airy function Ai’(z) with accuracy status.
airyprime_scaled
Scaled derivative of the Airy function, exp(ζ) · Ai’(z), where ζ = (2/3) z√z.
besseli
Modified Bessel function of the first kind, Iν(z).
besseli_scaled
Scaled modified Bessel function of the first kind, exp(−|Re(z)|) · Iν(z).
besseli_seq
Compute Iν+j(z) for j = 0, 1, …, n−1 in a single call.
besselj
Bessel function of the first kind, Jν(z).
besselj_scaled
Scaled Bessel function of the first kind, exp(−|Im(z)|) · Jν(z).
besselj_seq
Compute Jν+j(z) for j = 0, 1, …, n−1 in a single call.
besselk
Modified Bessel function of the second kind, Kν(z).
besselk_scaled
Scaled modified Bessel function of the second kind, exp(z) · Kν(z).
besselk_seq
Compute Kν+j(z) for j = 0, 1, …, n−1 in a single call.
bessely
Bessel function of the second kind, Yν(z).
bessely_scaled
Scaled Bessel function of the second kind, exp(−|Im(z)|) · Yν(z).
bessely_seq
Compute Yν+j(z) for j = 0, 1, …, n−1 in a single call.
biry
Airy function of the second kind, Bi(z).
biry_raw
Airy function of the second kind Bi(z) with accuracy status.
biry_scaled
Scaled Airy function of the second kind, exp(−|Re(ζ)|) · Bi(z), where ζ = (2/3) z√z.
biryprime
Derivative of the Airy function of the second kind, Bi’(z).
biryprime_raw
Derivative of the Airy function of the second kind Bi’(z) with accuracy status.
biryprime_scaled
Scaled derivative of the Airy function of the second kind, exp(−|Re(ζ)|) · Bi’(z), where ζ = (2/3) z√z.
hankel1
Hankel function of the first kind, Hν(1)(z).
hankel2
Hankel function of the second kind, Hν(2)(z).
hankel1_scaled
Scaled Hankel function of the first kind, exp(−iz) · Hν(1)(z).
hankel1_seq
Compute Hν+j(1)(z) for j = 0, 1, …, n−1 in a single call.
hankel2_scaled
Scaled Hankel function of the second kind, exp(iz) · Hν(2)(z).
hankel2_seq
Compute Hν+j(2)(z) for j = 0, 1, …, n−1 in a single call.