API Reference
SpheroidalWaveFunctions.Fmn — Method
Fmn(m, n, dr)Compute the normalization factor Fₘₙ for radial spheroidal wave functions.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberdr::Vector: expansion coefficients
Returns
- Normalization factor Fₘₙ
SpheroidalWaveFunctions.Nmn — Method
Nmn(m, n, dr)Compute the normalization constant Nₘₙ(c) for spheroidal wave functions.
Calculates the normalization constant using the expansion coefficients and factorial recursion relations.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberdr::Vector: expansion coefficients
Returns
- Normalization constant Nₘₙ(c)
Note
The normalization ensures proper orthogonality and norm conditions for the spheroidal wave functions.
SpheroidalWaveFunctions.Nᵣᵐ1 — Method
Nᵣᵐ1(m, n, c, λ, Nmax)Compute the backward continued fraction for the characteristic value equation.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterλ::Number: trial characteristic valueNmax::Integer: maximum number of terms in continued fraction
Returns
- Value of the backward continued fraction
SpheroidalWaveFunctions.Nᵣᵐ1_prime — Method
Nᵣᵐ1_prime(m, n, c, λ, Nmax)Compute the derivative of the backward continued fraction with respect to λ.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterλ::Number: trial characteristic valueNmax::Integer: maximum number of terms
Returns
- Derivative dNᵣᵐ1/dλ
SpheroidalWaveFunctions.Nᵣᵐ2 — Method
Nᵣᵐ2(m, n, c, λ, Nmax)Compute the forward continued fraction for the characteristic value equation.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterλ::Number: trial characteristic valueNmax::Integer: maximum number of terms in continued fraction
Returns
- Value of the forward continued fraction
SpheroidalWaveFunctions.Nᵣᵐ2_prime — Method
Nᵣᵐ2_prime(m, n, c, λ, Nmax)Compute the derivative of the forward continued fraction with respect to λ.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterλ::Number: trial characteristic valueNmax::Integer: maximum number of terms
Returns
- Derivative dNᵣᵐ2/dλ
SpheroidalWaveFunctions.assoc_legendre_Pm — Method
assoc_legendre_Pm(m, N, x)Compute associated Legendre functions Pₙᵐ(x) and their derivatives for n = m, m+1, ..., N.
Uses upward recurrence relations from DLMF for numerical stability.
Arguments
m::Integer: order (m ≥ 0)N::Integer: maximum degree (N ≥ m)x::Number: argument in [-1, 1]
Returns
P::Vector: values [Pₘᵐ(x), Pₘ₊₁ᵐ(x), ..., Pₙᵐ(x)]dP::Vector: derivatives [dPₘᵐ/dx, dPₘ₊₁ᵐ/dx, ..., dPₙᵐ/dx]
Method
- DLMF 14.10.1: Initial value Pₘᵐ(x) = (-1)ᵐ (2m-1)!! (1-x²)^(m/2)
- DLMF 14.10.2: Pₘ₊₁ᵐ(x) = x(2m+1)Pₘᵐ(x)
- DLMF 14.10.3: Upward recurrence for n ≥ m+1
- DLMF 14.10.5: Derivative formula
Note
For large m, numerical overflow may occur. Consider using a specialized library.
SpheroidalWaveFunctions.asymp_λ_pro — Method
asymp_λ_pro(m, n, c)Compute the characteristic value using asymptotic expansion for large c (oblate case).
Uses the asymptotic expansion for large values of the spheroidal parameter c in the prolate case.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameter (should be large)
Returns
- Approximate characteristic value λₘₙ(c) for large c
Reference
Page 754, Abramowitz and Stegun, "Handbook of Mathematical Functions"
SpheroidalWaveFunctions.build_tridiagonal — Method
build_tridiagonal(m, n, c, N)Construct the tridiagonal matrix for the three-term recurrence relation of spheroidal wave functions.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterN::Integer: matrix dimension
Returns
- Tridiagonal matrix with elements αᵣ, βᵣ, γᵣ
SpheroidalWaveFunctions.compute_c2k — Method
compute_c2k(m, n, dr)Compute the power series expansion coefficients c₂ₖ from the Legendre expansion coefficients dᵣ.
Transforms the expansion coefficients from the Legendre function basis to the power series basis (1-x²)ᵏ, which is useful for certain computational approaches.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberdr::Vector: Legendre expansion coefficients
Returns
- Vector of power series coefficients c₂ₖ
Notes
The recurrence relations used are:
- (2m+r)!/r! → (2m+r+1)(2m+r+2)/((r+1)(r+2))
- (-r/2)ₖ → (-r/2 + 1) / (-r/2 + k - 1)
- (m+r/2+0.5)ₖ → (m + r/2 + 0.5 + k)/(m + r/2 + 0.5)
SpheroidalWaveFunctions.compute_dr2_mix — Function
compute_dr2_mix(m, n, c, λ, max_terms=25 + div(n - m, 2) + round(Int, abs(c)))Compute the expansion coefficients dᵣ for the spheroidal wave function using a mixed method.
Combines backward and forward recursion to compute the expansion coefficients in the Legendre function expansion. The coefficients satisfy a three-term recurrence relation.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterλ::Number: characteristic valuemax_terms::Int: maximum number of expansion terms (default: adaptive)
Returns
- Vector of normalized expansion coefficients dᵣ
Method
Uses backward recursion from the tail and forward recursion from the start, joining at an optimal point to maintain numerical stability.
SpheroidalWaveFunctions.compute_normalization_constant — Method
compute_normalization_constant(m, n, is_even)Compute the normalization constant based on the parity of (n-m).
Arguments
m::Integer: azimuthal ordern::Integer: mode numberis_even::Bool: true if (n-m) is even
Returns
- Normalization constant
SpheroidalWaveFunctions.compute_normalization_sum — Method
compute_normalization_sum(dr, r_values, m, is_even)Compute the normalization sum for the expansion coefficients.
Arguments
dr::Vector: expansion coefficientsr_values::Range: range of r indicesm::Integer: azimuthal orderis_even::Bool: true if (n-m) is even
Returns
- Normalization sum value
SpheroidalWaveFunctions.cv_matrix — Function
cv_matrix(m, n, c, N=10 + div(n - m, 2) + round(Int, abs(c)))Compute the characteristic value using the matrix eigenvalue method.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterN::Integer: matrix dimension (default: adaptive based on n, m, c)
Returns
- Characteristic value λₘₙ(c)
SpheroidalWaveFunctions.cv_matrix_seq — Function
cv_matrix_seq(m, n, c, N=10 + div(n - m, 2) + round(Int, abs(c)))Compute a sequence of characteristic values using the matrix eigenvalue method.
Arguments
m::Integer: azimuthal ordern::Integer: maximum mode numberc::Number: spheroidal parameterN::Integer: matrix dimension (default: adaptive based on n, m, c)
Returns
- Vector of characteristic values
SpheroidalWaveFunctions.expr_1 — Method
expr_1(n)Compute (2n)!/n! = (n+1)(n+2)⋯(2n).
Arguments
n::Integer: input value
Returns
- Value of (2n)!/n!
SpheroidalWaveFunctions.expr_2 — Method
expr_2(n)Compute (2n+2)!/(2(n+1)!).
Arguments
n::Integer: input value
Returns
- Value of (2n+2)!/(2(n+1)!)
SpheroidalWaveFunctions.expr_3 — Method
expr_3(m, n)Compute the normalization factor for even (n-m) case.
Calculates (n+m)! / ((n-m)/2)! / 2^(n-m).
Arguments
m::Integer: azimuthal ordern::Integer: mode number
Returns
- Normalization factor for even (n-m)
SpheroidalWaveFunctions.expr_4 — Method
expr_4(m, n)Compute the normalization factor for odd (n-m) case.
Calculates (n+m+1)! / ((n-m-1)/2)! / 2^(n-m).
Arguments
m::Integer: azimuthal ordern::Integer: mode number
Returns
- Normalization factor for odd (n-m)
SpheroidalWaveFunctions.find_eigenvalue — Method
find_eigenvalue(m, n, c; tol=1e-12, max_iter=25, num_terms=100)Compute the characteristic value using Newton-Raphson iteration with continued fractions.
Uses the matrix method to obtain an initial guess, then refines it using Newton-Raphson iteration on the continued fraction equation.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parametertol::Float64: convergence tolerance (default: 1e-12)max_iter::Int: maximum number of iterations (default: 25)num_terms::Int: number of terms in continued fraction (default: 100)
Returns
- Characteristic value λₘₙ(c)
Examples
λ = find_eigenvalue(1, 2, 0.5)
λ = find_eigenvalue(1, 2, 0.5, tol=1e-14, max_iter=50)SpheroidalWaveFunctions.find_eigenvalue_guess — Method
find_eigenvalue_guess(m, n, c, λ_guess; tol=1e-12, max_iter=25, num_terms=40)Refine a characteristic value using Newton-Raphson iteration with a given initial guess.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterλ_guess::Number: initial guess for the characteristic valuetol::Float64: convergence tolerance (default: 1e-12)max_iter::Int: maximum number of iterations (default: 25)num_terms::Int: number of terms in continued fraction (default: 40)
Returns
- Refined characteristic value λₘₙ(c)
SpheroidalWaveFunctions.find_eigenvalue_seq — Method
find_eigenvalue_seq(m, n, c; max_iter=20, num_terms=30)Compute a sequence of characteristic values for modes m through n.
Arguments
m::Integer: azimuthal ordern::Integer: maximum mode numberc::Number: spheroidal parametermax_iter::Int: maximum iterations per eigenvalue (default: 20)num_terms::Int: number of terms in continued fraction (default: 30)
Returns
- Vector of characteristic values [λₘₘ, λₘ₊₁,ₘ, ..., λₙₘ]
SpheroidalWaveFunctions.kmn1 — Method
kmn1(m, n, c)Compute the normalization constant kₘₙ for the first radial function.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameter
Returns
- Normalization constant kₘₙ
Reference
DLMF: (30.11.10) / (30.11.11)
SpheroidalWaveFunctions.oblate_angular_leg — Method
oblate_angular_leg(m, n, c, x)
oblate_angular_leg(m, n, c, λ, x)
oblate_angular_leg(m, n, c, λ, c2k, x)Compute the oblate angular spheroidal wave function and its derivative using the Legendre expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)c2k::Vector: expansion coefficients (optional, computed if not provided)x::Number: argument in [-1, 1]
Returns
S: value of the angular function Sₘₙ(-ic, x)dS: derivative dSₘₙ/dx
Examples
S, dS = oblate_angular_leg(1, 2, 0.5, 0.3)SpheroidalWaveFunctions.oblate_angular_ps — Method
oblate_angular_ps(m, n, c, x)
oblate_angular_ps(m, n, c, λ, x)
oblate_angular_ps(m, n, c, λ, c2k, x)Compute the oblate angular spheroidal wave function and its derivative using the power series expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)c2k::Vector: expansion coefficients (optional, computed if not provided)x::Number: argument in [-1, 1]
Returns
S: value of the angular function Sₘₙ(-ic, x)dS: derivative dSₘₙ/dx
Examples
S, dS = oblate_angular_ps(1, 2, 0.5, 0.3)SpheroidalWaveFunctions.oblate_cv — Method
oblate_cv(m, n, c)Compute the characteristic value (eigenvalue) for oblate spheroidal wave functions.
The characteristic value λₘₙ(-ic) is the eigenvalue of the differential equation for the angular spheroidal wave function in the oblate case.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameter
Returns
λ: characteristic value λₘₙ(-ic)
Examples
λ = oblate_cv(1, 2, 0.5)SpheroidalWaveFunctions.oblate_cv_seq — Method
oblate_cv_seq(m, n, c)Compute a sequence of characteristic values for oblate spheroidal wave functions.
Returns the characteristic values λₘₘ(-ic), λₘ₊₁,ₘ(-ic), ..., λₙₘ(-ic) for all modes from m to n.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: maximum mode number (n ≥ m)c::Number: spheroidal parameter
Returns
- Vector of characteristic values [λₘₘ(-ic), λₘ₊₁,ₘ(-ic), ..., λₙₘ(-ic)]
Examples
λ_seq = oblate_cv_seq(1, 5, 0.5) # Returns [λ₁₁, λ₂₁, λ₃₁, λ₄₁, λ₅₁]SpheroidalWaveFunctions.oblate_radial1 — Method
oblate_radial1(m, n, c, ξ)
oblate_radial1(m, n, c, λ, ξ)
oblate_radial1(m, n, c, λ, dr, ξ)Compute the oblate radial spheroidal wave function of the first kind and its derivative.
The radial function is regular at the origin and corresponds to the modified spherical Bessel function expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)dr::Vector: expansion coefficients (optional, computed if not provided)ξ::Number: radial argument (ξ ≥ 1)
Returns
R: value of the radial function Rₘₙ⁽¹⁾(-ic, iξ)dR: derivative dRₘₙ⁽¹⁾/dξ
Examples
R, dR = oblate_radial1(1, 2, 0.5, 1.5)SpheroidalWaveFunctions.oblate_radial2 — Method
oblate_radial2(m, n, c, ξ)
oblate_radial2(m, n, c, λ, ξ)
oblate_radial2(m, n, c, λ, dr, ξ)Compute the oblate radial spheroidal wave function of the second kind and its derivative.
The radial function is irregular at the origin and corresponds to the modified spherical Neumann function expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)dr::Vector: expansion coefficients (optional, computed if not provided)ξ::Number: radial argument (ξ ≥ 1)
Returns
R: value of the radial function Rₘₙ⁽²⁾(-ic, iξ)dR: derivative dRₘₙ⁽²⁾/dξ
Examples
R, dR = oblate_radial2(1, 2, 0.5, 1.5)SpheroidalWaveFunctions.pochhammer — Method
pochhammer(a, k)Compute the Pochhammer symbol (rising factorial) (a)ₖ = a(a+1)(a+2)⋯(a+k-1).
Arguments
a::Number: base valuek::Integer: number of terms
Returns
- Pochhammer symbol (a)ₖ
SpheroidalWaveFunctions.prolate_angular_leg — Method
prolate_angular_leg(m, n, c, x)
prolate_angular_leg(m, n, c, λ, x)
prolate_angular_leg(m, n, c, λ, dr, x)Compute the prolate angular spheroidal wave function and its derivative using the Legendre expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)dr::Vector: expansion coefficients (optional, computed if not provided)x::Number: argument in [-1, 1]
Returns
S: value of the angular function Sₘₙ(c, x)dS: derivative dSₘₙ/dx
Examples
S, dS = prolate_angular_leg(1, 2, 0.5, 0.3)SpheroidalWaveFunctions.prolate_angular_ps — Method
prolate_angular_ps(m, n, c, x)
prolate_angular_ps(m, n, c, λ, x)
prolate_angular_ps(m, n, c, λ, c2k, x)Compute the prolate angular spheroidal wave function and its derivative using the power series expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)c2k::Vector: expansion coefficients (optional, computed if not provided)x::Number: argument in [-1, 1]
Returns
S: value of the angular function Sₘₙ(c, x)dS: derivative dSₘₙ/dx
Examples
S, dS = prolate_angular_ps(1, 2, 0.5, 0.3)SpheroidalWaveFunctions.prolate_cv — Method
prolate_cv(m, n, c)Compute the characteristic value (eigenvalue) for prolate spheroidal wave functions.
The characteristic value λₘₙ(c) is the eigenvalue of the differential equation for the angular spheroidal wave function.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameter
Returns
λ: characteristic value λₘₙ(c)
Examples
λ = prolate_cv(1, 2, 0.5)SpheroidalWaveFunctions.prolate_cv_seq — Method
prolate_cv_seq(m, n, c)Compute a sequence of characteristic values for prolate spheroidal wave functions.
Returns the characteristic values λₘₘ(c), λₘ₊₁,ₘ(c), ..., λₙₘ(c) for all modes from m to n.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: maximum mode number (n ≥ m)c::Number: spheroidal parameter
Returns
- Vector of characteristic values [λₘₘ(c), λₘ₊₁,ₘ(c), ..., λₙₘ(c)]
Examples
λ_seq = prolate_cv_seq(1, 5, 0.5) # Returns [λ₁₁, λ₂₁, λ₃₁, λ₄₁, λ₅₁]SpheroidalWaveFunctions.prolate_radial1 — Method
prolate_radial1(m, n, c, ξ)
prolate_radial1(m, n, c, λ, ξ)
prolate_radial1(m, n, c, λ, dr, ξ)Compute the prolate radial spheroidal wave function of the first kind and its derivative.
The radial function is regular at the origin and corresponds to the spherical Bessel function expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)dr::Vector: expansion coefficients (optional, computed if not provided)ξ::Number: radial argument (ξ ≥ 1)
Returns
R: value of the radial function Rₘₙ⁽¹⁾(c, ξ)dR: derivative dRₘₙ⁽¹⁾/dξ
Examples
R, dR = prolate_radial1(1, 2, 0.5, 1.5)SpheroidalWaveFunctions.prolate_radial2 — Method
prolate_radial2(m, n, c, ξ)
prolate_radial2(m, n, c, λ, ξ)
prolate_radial2(m, n, c, λ, dr, ξ)Compute the prolate radial spheroidal wave function of the second kind and its derivative.
The radial function is irregular at the origin and corresponds to the spherical Neumann function expansion.
Arguments
m::Integer: azimuthal order (m ≥ 0)n::Integer: mode number (n ≥ m)c::Number: spheroidal parameterλ::Number: characteristic value (optional, computed if not provided)dr::Vector: expansion coefficients (optional, computed if not provided)ξ::Number: radial argument (ξ ≥ 1)
Returns
R: value of the radial function Rₘₙ⁽²⁾(c, ξ)dR: derivative dRₘₙ⁽²⁾/dξ
Examples
R, dR = prolate_radial2(1, 2, 0.5, 1.5)SpheroidalWaveFunctions.ps_λ — Method
ps_λ(m, n, c)Compute the characteristic value using power series expansion for small c.
Uses the asymptotic expansion λₘₙ(c) ≈ n(n+1) + ℓ₂c² + ℓ₄c⁴ + ℓ₆c⁶ for small values of the spheroidal parameter c.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameter (should be small)
Returns
- Approximate characteristic value λₘₙ(c)
Reference
Page 754, Abramowitz and Stegun, "Handbook of Mathematical Functions"
SpheroidalWaveFunctions.sphericalbesselj_derivative — Method
sphericalbesselj_derivative(ν, z)Compute the derivative of the spherical Bessel function jᵥ(z) with respect to z.
Uses the recurrence relation: jᵥ'(z) = jᵥ₋₁(z) - (ν+1)/z jᵥ(z).
Arguments
ν::Number: orderz::Number: argument
Returns
- Derivative djᵥ/dz
SpheroidalWaveFunctions.sphericalbessely_derivative — Method
sphericalbessely_derivative(ν, z)Compute the derivative of the spherical Neumann function yᵥ(z) with respect to z.
Uses the recurrence relation: yᵥ'(z) = yᵥ₋₁(z) - (ν+1)/z yᵥ(z).
Arguments
ν::Number: orderz::Number: argument
Returns
- Derivative dyᵥ/dz
SpheroidalWaveFunctions.spheroidal_ang_1 — Method
spheroidal_ang_1(m, n, dr, x)Compute the angular spheroidal wave function using Legendre function expansion (method 1).
Evaluates the spheroidal function as:
- Sₘₙ(x, γ²) = ∑ₖ d₂ₖ(c) Pₘ₊₂ₖᵐ(x) for even (n-m)
- Sₘₙ(x, γ²) = ∑ₖ d₂ₖ₊₁(c) Pₘ₊₂ₖ₊₁ᵐ(x) for odd (n-m)
where Pₙᵐ are associated Legendre functions.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberdr::Vector: expansion coefficients dᵣx::Number: argument in [-1, 1]
Returns
S: value of the angular function Sₘₙ(x)dS: derivative dSₘₙ/dx
SpheroidalWaveFunctions.spheroidal_ang_2 — Method
spheroidal_ang_2(m, n, c2k, x)Compute the angular spheroidal wave function using power series expansion (method 2).
Evaluates the spheroidal function as:
- Sₘₙ(x, γ²) = (-1)ᵐ (1 - x²)^(m/2) ∑ₖ c₂ₖ (1 - x²)ᵏ for even (n-m)
- Sₘₙ(x, γ²) = (-1)ᵐ x(1 - x²)^(m/2) ∑ₖ c₂ₖ (1 - x²)ᵏ for odd (n-m)
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc2k::Vector: power series coefficients c₂ₖx::Number: argument in [-1, 1]
Returns
S: value of the angular function Sₘₙ(x)dS: derivative dSₘₙ/dx
SpheroidalWaveFunctions.spheroidal_rad_1 — Method
spheroidal_rad_1(m, n, c, dr, ξ)Compute the radial spheroidal wave function of the first kind using spherical Bessel functions.
The first kind radial function is regular at ξ = 1 and is expressed as a series of spherical Bessel functions jₙ.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterdr::Vector: expansion coefficientsξ::Number: radial argument (ξ ≥ 1)
Returns
R: value of the radial function Rₘₙ⁽¹⁾(c, ξ)dR: derivative dRₘₙ⁽¹⁾/dξ
SpheroidalWaveFunctions.spheroidal_rad_2 — Method
spheroidal_rad_2(m, n, c, dr, ξ)Compute the radial spheroidal wave function of the second kind using spherical Neumann functions.
The second kind radial function is irregular at ξ = 1 and is expressed as a series of spherical Neumann functions yₙ.
Arguments
m::Integer: azimuthal ordern::Integer: mode numberc::Number: spheroidal parameterdr::Vector: expansion coefficientsξ::Number: radial argument (ξ ≥ 1)
Returns
R: value of the radial function Rₘₙ⁽²⁾(c, ξ)dR: derivative dRₘₙ⁽²⁾/dξ
SpheroidalWaveFunctions.symmetrize_tridiagonal — Method
symmetrize_tridiagonal(A)Convert a tridiagonal matrix to a symmetric tridiagonal form for eigenvalue computation.
Arguments
A::Tridiagonal: input tridiagonal matrix
Returns
- Symmetric tridiagonal matrix with off-diagonal elements √(Aᵢⱼ * Aⱼᵢ)
SpheroidalWaveFunctions.αᵣ — Method
αᵣ(m, r, γ²)Compute the upper diagonal element of the three-term recurrence relation.
Arguments
m::Integer: azimuthal orderr::Integer: indexγ²::Number: γ² = c²
Returns
- Coefficient αᵣ
SpheroidalWaveFunctions.βᵐ — Method
βᵐ(m, r, γ²)Compute the product γᵣ(m, r, γ²) * αᵣ(m, r-2, γ²) used in the continued fraction.
SpheroidalWaveFunctions.βᵣ — Method
βᵣ(m, r, γ²)Compute the diagonal element of the three-term recurrence relation.
Arguments
m::Integer: azimuthal orderr::Integer: indexγ²::Number: γ² = c²
Returns
- Coefficient βᵣ
SpheroidalWaveFunctions.γᵐ — Method
γᵐ(m, r, γ²)Alias for βᵣ used in the continued fraction formulation.
SpheroidalWaveFunctions.γᵣ — Method
γᵣ(m, r, γ²)Compute the lower diagonal element of the three-term recurrence relation.
Arguments
m::Integer: azimuthal orderr::Integer: indexγ²::Number: γ² = c²
Returns
- Coefficient γᵣ
Notes
All angular and radial functions return tuples (value, derivative):
- Angular functions:
(S, dS/dx) - Radial functions:
(R, dR/dξ)
Parameters
m::Integer: Azimuthal order (m ≥ 0)n::Integer: Mode number (n ≥ m)c::Number: Spheroidal parameterλ::Number: Characteristic value (optional, computed if not provided)dr::Vector: Expansion coefficients (optional, computed if not provided)c2k::Vector: Power series coefficients (optional, computed if not provided)x::Number: Angular argument in [-1, 1]ξ::Number: Radial argument (ξ ≥ 1)
Coordinate Systems
Prolate spheroidal coordinates: Used for elongated geometries.
Oblate spheroidal coordinates: Used for flattened geometries.