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>orComplex<f32> - Complete function set — J, Y, I, K, H(1), H(2), Ai, Bi
- Consecutive orders —
_seqvariants return ν, ν+1, …, ν+n−1 in one call - Exponential scaling —
_scaledvariants prevent overflow/underflow - Negative orders — supports ν < 0 via DLMF reflection formulas (not in Amos)
no_stdsupport — 3-tier: bareno_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
| Function | Description | Scaled variant returns |
|---|---|---|
besselj(ν, z)besselj_scaled(ν, z)besselj_seq(ν, z, n, scaling) | Jν(z), Bessel first kind | exp(−|Im(z)|) · Jν(z) |
bessely(ν, z)bessely_scaled(ν, z)bessely_seq(ν, z, n, scaling) | Yν(z), Bessel second kind | exp(−|Im(z)|) · Yν(z) |
besseli(ν, z)besseli_scaled(ν, z)besseli_seq(ν, z, n, scaling) | Iν(z), modified first kind | exp(−|Re(z)|) · Iν(z) |
besselk(ν, z)besselk_scaled(ν, z)besselk_seq(ν, z, n, scaling) | Kν(z), modified second kind | exp(z) · Kν(z) |
hankel1(ν, z)hankel1_scaled(ν, z)hankel1_seq(ν, z, n, scaling) | Hν(1)(z), Hankel first kind | exp(−iz) · Hν(1)(z) |
hankel2(ν, z)hankel2_scaled(ν, z)hankel2_seq(ν, z, n, scaling) | Hν(2)(z), Hankel second kind | exp(iz) · Hν(2)(z) |
airy(z)airy_scaled(z)airy_raw(z, scaling) | Ai(z), Airy first kind | exp(ζ) · Ai(z) |
airyprime(z)airyprime_scaled(z)airyprime_raw(z, scaling) | Ai′(z), derivative of Airy first kind | exp(ζ) · Ai′(z) |
biry(z)biry_scaled(z)biry_raw(z, scaling) | Bi(z), Airy second kind | exp(−|Re(ζ)|) · Bi(z) |
biryprime(z)biryprime_scaled(z)biryprime_raw(z, scaling) | Bi′(z), derivative of Airy second kind | exp(−|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.
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:
| Variant | Cause |
|---|---|
InvalidInput | z = 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 |
ConvergenceFailure | Internal algorithm did not converge |
§no_std support
| Cargo features | Available API |
|---|---|
default-features = false | 24 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
BesselFloattrait. - 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.