1#![allow(clippy::similar_names)]
10#![allow(clippy::many_single_char_names)]
11#![allow(clippy::unreadable_literal)]
12
13use crate::error::{check_coordinates, check_elevation_angle};
14use crate::math::{
15 acos, asin, atan, atan2, cos, degrees_to_radians, mul_add, normalize_degrees_0_to_360,
16 polynomial, powi, radians_to_degrees, rem_euclid, sin, tan,
17};
18use crate::time::JulianDate;
19use crate::{Horizon, RefractionCorrection, Result, SolarPosition};
20
21pub mod coefficients;
22use coefficients::{
23 NUTATION_COEFFS, OBLIQUITY_COEFFS, TERMS_B, TERMS_L, TERMS_PE, TERMS_R, TERMS_Y,
24};
25
26#[cfg(feature = "chrono")]
27use chrono::{offset::Offset, DateTime, Datelike, NaiveDate, TimeZone};
28
29const ABERRATION_CONSTANT: f64 = -20.4898;
31
32const EARTH_FLATTENING_FACTOR: f64 = 0.99664719;
34
35const EARTH_RADIUS_METERS: f64 = 6378140.0;
37
38const SECONDS_PER_HOUR: f64 = 3600.0;
40
41#[cfg(feature = "chrono")]
88#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
89#[allow(clippy::needless_pass_by_value)]
90pub fn solar_position<Tz: TimeZone>(
91 datetime: DateTime<Tz>,
92 latitude: f64,
93 longitude: f64,
94 elevation: f64,
95 delta_t: f64,
96 refraction: Option<RefractionCorrection>,
97) -> Result<SolarPosition> {
98 let jd = JulianDate::from_datetime(&datetime, delta_t)?;
99 solar_position_from_julian(jd, latitude, longitude, elevation, refraction)
100}
101
102pub fn solar_position_from_julian(
138 jd: JulianDate,
139 latitude: f64,
140 longitude: f64,
141 elevation: f64,
142 refraction: Option<RefractionCorrection>,
143) -> Result<SolarPosition> {
144 let time_dependent = spa_time_dependent_from_julian(jd)?;
145 spa_with_time_dependent_parts(latitude, longitude, elevation, refraction, &time_dependent)
146}
147
148#[derive(Debug, Clone)]
153pub struct SpaTimeDependent {
154 pub(crate) r: f64,
156 pub(crate) nu_degrees: f64,
158 pub(crate) alpha_degrees: f64,
160 pub(crate) delta_degrees: f64,
162}
163
164impl SpaTimeDependent {
165 #[must_use]
167 pub const fn earth_radius_vector(&self) -> f64 {
168 self.r
169 }
170
171 #[must_use]
173 pub const fn right_ascension(&self) -> f64 {
174 self.alpha_degrees
175 }
176
177 #[must_use]
179 pub const fn declination(&self) -> f64 {
180 self.delta_degrees
181 }
182}
183
184#[derive(Debug, Clone, Copy)]
185struct DeltaPsiEpsilon {
186 delta_psi: f64,
187 delta_epsilon: f64,
188}
189
190fn calculate_lbr_terms(jme: f64, term_coeffs: &[&[&[f64; 3]]]) -> [f64; 6] {
192 let mut lbr_terms = [0.0; 6];
195
196 for (i, term_set) in term_coeffs.iter().enumerate().take(6) {
197 let mut lbr_sum = 0.0;
198 for term in *term_set {
199 lbr_sum += term[0] * cos(mul_add(term[2], jme, term[1]));
200 }
201 lbr_terms[i] = lbr_sum;
202 }
203
204 lbr_terms
205}
206
207fn calculate_lbr_polynomial(jme: f64, terms: &[f64], num_terms: usize) -> f64 {
209 polynomial(&terms[..num_terms], jme) / 1e8
210}
211
212fn lbr_to_normalized_degrees(jme: f64, terms: &[f64], num_terms: usize) -> f64 {
214 normalize_degrees_0_to_360(radians_to_degrees(calculate_lbr_polynomial(
215 jme, terms, num_terms,
216 )))
217}
218
219fn calculate_nutation_terms(jce: f64) -> [f64; 5] {
221 [
224 polynomial(NUTATION_COEFFS[0], jce),
225 polynomial(NUTATION_COEFFS[1], jce),
226 polynomial(NUTATION_COEFFS[2], jce),
227 polynomial(NUTATION_COEFFS[3], jce),
228 polynomial(NUTATION_COEFFS[4], jce),
229 ]
230}
231
232fn calculate_delta_psi_epsilon(jce: f64, x: &[f64]) -> DeltaPsiEpsilon {
234 let mut delta_psi = 0.0;
235 let mut delta_epsilon = 0.0;
236
237 for (i, pe_term) in TERMS_PE.iter().enumerate() {
238 let mut xj_yterm_sum = 0.0;
239 for (j, &x_val) in x.iter().enumerate() {
240 xj_yterm_sum += x_val * f64::from(TERMS_Y[i][j]);
241 }
242 let xj_yterm_sum = degrees_to_radians(xj_yterm_sum);
243
244 let delta_psi_contrib = mul_add(pe_term[1], jce, pe_term[0]) * sin(xj_yterm_sum);
246 let delta_epsilon_contrib = mul_add(pe_term[3], jce, pe_term[2]) * cos(xj_yterm_sum);
247
248 delta_psi += delta_psi_contrib;
249 delta_epsilon += delta_epsilon_contrib;
250 }
251
252 DeltaPsiEpsilon {
253 delta_psi: delta_psi / 36_000_000.0,
254 delta_epsilon: delta_epsilon / 36_000_000.0,
255 }
256}
257
258fn calculate_true_obliquity_of_ecliptic(jd: &JulianDate, delta_epsilon: f64) -> f64 {
260 let epsilon0 = polynomial(OBLIQUITY_COEFFS, jd.julian_ephemeris_millennium() / 10.0);
261 epsilon0 / 3600.0 + delta_epsilon
262}
263
264fn calculate_apparent_sidereal_time_at_greenwich(
266 jd: &JulianDate,
267 delta_psi: f64,
268 epsilon_degrees: f64,
269) -> f64 {
270 let nu0_degrees = normalize_degrees_0_to_360(mul_add(
271 powi(jd.julian_century(), 2),
272 0.000387933 - jd.julian_century() / 38710000.0,
273 mul_add(
274 360.98564736629f64,
275 jd.julian_date() - 2451545.0,
276 280.46061837,
277 ),
278 ));
279
280 mul_add(
281 delta_psi,
282 cos(degrees_to_radians(epsilon_degrees)),
283 nu0_degrees,
284 )
285}
286
287fn calculate_geocentric_sun_right_ascension(
289 beta_rad: f64,
290 epsilon_rad: f64,
291 lambda_rad: f64,
292) -> f64 {
293 let alpha = atan2(
294 mul_add(
295 sin(lambda_rad),
296 cos(epsilon_rad),
297 -(tan(beta_rad) * sin(epsilon_rad)),
298 ),
299 cos(lambda_rad),
300 );
301 normalize_degrees_0_to_360(radians_to_degrees(alpha))
302}
303
304fn calculate_geocentric_sun_declination(beta_rad: f64, epsilon_rad: f64, lambda_rad: f64) -> f64 {
306 asin(mul_add(
307 sin(beta_rad),
308 cos(epsilon_rad),
309 cos(beta_rad) * sin(epsilon_rad) * sin(lambda_rad),
310 ))
311}
312
313pub fn sunrise_sunset_utc(
354 year: i32,
355 month: u32,
356 day: u32,
357 latitude: f64,
358 longitude: f64,
359 delta_t: f64,
360 elevation_angle: f64,
361) -> Result<crate::SunriseResult<crate::HoursUtc>> {
362 check_coordinates(latitude, longitude)?;
363 check_elevation_angle(elevation_angle)?;
364
365 let jd_midnight = JulianDate::from_utc(year, month, day, 0, 0, 0.0, delta_t)?;
367
368 Ok(calculate_sunrise_sunset_core(
370 jd_midnight,
371 latitude,
372 longitude,
373 delta_t,
374 elevation_angle,
375 ))
376}
377
378pub fn sunrise_sunset_utc_for_horizon(
420 year: i32,
421 month: u32,
422 day: u32,
423 latitude: f64,
424 longitude: f64,
425 delta_t: f64,
426 horizon: crate::Horizon,
427) -> Result<crate::SunriseResult<crate::HoursUtc>> {
428 sunrise_sunset_utc(
429 year,
430 month,
431 day,
432 latitude,
433 longitude,
434 delta_t,
435 horizon.elevation_angle(),
436 )
437}
438
439#[cfg(feature = "chrono")]
484#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
485#[allow(clippy::needless_pass_by_value)]
486pub fn sunrise_sunset<Tz: TimeZone>(
487 date: DateTime<Tz>,
488 latitude: f64,
489 longitude: f64,
490 delta_t: f64,
491 elevation_angle: f64,
492) -> Result<crate::SunriseResult<DateTime<Tz>>> {
493 check_coordinates(latitude, longitude)?;
494
495 let tz = date.timezone();
496 let local_date = date.date_naive();
497 let base_utc_date_guess = local_date;
501 let (_base_utc_date, result) =
502 select_utc_date_by_transit(local_date, base_utc_date_guess, |d| {
503 let converted = match sunrise_sunset_utc(
504 d.year(),
505 d.month(),
506 d.day(),
507 latitude,
508 longitude,
509 delta_t,
510 elevation_angle,
511 )? {
512 crate::SunriseResult::RegularDay {
513 sunrise,
514 transit,
515 sunset,
516 } => crate::SunriseResult::RegularDay {
517 sunrise: hours_utc_to_datetime(&tz, d, sunrise),
518 transit: hours_utc_to_datetime(&tz, d, transit),
519 sunset: hours_utc_to_datetime(&tz, d, sunset),
520 },
521 crate::SunriseResult::AllDay { transit } => crate::SunriseResult::AllDay {
522 transit: hours_utc_to_datetime(&tz, d, transit),
523 },
524 crate::SunriseResult::AllNight { transit } => crate::SunriseResult::AllNight {
525 transit: hours_utc_to_datetime(&tz, d, transit),
526 },
527 };
528
529 let transit_local_date = converted.transit().date_naive();
530
531 Ok((transit_local_date, converted))
532 })?;
533
534 Ok(ensure_events_bracket_transit(result))
535}
536
537fn precompute_sunrise_sunset_for_jd_midnight(jd_midnight: JulianDate) -> (f64, [AlphaDelta; 3]) {
539 let jce_day = jd_midnight.julian_ephemeris_century();
541 let x_terms = calculate_nutation_terms(jce_day);
542 let delta_psi_epsilon = calculate_delta_psi_epsilon(jce_day, &x_terms);
543 let epsilon_degrees =
544 calculate_true_obliquity_of_ecliptic(&jd_midnight, delta_psi_epsilon.delta_epsilon);
545 let nu_degrees = calculate_apparent_sidereal_time_at_greenwich(
546 &jd_midnight,
547 delta_psi_epsilon.delta_psi,
548 epsilon_degrees,
549 );
550
551 let mut alpha_deltas = [AlphaDelta {
553 alpha: 0.0,
554 delta: 0.0,
555 }; 3];
556 for (i, alpha_delta) in alpha_deltas.iter_mut().enumerate() {
557 let current_jd = jd_midnight.add_days((i as f64) - 1.0);
558 let current_jme = current_jd.julian_ephemeris_millennium();
559 let current_jce = current_jd.julian_ephemeris_century();
560 let current_x_terms = calculate_nutation_terms(current_jce);
561 let current_delta_psi_epsilon = calculate_delta_psi_epsilon(current_jce, ¤t_x_terms);
562 let current_epsilon_degrees = calculate_true_obliquity_of_ecliptic(
563 ¤t_jd,
564 current_delta_psi_epsilon.delta_epsilon,
565 );
566 *alpha_delta = calculate_alpha_delta(
567 current_jme,
568 current_delta_psi_epsilon.delta_psi,
569 current_epsilon_degrees,
570 );
571 }
572
573 (nu_degrees, alpha_deltas)
574}
575
576fn calculate_sunrise_sunset_hours_with_precomputed(
577 latitude: f64,
578 longitude: f64,
579 delta_t: f64,
580 elevation_angle: f64,
581 nu_degrees: f64,
582 alpha_deltas: [AlphaDelta; 3],
583) -> crate::SunriseResult<crate::HoursUtc> {
584 let m0 = (alpha_deltas[1].alpha - longitude - nu_degrees) / 360.0;
585 let transit_m = rem_euclid(m0, 1.0);
586 let phi = degrees_to_radians(latitude);
587 let delta1_rad = degrees_to_radians(alpha_deltas[1].delta);
588 let elevation_rad = degrees_to_radians(elevation_angle);
589 let acos_arg =
590 mul_add(sin(phi), -sin(delta1_rad), sin(elevation_rad)) / (cos(phi) * cos(delta1_rad));
591
592 let polar_transit_hours =
593 calculate_transit_hours(transit_m, longitude, delta_t, nu_degrees, &alpha_deltas);
594
595 if acos_arg < -1.0 {
596 return crate::SunriseResult::AllDay {
597 transit: polar_transit_hours,
598 };
599 }
600 if acos_arg > 1.0 {
601 return crate::SunriseResult::AllNight {
602 transit: polar_transit_hours,
603 };
604 }
605
606 let h0_degrees = radians_to_degrees(acos(acos_arg));
607 let m_values = [
608 transit_m,
609 rem_euclid(m0 - h0_degrees / 360.0, 1.0),
610 rem_euclid(m0 + h0_degrees / 360.0, 1.0),
611 ];
612
613 let (t_frac, r_frac, s_frac) = calculate_final_time_fractions(
614 m_values,
615 nu_degrees,
616 delta_t,
617 latitude,
618 longitude,
619 elevation_angle,
620 alpha_deltas,
621 );
622
623 let transit_hours = crate::HoursUtc::from_hours(t_frac * 24.0);
624 let sunrise_hours = crate::HoursUtc::from_hours(r_frac * 24.0);
625 let sunset_hours = crate::HoursUtc::from_hours(s_frac * 24.0);
626
627 crate::SunriseResult::RegularDay {
628 sunrise: sunrise_hours,
629 transit: transit_hours,
630 sunset: sunset_hours,
631 }
632}
633
634fn calculate_sunrise_sunset_core(
638 jd_midnight: JulianDate,
639 latitude: f64,
640 longitude: f64,
641 delta_t: f64,
642 elevation_angle: f64,
643) -> crate::SunriseResult<crate::HoursUtc> {
644 let (nu_degrees, alpha_deltas) = precompute_sunrise_sunset_for_jd_midnight(jd_midnight);
645 calculate_sunrise_sunset_hours_with_precomputed(
646 latitude,
647 longitude,
648 delta_t,
649 elevation_angle,
650 nu_degrees,
651 alpha_deltas,
652 )
653}
654
655fn calculate_transit_hours(
661 transit_m: f64,
662 longitude: f64,
663 delta_t: f64,
664 nu_degrees: f64,
665 alpha_deltas: &[AlphaDelta; 3],
666) -> crate::HoursUtc {
667 let transit_nu = mul_add(360.985647f64, transit_m, nu_degrees);
668 let transit_n = [transit_m + delta_t / 86400.0, 0.0, 0.0];
669 let transit_alpha_delta = calculate_interpolated_alpha_deltas(alpha_deltas, &transit_n)[0];
670 let transit_h_prime = limit_h_prime(transit_nu + longitude - transit_alpha_delta.alpha);
671 crate::HoursUtc::from_hours((transit_m - transit_h_prime / 360.0) * 24.0)
672}
673
674fn calculate_final_time_fractions(
677 m_values: [f64; 3],
678 nu_degrees: f64,
679 delta_t: f64,
680 latitude: f64,
681 longitude: f64,
682 elevation_angle: f64,
683 alpha_deltas: [AlphaDelta; 3],
684) -> (f64, f64, f64) {
685 let mut nu = [0.0; 3];
687 for (i, nu_item) in nu.iter_mut().enumerate() {
688 *nu_item = mul_add(360.985647f64, m_values[i], nu_degrees);
689 }
690
691 let mut n = [0.0; 3];
693 for (i, n_item) in n.iter_mut().enumerate() {
694 *n_item = m_values[i] + delta_t / 86400.0;
695 }
696
697 let alpha_delta_primes = calculate_interpolated_alpha_deltas(&alpha_deltas, &n);
699
700 let mut h_prime = [0.0; 3];
702 for i in 0..3 {
703 let h_prime_i = nu[i] + longitude - alpha_delta_primes[i].alpha;
704 h_prime[i] = limit_h_prime(h_prime_i);
705 }
706
707 let phi = degrees_to_radians(latitude);
709 let mut h = [0.0; 3];
710 for i in 0..3 {
711 let delta_prime_rad = degrees_to_radians(alpha_delta_primes[i].delta);
712 h[i] = radians_to_degrees(asin(mul_add(
713 sin(phi),
714 sin(delta_prime_rad),
715 cos(phi) * cos(delta_prime_rad) * cos(degrees_to_radians(h_prime[i])),
716 )));
717 }
718
719 let t = m_values[0] - h_prime[0] / 360.0;
721 let r = m_values[1]
722 + (h[1] - elevation_angle)
723 / (360.0
724 * cos(degrees_to_radians(alpha_delta_primes[1].delta))
725 * cos(phi)
726 * sin(degrees_to_radians(h_prime[1])));
727 let s = m_values[2]
728 + (h[2] - elevation_angle)
729 / (360.0
730 * cos(degrees_to_radians(alpha_delta_primes[2].delta))
731 * cos(phi)
732 * sin(degrees_to_radians(h_prime[2])));
733
734 (t, r, s)
735}
736
737fn calculate_interpolated_alpha_deltas(
739 alpha_deltas: &[AlphaDelta; 3],
740 n: &[f64; 3],
741) -> [AlphaDelta; 3] {
742 let a = limit_if_necessary(alpha_deltas[1].alpha - alpha_deltas[0].alpha);
743 let a_prime = limit_if_necessary(alpha_deltas[1].delta - alpha_deltas[0].delta);
744
745 let b = limit_if_necessary(alpha_deltas[2].alpha - alpha_deltas[1].alpha);
746 let b_prime = limit_if_necessary(alpha_deltas[2].delta - alpha_deltas[1].delta);
747
748 let c = b - a;
749 let c_prime = b_prime - a_prime;
750
751 let mut alpha_delta_primes = [AlphaDelta {
752 alpha: 0.0,
753 delta: 0.0,
754 }; 3];
755 for i in 0..3 {
756 alpha_delta_primes[i].alpha =
757 alpha_deltas[1].alpha + (n[i] * (mul_add(c, n[i], a + b))) / 2.0;
758 alpha_delta_primes[i].delta =
759 alpha_deltas[1].delta + (n[i] * (mul_add(c_prime, n[i], a_prime + b_prime))) / 2.0;
760 }
761 alpha_delta_primes
762}
763
764#[derive(Debug, Clone, Copy)]
765struct AlphaDelta {
766 alpha: f64,
767 delta: f64,
768}
769
770#[cfg(feature = "chrono")]
816#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
817pub fn sunrise_sunset_for_horizon<Tz: TimeZone>(
818 date: DateTime<Tz>,
819 latitude: f64,
820 longitude: f64,
821 delta_t: f64,
822 horizon: Horizon,
823) -> Result<crate::SunriseResult<DateTime<Tz>>> {
824 sunrise_sunset(
825 date,
826 latitude,
827 longitude,
828 delta_t,
829 horizon.elevation_angle(),
830 )
831}
832
833fn calculate_alpha_delta(jme: f64, delta_psi: f64, epsilon_degrees: f64) -> AlphaDelta {
836 let b_terms = calculate_lbr_terms(jme, TERMS_B);
840 let b_degrees = lbr_to_normalized_degrees(jme, &b_terms, TERMS_B.len());
841
842 let r_terms = calculate_lbr_terms(jme, TERMS_R);
844 let r = calculate_lbr_polynomial(jme, &r_terms, TERMS_R.len());
845 assert!(
846 r != 0.0,
847 "Earth radius vector is zero - astronomical impossibility"
848 );
849
850 let l_terms = calculate_lbr_terms(jme, TERMS_L);
852 let l_degrees = lbr_to_normalized_degrees(jme, &l_terms, TERMS_L.len());
853
854 let theta_degrees = normalize_degrees_0_to_360(l_degrees + 180.0);
856
857 let beta_degrees = -b_degrees;
859 let beta = degrees_to_radians(beta_degrees);
860 let epsilon = degrees_to_radians(epsilon_degrees);
861
862 let delta_tau = ABERRATION_CONSTANT / (SECONDS_PER_HOUR * r);
864
865 let lambda_degrees = theta_degrees + delta_psi + delta_tau;
867 let lambda = degrees_to_radians(lambda_degrees);
868
869 let alpha_degrees = calculate_geocentric_sun_right_ascension(beta, epsilon, lambda);
871 let delta_degrees =
872 radians_to_degrees(calculate_geocentric_sun_declination(beta, epsilon, lambda));
873
874 AlphaDelta {
875 alpha: alpha_degrees,
876 delta: delta_degrees,
877 }
878}
879
880#[cfg(feature = "chrono")]
881fn select_utc_date_by_transit<V, F>(
882 local_date: NaiveDate,
883 mut utc_date: NaiveDate,
884 mut compute: F,
885) -> Result<(NaiveDate, V)>
886where
887 F: FnMut(NaiveDate) -> Result<(NaiveDate, V)>,
888{
889 let (transit_local_date, value) = compute(utc_date)?;
890 if transit_local_date == local_date {
891 return Ok((utc_date, value));
892 }
893
894 utc_date = if transit_local_date > local_date {
895 utc_date.pred_opt().unwrap_or(utc_date)
896 } else {
897 utc_date.succ_opt().unwrap_or(utc_date)
898 };
899
900 let (transit_local_date, value) = compute(utc_date)?;
901 debug_assert_eq!(transit_local_date, local_date);
902 Ok((utc_date, value))
903}
904
905#[cfg(feature = "chrono")]
906fn hours_utc_to_datetime<Tz: TimeZone>(
907 tz: &Tz,
908 base_utc_date: NaiveDate,
909 hours: crate::HoursUtc,
910) -> DateTime<Tz> {
911 let base_utc_midnight = base_utc_date
912 .and_hms_opt(0, 0, 0)
913 .expect("midnight is always valid")
914 .and_utc();
915
916 let millis_plus = (hours.hours() * 3_600_000.0) as i64;
919 let utc_dt = base_utc_midnight + chrono::Duration::milliseconds(millis_plus);
920
921 tz.from_utc_datetime(&utc_dt.naive_utc())
922}
923
924#[cfg(feature = "chrono")]
925fn ensure_events_bracket_transit<Tz: TimeZone>(
926 result: crate::SunriseResult<DateTime<Tz>>,
927) -> crate::SunriseResult<DateTime<Tz>> {
928 let crate::SunriseResult::RegularDay {
929 mut sunrise,
930 transit,
931 mut sunset,
932 } = result
933 else {
934 return result;
935 };
936
937 if transit.offset().fix().local_minus_utc() == 0 {
941 return crate::SunriseResult::RegularDay {
942 sunrise,
943 transit,
944 sunset,
945 };
946 }
947
948 if sunrise > transit {
950 sunrise -= chrono::Duration::days(1);
951 }
952
953 if sunset < transit {
954 sunset += chrono::Duration::days(1);
955 }
956
957 crate::SunriseResult::RegularDay {
958 sunrise,
959 transit,
960 sunset,
961 }
962}
963
964fn limit_if_necessary(val: f64) -> f64 {
966 if val.abs() > 2.0 {
967 rem_euclid(val, 1.0)
968 } else {
969 val
970 }
971}
972
973fn limit_h_prime(h_prime: f64) -> f64 {
975 let limited = rem_euclid(h_prime, 360.0);
976 if limited > 180.0 {
977 limited - 360.0
978 } else {
979 limited
980 }
981}
982
983#[cfg(feature = "chrono")]
1038#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
1039#[allow(clippy::needless_pass_by_value)]
1040pub fn sunrise_sunset_multiple<Tz, H>(
1041 date: DateTime<Tz>,
1042 latitude: f64,
1043 longitude: f64,
1044 delta_t: f64,
1045 horizons: H,
1046) -> impl Iterator<Item = Result<(Horizon, crate::SunriseResult<DateTime<Tz>>)>>
1047where
1048 Tz: TimeZone,
1049 H: IntoIterator<Item = Horizon>,
1050{
1051 let tz = date.timezone();
1052 let local_date = date.date_naive();
1053 let base_utc_date_guess = local_date;
1057
1058 let precomputed = (|| -> Result<_> {
1060 check_coordinates(latitude, longitude)?;
1061 let (base_utc_date, (nu_degrees, alpha_deltas)) =
1062 select_utc_date_by_transit(local_date, base_utc_date_guess, |d| {
1063 let jd_midnight =
1064 JulianDate::from_utc(d.year(), d.month(), d.day(), 0, 0, 0.0, delta_t)?;
1065
1066 let (nu_degrees, alpha_deltas) =
1067 precompute_sunrise_sunset_for_jd_midnight(jd_midnight);
1068 let transit_m = rem_euclid(
1069 (alpha_deltas[1].alpha - longitude - nu_degrees) / 360.0,
1070 1.0,
1071 );
1072 let transit_hours = calculate_transit_hours(
1073 transit_m,
1074 longitude,
1075 delta_t,
1076 nu_degrees,
1077 &alpha_deltas,
1078 );
1079
1080 let transit_local_date = hours_utc_to_datetime(&tz, d, transit_hours).date_naive();
1081 Ok((transit_local_date, (nu_degrees, alpha_deltas)))
1082 })?;
1083
1084 Ok((base_utc_date, nu_degrees, alpha_deltas))
1085 })();
1086
1087 horizons.into_iter().map(move |horizon| {
1088 let (base_utc_date, nu_degrees, alpha_deltas) = precomputed.clone()?;
1089 let hours_result = calculate_sunrise_sunset_hours_with_precomputed(
1090 latitude,
1091 longitude,
1092 delta_t,
1093 horizon.elevation_angle(),
1094 nu_degrees,
1095 alpha_deltas,
1096 );
1097
1098 let result = ensure_events_bracket_transit(match hours_result {
1099 crate::SunriseResult::RegularDay {
1100 sunrise,
1101 transit,
1102 sunset,
1103 } => crate::SunriseResult::RegularDay {
1104 sunrise: hours_utc_to_datetime(&tz, base_utc_date, sunrise),
1105 transit: hours_utc_to_datetime(&tz, base_utc_date, transit),
1106 sunset: hours_utc_to_datetime(&tz, base_utc_date, sunset),
1107 },
1108 crate::SunriseResult::AllDay { transit } => crate::SunriseResult::AllDay {
1109 transit: hours_utc_to_datetime(&tz, base_utc_date, transit),
1110 },
1111 crate::SunriseResult::AllNight { transit } => crate::SunriseResult::AllNight {
1112 transit: hours_utc_to_datetime(&tz, base_utc_date, transit),
1113 },
1114 });
1115
1116 Ok((horizon, result))
1117 })
1118}
1119
1120#[cfg(feature = "chrono")]
1158#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
1159#[allow(clippy::needless_pass_by_value)]
1160pub fn spa_time_dependent_parts<Tz: TimeZone>(
1161 datetime: DateTime<Tz>,
1162 delta_t: f64,
1163) -> Result<SpaTimeDependent> {
1164 let jd = JulianDate::from_datetime(&datetime, delta_t)?;
1165 spa_time_dependent_from_julian(jd)
1166}
1167
1168pub fn spa_time_dependent_from_julian(jd: JulianDate) -> Result<SpaTimeDependent> {
1178 let jme = jd.julian_ephemeris_millennium();
1179 let jce = jd.julian_ephemeris_century();
1180
1181 let l_terms = calculate_lbr_terms(jme, TERMS_L);
1183 let l_degrees = lbr_to_normalized_degrees(jme, &l_terms, TERMS_L.len());
1184
1185 let b_terms = calculate_lbr_terms(jme, TERMS_B);
1187 let b_degrees = lbr_to_normalized_degrees(jme, &b_terms, TERMS_B.len());
1188
1189 let r_terms = calculate_lbr_terms(jme, TERMS_R);
1191 let r = calculate_lbr_polynomial(jme, &r_terms, TERMS_R.len());
1192
1193 assert!(
1195 r != 0.0,
1196 "Earth radius vector is zero - astronomical impossibility"
1197 );
1198
1199 let theta_degrees = normalize_degrees_0_to_360(l_degrees + 180.0);
1201 let beta_degrees = -b_degrees;
1203
1204 let x_terms = calculate_nutation_terms(jce);
1206 let delta_psi_epsilon = calculate_delta_psi_epsilon(jce, &x_terms);
1207
1208 let epsilon_degrees =
1210 calculate_true_obliquity_of_ecliptic(&jd, delta_psi_epsilon.delta_epsilon);
1211
1212 let delta_tau = ABERRATION_CONSTANT / (SECONDS_PER_HOUR * r);
1214
1215 let lambda_degrees = theta_degrees + delta_psi_epsilon.delta_psi + delta_tau;
1217
1218 let nu_degrees = calculate_apparent_sidereal_time_at_greenwich(
1220 &jd,
1221 delta_psi_epsilon.delta_psi,
1222 epsilon_degrees,
1223 );
1224
1225 let beta = degrees_to_radians(beta_degrees);
1227 let epsilon = degrees_to_radians(epsilon_degrees);
1228 let lambda = degrees_to_radians(lambda_degrees);
1229 let alpha_degrees = calculate_geocentric_sun_right_ascension(beta, epsilon, lambda);
1230
1231 let delta_degrees =
1233 radians_to_degrees(calculate_geocentric_sun_declination(beta, epsilon, lambda));
1234
1235 Ok(SpaTimeDependent {
1236 r,
1237 nu_degrees,
1238 alpha_degrees,
1239 delta_degrees,
1240 })
1241}
1242
1243pub fn spa_with_time_dependent_parts(
1280 latitude: f64,
1281 longitude: f64,
1282 elevation: f64,
1283 refraction: Option<RefractionCorrection>,
1284 time_dependent: &SpaTimeDependent,
1285) -> Result<SolarPosition> {
1286 check_coordinates(latitude, longitude)?;
1287
1288 let nu_degrees = time_dependent.nu_degrees;
1291
1292 let h_degrees =
1294 normalize_degrees_0_to_360(nu_degrees + longitude - time_dependent.alpha_degrees);
1295 let h = degrees_to_radians(h_degrees);
1296
1297 let xi_degrees = 8.794 / (3600.0 * time_dependent.r);
1299 let xi = degrees_to_radians(xi_degrees);
1300 let phi = degrees_to_radians(latitude);
1301 let delta = degrees_to_radians(time_dependent.delta_degrees);
1302
1303 let u = atan(EARTH_FLATTENING_FACTOR * tan(phi));
1304 let y = mul_add(
1305 EARTH_FLATTENING_FACTOR,
1306 sin(u),
1307 (elevation / EARTH_RADIUS_METERS) * sin(phi),
1308 );
1309 let x = mul_add(elevation / EARTH_RADIUS_METERS, cos(phi), cos(u));
1310
1311 let delta_alpha_prime_degrees = radians_to_degrees(atan2(
1312 -x * sin(xi) * sin(h),
1313 mul_add(x * sin(xi), -cos(h), cos(delta)),
1314 ));
1315
1316 let delta_prime = radians_to_degrees(atan2(
1317 mul_add(y, -sin(xi), sin(delta)) * cos(degrees_to_radians(delta_alpha_prime_degrees)),
1318 mul_add(x * sin(xi), -cos(h), cos(delta)),
1319 ));
1320
1321 let h_prime_degrees = h_degrees - delta_alpha_prime_degrees;
1323
1324 let zenith_angle = radians_to_degrees(acos(mul_add(
1326 sin(degrees_to_radians(latitude)),
1327 sin(degrees_to_radians(delta_prime)),
1328 cos(degrees_to_radians(latitude))
1329 * cos(degrees_to_radians(delta_prime))
1330 * cos(degrees_to_radians(h_prime_degrees)),
1331 )));
1332
1333 let azimuth = normalize_degrees_0_to_360(
1335 180.0
1336 + radians_to_degrees(atan2(
1337 sin(degrees_to_radians(h_prime_degrees)),
1338 cos(degrees_to_radians(h_prime_degrees)) * sin(degrees_to_radians(latitude))
1339 - tan(degrees_to_radians(delta_prime)) * cos(degrees_to_radians(latitude)),
1340 )),
1341 );
1342
1343 let elevation_angle = 90.0 - zenith_angle;
1345 let final_zenith = refraction.map_or(zenith_angle, |correction| {
1346 if elevation_angle > Horizon::SunriseSunset.elevation_angle() {
1347 let pressure = correction.pressure();
1348 let temperature = correction.temperature();
1349 zenith_angle
1351 - (pressure / 1010.0) * (283.0 / (273.0 + temperature)) * 1.02
1352 / (60.0
1353 * tan(degrees_to_radians(
1354 elevation_angle + 10.3 / (elevation_angle + 5.11),
1355 )))
1356 } else {
1357 zenith_angle
1358 }
1359 });
1360
1361 SolarPosition::new(azimuth, final_zenith)
1362}
1363
1364#[cfg(all(test, feature = "chrono", feature = "std"))]
1365mod tests {
1366 use super::*;
1367 use chrono::{DateTime, FixedOffset};
1368
1369 fn angular_distance(a: f64, b: f64) -> f64 {
1370 let diff = (a - b).abs();
1371 diff.min(360.0 - diff)
1372 }
1373
1374 #[test]
1375 fn test_spa_basic_functionality() {
1376 let datetime = "2023-06-21T12:00:00Z"
1377 .parse::<DateTime<FixedOffset>>()
1378 .unwrap();
1379
1380 let result = solar_position(
1381 datetime,
1382 37.7749, -122.4194,
1384 0.0,
1385 69.0,
1386 Some(RefractionCorrection::new(1013.25, 15.0).unwrap()),
1387 );
1388
1389 assert!(result.is_ok());
1390 let position = result.unwrap();
1391 assert!(position.azimuth() >= 0.0 && position.azimuth() <= 360.0);
1392 assert!(position.zenith_angle() >= 0.0 && position.zenith_angle() <= 180.0);
1393 }
1394
1395 #[test]
1396 fn test_time_dependent_tracks_seasonal_geometry() {
1397 let june_solstice = spa_time_dependent_from_julian(
1398 JulianDate::from_utc(2023, 6, 21, 12, 0, 0.0, 69.0).unwrap(),
1399 )
1400 .unwrap();
1401 let december_solstice = spa_time_dependent_from_julian(
1402 JulianDate::from_utc(2023, 12, 22, 12, 0, 0.0, 69.0).unwrap(),
1403 )
1404 .unwrap();
1405 let march_equinox = spa_time_dependent_from_julian(
1406 JulianDate::from_utc(2023, 3, 20, 12, 0, 0.0, 69.0).unwrap(),
1407 )
1408 .unwrap();
1409
1410 assert!(june_solstice.declination() > 23.0);
1411 assert!(june_solstice.declination() < 24.0);
1412 assert!(angular_distance(june_solstice.right_ascension(), 90.0) < 2.0);
1413
1414 assert!(december_solstice.declination() < -23.0);
1415 assert!(december_solstice.declination() > -24.0);
1416 assert!(angular_distance(december_solstice.right_ascension(), 270.0) < 2.0);
1417
1418 assert!(march_equinox.declination().abs() < 1.0);
1419 }
1420
1421 #[test]
1422 fn test_time_dependent_earth_radius_vector_changes_over_year() {
1423 let near_perihelion = spa_time_dependent_from_julian(
1424 JulianDate::from_utc(2023, 1, 4, 12, 0, 0.0, 69.0).unwrap(),
1425 )
1426 .unwrap();
1427 let near_aphelion = spa_time_dependent_from_julian(
1428 JulianDate::from_utc(2023, 7, 4, 12, 0, 0.0, 69.0).unwrap(),
1429 )
1430 .unwrap();
1431
1432 assert!(near_perihelion.earth_radius_vector() > 0.98);
1433 assert!(near_perihelion.earth_radius_vector() < 0.99);
1434 assert!(near_aphelion.earth_radius_vector() > 1.01);
1435 assert!(near_aphelion.earth_radius_vector() < 1.02);
1436 assert!(near_aphelion.earth_radius_vector() > near_perihelion.earth_radius_vector());
1437 }
1438
1439 #[test]
1440 fn test_sunrise_sunset_multiple() {
1441 let datetime = "2023-06-21T12:00:00Z"
1442 .parse::<DateTime<FixedOffset>>()
1443 .unwrap();
1444 let horizons = [
1445 Horizon::SunriseSunset,
1446 Horizon::CivilTwilight,
1447 Horizon::NauticalTwilight,
1448 ];
1449
1450 let results: Result<Vec<_>> = sunrise_sunset_multiple(
1451 datetime,
1452 37.7749, -122.4194, 69.0, horizons.iter().copied(),
1456 )
1457 .collect();
1458
1459 let results = results.unwrap();
1460
1461 assert_eq!(results.len(), 3);
1463
1464 for expected_horizon in horizons {
1466 assert!(results.iter().any(|(h, _)| *h == expected_horizon));
1467 }
1468
1469 for (horizon, bulk_result) in &results {
1471 let individual_result =
1472 sunrise_sunset_for_horizon(datetime, 37.7749, -122.4194, 69.0, *horizon).unwrap();
1473
1474 match (&individual_result, bulk_result) {
1476 (
1477 crate::SunriseResult::RegularDay {
1478 sunrise: s1,
1479 transit: t1,
1480 sunset: ss1,
1481 },
1482 crate::SunriseResult::RegularDay {
1483 sunrise: s2,
1484 transit: t2,
1485 sunset: ss2,
1486 },
1487 ) => {
1488 assert_eq!(s1, s2);
1489 assert_eq!(t1, t2);
1490 assert_eq!(ss1, ss2);
1491 }
1492 (
1493 crate::SunriseResult::AllDay { transit: t1 },
1494 crate::SunriseResult::AllDay { transit: t2 },
1495 )
1496 | (
1497 crate::SunriseResult::AllNight { transit: t1 },
1498 crate::SunriseResult::AllNight { transit: t2 },
1499 ) => {
1500 assert_eq!(t1, t2);
1501 }
1502 _ => panic!("Bulk and individual results differ in type for {horizon:?}"),
1503 }
1504 }
1505 }
1506
1507 #[test]
1508 fn test_sunrise_sunset_multiple_polar_consistency() {
1509 let datetime = "2023-06-21T12:00:00Z"
1510 .parse::<DateTime<FixedOffset>>()
1511 .unwrap();
1512
1513 let individual = sunrise_sunset_for_horizon(
1514 datetime,
1515 80.0, 0.0,
1517 69.0,
1518 Horizon::SunriseSunset,
1519 )
1520 .unwrap();
1521
1522 let bulk_results: Result<Vec<_>> =
1523 sunrise_sunset_multiple(datetime, 80.0, 0.0, 69.0, [Horizon::SunriseSunset]).collect();
1524
1525 let (_, bulk) = bulk_results.unwrap().into_iter().next().unwrap();
1526
1527 match (bulk, individual) {
1528 (
1529 crate::SunriseResult::AllDay { transit: t1 },
1530 crate::SunriseResult::AllDay { transit: t2 },
1531 )
1532 | (
1533 crate::SunriseResult::AllNight { transit: t1 },
1534 crate::SunriseResult::AllNight { transit: t2 },
1535 ) => assert_eq!(t1, t2),
1536 _ => panic!("expected matching polar-day/night results between bulk and individual"),
1537 }
1538 }
1539
1540 #[test]
1541 fn test_spa_no_refraction() {
1542 let datetime = "2023-06-21T12:00:00Z"
1543 .parse::<DateTime<FixedOffset>>()
1544 .unwrap();
1545
1546 let result = solar_position(datetime, 37.7749, -122.4194, 0.0, 69.0, None);
1547
1548 assert!(result.is_ok());
1549 let position = result.unwrap();
1550 assert!(position.azimuth() >= 0.0 && position.azimuth() <= 360.0);
1551 assert!(position.zenith_angle() >= 0.0 && position.zenith_angle() <= 180.0);
1552 }
1553
1554 #[test]
1555 fn test_spa_coordinate_validation() {
1556 let datetime = "2023-06-21T12:00:00Z"
1557 .parse::<DateTime<FixedOffset>>()
1558 .unwrap();
1559
1560 assert!(solar_position(
1562 datetime,
1563 95.0,
1564 0.0,
1565 0.0,
1566 0.0,
1567 Some(RefractionCorrection::new(1013.25, 15.0).unwrap())
1568 )
1569 .is_err());
1570
1571 assert!(solar_position(
1573 datetime,
1574 0.0,
1575 185.0,
1576 0.0,
1577 0.0,
1578 Some(RefractionCorrection::new(1013.25, 15.0).unwrap())
1579 )
1580 .is_err());
1581 }
1582
1583 #[test]
1584 fn test_sunrise_sunset_basic() {
1585 let date = "2023-06-21T00:00:00Z"
1586 .parse::<DateTime<FixedOffset>>()
1587 .unwrap();
1588
1589 let result = sunrise_sunset(date, 37.7749, -122.4194, 69.0, -0.833);
1590 assert!(result.is_ok());
1591
1592 let result =
1593 sunrise_sunset_for_horizon(date, 37.7749, -122.4194, 69.0, Horizon::SunriseSunset);
1594 assert!(result.is_ok());
1595 }
1596
1597 #[test]
1598 fn test_horizon_enum() {
1599 assert_eq!(Horizon::SunriseSunset.elevation_angle(), -0.83337);
1600 assert_eq!(Horizon::CivilTwilight.elevation_angle(), -6.0);
1601 assert_eq!(Horizon::Custom(-10.5).elevation_angle(), -10.5);
1602 }
1603}