use num_traits::Float;
use crate::{carlson::elliprj_unchecked, crate_util::check, ellipf, ellipk, StrErr};
#[numeric_literals::replace_float_literals(T::from(literal).unwrap())]
pub fn heuman_lambda<T: Float>(phi: T, m: T) -> Result<T, StrErr> {
let ans = heuman_lambda_unchecked(phi, m);
#[cfg(not(feature = "test_force_fail"))]
if ans.is_finite() {
return Ok(ans);
}
check!(@nan, heuman_lambda, [phi, m]);
if m < 0.0 || m >= 1.0 {
return Err("heuman_lambda: m must satisfy 0.0 ≤ m < 1.0.");
}
check!(@inf, heuman_lambda, [phi]);
Err("heuman_lambda: Unexpected error.")
}
#[numeric_literals::replace_float_literals(T::from(literal).unwrap())]
pub fn heuman_lambda_unchecked<T: Float>(phi: T, m: T) -> T {
if m <= 0.0 {
if m == 0.0 {
return phi.sin();
}
return nan!();
}
let n = (phi / pi_2!()).round();
if (phi - n * pi_2!()).abs() < epsilon!() {
return n;
}
let mc = 1.0 - m;
let f = ellipf(phi, mc).unwrap_or(nan!());
let k_m = ellipk(m).unwrap_or(nan!());
let k_mc = ellipk(mc).unwrap_or(nan!());
let zeta = jacobi_zeta_unchecked_k(phi, mc, k_mc);
f / k_mc + k_m * zeta / pi_2!()
}
#[inline]
#[numeric_literals::replace_float_literals(T::from(literal).unwrap())]
fn jacobi_zeta_unchecked_k<T: Float>(phi: T, m: T, k: T) -> T {
let sign = phi.signum();
let phi = phi.abs();
let sinp = phi.sin();
let cosp = phi.cos();
let mc = 1.0 - m;
let c2p = cosp * cosp;
let one_m_ms2p = mc + m * c2p;
sign * m * sinp * cosp * one_m_ms2p.sqrt() * elliprj_unchecked(0.0, mc, 1.0, one_m_ms2p)
/ (3.0 * k)
}
#[cfg(not(feature = "test_force_fail"))]
#[cfg(test)]
mod tests {
use super::*;
use crate::compare_test_data_wolfram;
#[test]
fn test_heuman_lambda() {
compare_test_data_wolfram!("heuman_lambda_data.csv", heuman_lambda, 2, 3e-15);
}
#[test]
fn test_heuman_lambda_special_cases() {
use std::f64::{
consts::{FRAC_PI_2, PI},
INFINITY, NAN,
};
assert_eq!(heuman_lambda(FRAC_PI_2, 0.5).unwrap(), 1.0);
assert_eq!(heuman_lambda(PI, 0.5).unwrap(), 2.0);
assert_eq!(heuman_lambda(3.0 * FRAC_PI_2, 0.5).unwrap(), 3.0);
assert_eq!(
heuman_lambda(1.0, 1.5),
Err("heuman_lambda: m must satisfy 0.0 ≤ m < 1.0.")
);
assert_eq!(heuman_lambda(1.0, 0.0).unwrap(), 1.0.sin());
assert_eq!(
heuman_lambda(1.0, -1.0),
Err("heuman_lambda: m must satisfy 0.0 ≤ m < 1.0.")
);
assert_eq!(
heuman_lambda(NAN, 1.0),
Err("heuman_lambda: Arguments cannot be NAN.")
);
assert_eq!(
heuman_lambda(1.0, NAN),
Err("heuman_lambda: Arguments cannot be NAN.")
);
assert_eq!(
heuman_lambda(INFINITY, 0.5),
Err("heuman_lambda: phi cannot be infinite.")
);
assert_eq!(
heuman_lambda(-INFINITY, 0.5),
Err("heuman_lambda: phi cannot be infinite.")
);
}
}
#[cfg(feature = "test_force_fail")]
crate::test_force_unreachable! {
assert_eq!(heuman_lambda(0.5, 0.5), Err("heuman_lambda: Unexpected error."));
}