Skip to main content

celestial_core/nutation/
iau2000b.rs

1//! IAU 2000B Nutation Model
2//!
3//! This module implements the IAU 2000B nutation model, a truncated version of the full
4//! IAU 2000A model designed for applications where sub-milliarcsecond precision is not required.
5//!
6//! # Model Description
7//!
8//! IAU 2000B reduces computational cost by:
9//! - Using only the first 77 lunisolar terms (out of 678 in IAU 2000A)
10//! - Omitting the 687 planetary terms entirely
11//! - Applying fixed bias corrections to approximate the omitted planetary effects
12//!
13//! The planetary bias corrections are:
14//! - Longitude (Δψ): -0.135 milliarcseconds
15//! - Obliquity (Δε): +0.388 milliarcseconds
16//!
17//! # Accuracy
18//!
19//! The IAU 2000B model achieves accuracy of approximately 1 milliarcsecond over the
20//! period 1995-2050. This is sufficient for many practical applications including
21//! amateur telescope pointing and general ephemeris work, but insufficient for
22//! high-precision astrometry, VLBI, or pulsar timing.
23//!
24//! # Fundamental Arguments
25//!
26//! The model uses five Delaunay arguments computed from polynomial expressions
27//! in Julian centuries from J2000.0 (TDB):
28//!
29//! - `l` (mean anomaly of the Moon)
30//! - `l'` (mean anomaly of the Sun)
31//! - `F` (mean argument of latitude of the Moon)
32//! - `D` (mean elongation of the Moon from the Sun)
33//! - `Ω` (mean longitude of the Moon's ascending node)
34//!
35//! # References
36//!
37//! - McCarthy, D. D. & Luzum, B. J., "An Abridged Model of the Precession-Nutation
38//!   of the Celestial Pole", Celestial Mechanics and Dynamical Astronomy, 2003
39//! - IERS Conventions (2003), Chapter 5
40//! - SOFA Library: `iauNut00b`
41
42use super::lunisolar_terms::LUNISOLAR_TERMS;
43use super::types::NutationResult;
44use crate::constants::{
45    ARCSEC_TO_RAD, CIRCULAR_ARCSECONDS, MICROARCSEC_TO_RAD, MILLIARCSEC_TO_RAD, TWOPI,
46};
47use crate::errors::AstroResult;
48use crate::math::fmod;
49
50/// IAU 2000B nutation calculator.
51///
52/// A simplified nutation model using 77 lunisolar terms plus fixed planetary bias
53/// corrections. Provides ~1 mas accuracy, suitable for applications not requiring
54/// the full precision of [`NutationIAU2000A`](super::NutationIAU2000A).
55///
56/// # Example
57///
58/// ```
59/// use celestial_core::nutation::NutationIAU2000B;
60///
61/// let nut = NutationIAU2000B::new();
62/// // Compute nutation for J2000.0 (two-part JD: 2451545.0 + 0.0)
63/// let result = nut.compute(2451545.0, 0.0).unwrap();
64///
65/// // delta_psi and delta_eps are in radians
66/// println!("Δψ = {} rad", result.delta_psi);
67/// println!("Δε = {} rad", result.delta_eps);
68/// ```
69pub struct NutationIAU2000B;
70
71impl Default for NutationIAU2000B {
72    fn default() -> Self {
73        Self::new()
74    }
75}
76
77impl NutationIAU2000B {
78    /// Creates a new IAU 2000B nutation calculator.
79    pub fn new() -> Self {
80        Self
81    }
82
83    /// Computes nutation angles for a given Julian Date.
84    ///
85    /// # Arguments
86    ///
87    /// * `jd1` - First part of two-part Julian Date (TDB). Typically the integer part
88    ///   or J2000 epoch (2451545.0).
89    /// * `jd2` - Second part of two-part Julian Date (TDB). Typically the fractional
90    ///   part or offset from `jd1`.
91    ///
92    /// The two-part representation preserves precision. The split is arbitrary;
93    /// `jd1 + jd2` must equal the desired Julian Date.
94    ///
95    /// # Returns
96    ///
97    /// Returns a [`NutationResult`] containing:
98    /// - `delta_psi`: Nutation in longitude (radians)
99    /// - `delta_eps`: Nutation in obliquity (radians)
100    ///
101    /// Both values are IAU 2000B approximations with ~1 mas accuracy.
102    pub fn compute(&self, jd1: f64, jd2: f64) -> AstroResult<NutationResult> {
103        let t = crate::utils::jd_to_centuries(jd1, jd2);
104
105        let (delta_psi_ls, delta_eps_ls) = self.compute_lunisolar(t);
106
107        const PLANETARY_BIAS_LONGITUDE: f64 = -0.135 * MILLIARCSEC_TO_RAD;
108        const PLANETARY_BIAS_OBLIQUITY: f64 = 0.388 * MILLIARCSEC_TO_RAD;
109
110        let delta_psi = delta_psi_ls + PLANETARY_BIAS_LONGITUDE;
111        let delta_eps = delta_eps_ls + PLANETARY_BIAS_OBLIQUITY;
112        Ok(NutationResult {
113            delta_psi,
114            delta_eps,
115        })
116    }
117
118    /// Computes the lunisolar nutation contribution.
119    ///
120    /// Evaluates the first 77 terms of the lunisolar nutation series using
121    /// the five Delaunay fundamental arguments. Each term contributes
122    /// sine and cosine components to both longitude and obliquity.
123    ///
124    /// # Arguments
125    ///
126    /// * `t` - Julian centuries from J2000.0 (TDB)
127    ///
128    /// # Returns
129    ///
130    /// Tuple of (Δψ, Δε) in radians representing the lunisolar contribution
131    /// to nutation in longitude and obliquity.
132    fn compute_lunisolar(&self, t: f64) -> (f64, f64) {
133        // Delaunay arguments (arcseconds, then converted to radians)
134        // l: Mean anomaly of the Moon
135        let el = fmod(485868.249036 + 1717915923.2178 * t, CIRCULAR_ARCSECONDS) * ARCSEC_TO_RAD;
136        // l': Mean anomaly of the Sun
137        let elp = fmod(1287104.79305 + 129596581.0481 * t, CIRCULAR_ARCSECONDS) * ARCSEC_TO_RAD;
138        // F: Mean argument of latitude of the Moon
139        let f = fmod(335779.526232 + 1739527262.8478 * t, CIRCULAR_ARCSECONDS) * ARCSEC_TO_RAD;
140        // D: Mean elongation of the Moon from the Sun
141        let d = fmod(1072260.70369 + 1602961601.2090 * t, CIRCULAR_ARCSECONDS) * ARCSEC_TO_RAD;
142        // Ω: Mean longitude of the Moon's ascending node
143        let om = fmod(450160.398036 + -6962890.5431 * t, CIRCULAR_ARCSECONDS) * ARCSEC_TO_RAD;
144
145        let mut dpsi = 0.0;
146        let mut deps = 0.0;
147
148        for &(nl, nlp, nf, nd, nom, sp, spt, cp, ce, cet, se) in
149            LUNISOLAR_TERMS.iter().take(77).rev()
150        {
151            let arg = fmod(
152                (nl as f64) * el
153                    + (nlp as f64) * elp
154                    + (nf as f64) * f
155                    + (nd as f64) * d
156                    + (nom as f64) * om,
157                TWOPI,
158            );
159
160            let (sarg, carg) = libm::sincos(arg);
161
162            dpsi += (sp + spt * t) * sarg + cp * carg;
163            deps += (ce + cet * t) * carg + se * sarg;
164        }
165
166        (dpsi * MICROARCSEC_TO_RAD, deps * MICROARCSEC_TO_RAD)
167    }
168}