celestial_core/nutation/iau2000a.rs
1//! IAU 2000A nutation model.
2//!
3//! Implements the IAU 2000A nutation model as defined by the International
4//! Astronomical Union. This model computes the nutation in longitude (delta_psi)
5//! and nutation in obliquity (delta_eps) for a given epoch.
6//!
7//! ## Model Specification
8//!
9//! IAU 2000A is a nutation model based on the MHB2000 (Mathews,
10//! Herring, Buffett 2002) rigid-Earth series with:
11//!
12//! - **678 lunisolar terms**: Trigonometric series based on 5 fundamental arguments
13//! (Moon's mean anomaly, Sun's mean anomaly, Moon's argument of latitude,
14//! mean elongation of Moon from Sun, longitude of Moon's ascending node)
15//! - **687 planetary terms**: Additional terms involving planetary mean longitudes
16//! (Mercury through Neptune) and general precession in longitude
17//!
18//! ## Precision
19//!
20//! - Formal precision: ~0.1 microarcsecond (μas) for epochs near J2000.0
21//! - Accuracy degrades for epochs far from J2000.0 due to polynomial approximations
22//! - Suitable for applications requiring sub-milliarcsecond precision
23//!
24//! ## Reference
25//!
26//! - IERS Conventions (2010), Chapter 5
27//! - Mathews, Herring & Buffett (2002), J. Geophys. Res. 107, B4
28
29use super::fundamental_args::{IERS2010FundamentalArgs, MHB2000FundamentalArgs};
30use super::lunisolar_terms::LUNISOLAR_TERMS;
31use super::planetary_terms::PLANETARY_TERMS;
32use super::types::NutationResult;
33use crate::constants::{MICROARCSEC_TO_RAD, TWOPI};
34use crate::errors::AstroResult;
35use crate::math::fmod;
36
37/// IAU 2000A nutation calculator.
38///
39/// Computes nutation angles using the full IAU 2000A model with 678 lunisolar
40/// terms and 687 planetary terms. The computation follows the MHB2000 formulation
41/// with coefficients expressed in microarcseconds.
42///
43/// # Example
44///
45/// ```
46/// use celestial_core::nutation::iau2000a::NutationIAU2000A;
47///
48/// let nut = NutationIAU2000A::new();
49///
50/// // J2000.0 epoch (two-part JD for precision)
51/// let jd1 = 2451545.0;
52/// let jd2 = 0.0;
53///
54/// let result = nut.compute(jd1, jd2).unwrap();
55/// // result.delta_psi: nutation in longitude (radians)
56/// // result.delta_eps: nutation in obliquity (radians)
57/// ```
58#[derive(Debug, Clone, Copy)]
59pub struct NutationIAU2000A;
60
61impl Default for NutationIAU2000A {
62 fn default() -> Self {
63 Self::new()
64 }
65}
66
67impl NutationIAU2000A {
68 /// Creates a new IAU 2000A nutation calculator.
69 pub fn new() -> Self {
70 Self
71 }
72
73 /// Computes nutation for the given epoch.
74 ///
75 /// Evaluates both lunisolar and planetary nutation series at the specified
76 /// Julian Date, returning nutation in longitude (delta_psi) and obliquity
77 /// (delta_eps) in radians.
78 ///
79 /// # Arguments
80 ///
81 /// * `jd1` - First part of two-part Julian Date (typically the integer day)
82 /// * `jd2` - Second part of two-part Julian Date (typically the fractional day)
83 ///
84 /// The epoch is computed as `jd1 + jd2`. The two-part representation preserves
85 /// precision when the epoch is far from J2000.0.
86 ///
87 /// # Returns
88 ///
89 /// [`NutationResult`] containing:
90 /// - `delta_psi`: Nutation in longitude (radians)
91 /// - `delta_eps`: Nutation in obliquity (radians)
92 pub fn compute(&self, jd1: f64, jd2: f64) -> AstroResult<NutationResult> {
93 let t = crate::utils::jd_to_centuries(jd1, jd2);
94
95 let lunisolar_args = [
96 t.moon_mean_anomaly(),
97 t.sun_mean_anomaly_mhb(),
98 t.mean_argument_of_latitude(),
99 t.mean_elongation_mhb(),
100 t.moon_ascending_node_longitude(),
101 ];
102 let (delta_psi_ls, delta_eps_ls) = self.compute_lunisolar(&lunisolar_args, t);
103
104 let (delta_psi_planetary, delta_eps_planetary) = self.compute_planetary(t);
105
106 Ok(NutationResult {
107 delta_psi: delta_psi_planetary + delta_psi_ls,
108 delta_eps: delta_eps_planetary + delta_eps_ls,
109 })
110 }
111
112 /// Computes the lunisolar nutation contribution.
113 ///
114 /// Evaluates 678 terms of the lunisolar nutation series. Each term is a
115 /// trigonometric function of a linear combination of the five fundamental
116 /// arguments of lunisolar motion.
117 ///
118 /// The series has the form:
119 /// ```text
120 /// delta_psi = sum_i (A_i + A'_i * t) * sin(arg_i) + A''_i * cos(arg_i)
121 /// delta_eps = sum_i (B_i + B'_i * t) * cos(arg_i) + B''_i * sin(arg_i)
122 /// ```
123 ///
124 /// where `arg_i = n_l * l + n_lp * l' + n_F * F + n_D * D + n_Om * Om` and
125 /// coefficients are in microarcseconds.
126 ///
127 /// # Arguments
128 ///
129 /// * `args` - Five fundamental arguments in radians: \[l, l', F, D, Om\]
130 /// - l: Moon's mean anomaly
131 /// - l': Sun's mean anomaly
132 /// - F: Moon's argument of latitude
133 /// - D: Mean elongation of Moon from Sun
134 /// - Om: Longitude of Moon's ascending node
135 /// * `t` - Julian centuries from J2000.0 (TT)
136 ///
137 /// # Returns
138 ///
139 /// Tuple of (delta_psi, delta_eps) in radians.
140 pub fn compute_lunisolar(&self, args: &[f64; 5], t: f64) -> (f64, f64) {
141 let mut dpsi = 0.0;
142 let mut deps = 0.0;
143
144 for term in LUNISOLAR_TERMS.iter().rev() {
145 let arg = fmod(
146 (term.0 as f64) * args[0]
147 + (term.1 as f64) * args[1]
148 + (term.2 as f64) * args[2]
149 + (term.3 as f64) * args[3]
150 + (term.4 as f64) * args[4],
151 TWOPI,
152 );
153
154 let (sarg, carg) = libm::sincos(arg);
155
156 dpsi += (term.5 + term.6 * t) * sarg + term.7 * carg;
157 deps += (term.8 + term.9 * t) * carg + term.10 * sarg;
158 }
159
160 (dpsi * MICROARCSEC_TO_RAD, deps * MICROARCSEC_TO_RAD)
161 }
162
163 /// Computes the planetary nutation contribution.
164 ///
165 /// Evaluates 687 terms of the planetary nutation series. Each term depends on
166 /// the mean longitudes of the planets (Mercury through Neptune) plus the
167 /// general precession in longitude.
168 ///
169 /// The planetary series is smaller in amplitude than the lunisolar series
170 /// but essential for sub-milliarcsecond accuracy. The largest planetary terms
171 /// arise from resonances between planetary and lunar orbital periods.
172 ///
173 /// # Arguments
174 ///
175 /// * `t` - Julian centuries from J2000.0 (TT)
176 ///
177 /// # Returns
178 ///
179 /// Tuple of (delta_psi, delta_eps) in radians.
180 pub fn compute_planetary(&self, t: f64) -> (f64, f64) {
181 let al = fmod(2.35555598 + 8328.6914269554 * t, TWOPI);
182 let af = fmod(1.627905234 + 8433.466158131 * t, TWOPI);
183 let ad = fmod(5.198466741 + 7771.3771468121 * t, TWOPI);
184 let aom = fmod(2.18243920 - 33.757045 * t, TWOPI);
185 let apa = t.precession();
186
187 let alme = t.mercury_lng();
188 let alve = t.venus_lng();
189 let alea = t.earth_lng();
190 let alma = t.mars_lng();
191 let alju = t.jupiter_lng();
192 let alsa = t.saturn_lng();
193 let alur = t.uranus_lng();
194 let alne = t.neptune_longitude_mhb();
195
196 let mut dpsi = 0.0;
197 let mut deps = 0.0;
198
199 for &(nl, nf, nd, nom, nme, nve, nea, nma, nju, nsa, nur, nne, npa, sp, cp, se, ce) in
200 PLANETARY_TERMS.iter().rev()
201 {
202 let arg = fmod(
203 (nl as f64) * al
204 + (nf as f64) * af
205 + (nd as f64) * ad
206 + (nom as f64) * aom
207 + (nme as f64) * alme
208 + (nve as f64) * alve
209 + (nea as f64) * alea
210 + (nma as f64) * alma
211 + (nju as f64) * alju
212 + (nsa as f64) * alsa
213 + (nur as f64) * alur
214 + (nne as f64) * alne
215 + (npa as f64) * apa,
216 TWOPI,
217 );
218
219 let (sarg, carg) = libm::sincos(arg);
220
221 dpsi += (sp as f64) * sarg + (cp as f64) * carg;
222 deps += (se as f64) * sarg + (ce as f64) * carg;
223 }
224
225 (dpsi * MICROARCSEC_TO_RAD, deps * MICROARCSEC_TO_RAD)
226 }
227}