1#![allow(clippy::similar_names)]
10#![allow(clippy::many_single_char_names)]
11#![allow(clippy::unreadable_literal)]
12
13use crate::error::check_coordinates;
14use crate::math::{
15 acos, asin, atan, atan2, cos, degrees_to_radians, floor, mul_add, normalize_degrees_0_to_360,
16 polynomial, powi, radians_to_degrees, sin, tan,
17};
18use crate::time::JulianDate;
19#[cfg(feature = "chrono")]
20use crate::Horizon;
21use crate::{RefractionCorrection, Result, SolarPosition};
22
23pub mod coefficients;
24use coefficients::{
25 NUTATION_COEFFS, OBLIQUITY_COEFFS, TERMS_B, TERMS_L, TERMS_PE, TERMS_R, TERMS_Y,
26};
27
28#[cfg(feature = "chrono")]
29use chrono::{DateTime, LocalResult, TimeZone, Utc};
30
31const SUNRISE_SUNSET_ANGLE: f64 = -0.83337;
33
34const ABERRATION_CONSTANT: f64 = -20.4898;
36
37const EARTH_FLATTENING_FACTOR: f64 = 0.99664719;
39
40const EARTH_RADIUS_METERS: f64 = 6378140.0;
42
43const SECONDS_PER_HOUR: f64 = 3600.0;
45
46#[cfg(feature = "chrono")]
93#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
94#[allow(clippy::needless_pass_by_value)]
95pub fn solar_position<Tz: TimeZone>(
96 datetime: DateTime<Tz>,
97 latitude: f64,
98 longitude: f64,
99 elevation: f64,
100 delta_t: f64,
101 refraction: Option<RefractionCorrection>,
102) -> Result<SolarPosition> {
103 let jd = JulianDate::from_datetime(&datetime, delta_t)?;
104 solar_position_from_julian(jd, latitude, longitude, elevation, refraction)
105}
106
107pub fn solar_position_from_julian(
143 jd: JulianDate,
144 latitude: f64,
145 longitude: f64,
146 elevation: f64,
147 refraction: Option<RefractionCorrection>,
148) -> Result<SolarPosition> {
149 let time_dependent = spa_time_dependent_from_julian(jd)?;
150 spa_with_time_dependent_parts(latitude, longitude, elevation, refraction, &time_dependent)
151}
152
153#[derive(Debug, Clone)]
158pub struct SpaTimeDependent {
159 pub(crate) r: f64,
161 pub(crate) nu_degrees: f64,
163 pub(crate) alpha_degrees: f64,
165 pub(crate) delta_degrees: f64,
167}
168
169#[derive(Debug, Clone, Copy)]
170struct DeltaPsiEpsilon {
171 delta_psi: f64,
172 delta_epsilon: f64,
173}
174
175fn calculate_lbr_terms(jme: f64, term_coeffs: &[&[&[f64; 3]]]) -> [f64; 6] {
177 let mut lbr_terms = [0.0; 6];
180
181 for (i, term_set) in term_coeffs.iter().enumerate().take(6) {
182 let mut lbr_sum = 0.0;
183 for term in *term_set {
184 lbr_sum += term[0] * cos(mul_add(term[2], jme, term[1]));
185 }
186 lbr_terms[i] = lbr_sum;
187 }
188
189 lbr_terms
190}
191
192fn calculate_lbr_polynomial(jme: f64, terms: &[f64], num_terms: usize) -> f64 {
194 polynomial(&terms[..num_terms], jme) / 1e8
195}
196
197fn lbr_to_normalized_degrees(jme: f64, terms: &[f64], num_terms: usize) -> f64 {
199 normalize_degrees_0_to_360(radians_to_degrees(calculate_lbr_polynomial(
200 jme, terms, num_terms,
201 )))
202}
203
204fn calculate_nutation_terms(jce: f64) -> [f64; 5] {
206 [
209 polynomial(NUTATION_COEFFS[0], jce),
210 polynomial(NUTATION_COEFFS[1], jce),
211 polynomial(NUTATION_COEFFS[2], jce),
212 polynomial(NUTATION_COEFFS[3], jce),
213 polynomial(NUTATION_COEFFS[4], jce),
214 ]
215}
216
217fn calculate_delta_psi_epsilon(jce: f64, x: &[f64]) -> DeltaPsiEpsilon {
219 let mut delta_psi = 0.0;
220 let mut delta_epsilon = 0.0;
221
222 for (i, pe_term) in TERMS_PE.iter().enumerate() {
223 let xj_yterm_sum = degrees_to_radians(calculate_xj_yterm_sum(i, x));
224
225 let delta_psi_contrib = mul_add(pe_term[1], jce, pe_term[0]) * sin(xj_yterm_sum);
227 let delta_epsilon_contrib = mul_add(pe_term[3], jce, pe_term[2]) * cos(xj_yterm_sum);
228
229 delta_psi += delta_psi_contrib;
230 delta_epsilon += delta_epsilon_contrib;
231 }
232
233 DeltaPsiEpsilon {
234 delta_psi: delta_psi / 36_000_000.0,
235 delta_epsilon: delta_epsilon / 36_000_000.0,
236 }
237}
238
239fn calculate_xj_yterm_sum(i: usize, x: &[f64]) -> f64 {
241 let mut sum = 0.0;
242 for (j, &x_val) in x.iter().enumerate() {
243 sum += x_val * f64::from(TERMS_Y[i][j]);
244 }
245 sum
246}
247
248fn calculate_true_obliquity_of_ecliptic(jd: &JulianDate, delta_epsilon: f64) -> f64 {
250 let epsilon0 = polynomial(OBLIQUITY_COEFFS, jd.julian_ephemeris_millennium() / 10.0);
251 epsilon0 / 3600.0 + delta_epsilon
252}
253
254fn calculate_apparent_sidereal_time_at_greenwich(
256 jd: &JulianDate,
257 delta_psi: f64,
258 epsilon_degrees: f64,
259) -> f64 {
260 let nu0_degrees = normalize_degrees_0_to_360(mul_add(
261 powi(jd.julian_century(), 2),
262 0.000387933 - jd.julian_century() / 38710000.0,
263 mul_add(
264 360.98564736629f64,
265 jd.julian_date() - 2451545.0,
266 280.46061837,
267 ),
268 ));
269
270 mul_add(
271 delta_psi,
272 cos(degrees_to_radians(epsilon_degrees)),
273 nu0_degrees,
274 )
275}
276
277fn calculate_geocentric_sun_right_ascension(
279 beta_rad: f64,
280 epsilon_rad: f64,
281 lambda_rad: f64,
282) -> f64 {
283 let alpha = atan2(
284 mul_add(
285 sin(lambda_rad),
286 cos(epsilon_rad),
287 -(tan(beta_rad) * sin(epsilon_rad)),
288 ),
289 cos(lambda_rad),
290 );
291 normalize_degrees_0_to_360(radians_to_degrees(alpha))
292}
293
294fn calculate_geocentric_sun_declination(beta_rad: f64, epsilon_rad: f64, lambda_rad: f64) -> f64 {
296 asin(mul_add(
297 sin(beta_rad),
298 cos(epsilon_rad),
299 cos(beta_rad) * sin(epsilon_rad) * sin(lambda_rad),
300 ))
301}
302
303pub fn sunrise_sunset_utc(
344 year: i32,
345 month: u32,
346 day: u32,
347 latitude: f64,
348 longitude: f64,
349 delta_t: f64,
350 elevation_angle: f64,
351) -> Result<crate::SunriseResult<crate::HoursUtc>> {
352 check_coordinates(latitude, longitude)?;
353
354 let jd_midnight = JulianDate::from_utc(year, month, day, 0, 0, 0.0, delta_t)?;
356
357 calculate_sunrise_sunset_core(jd_midnight, latitude, longitude, delta_t, elevation_angle)
359}
360
361pub fn sunrise_sunset_utc_for_horizon(
402 year: i32,
403 month: u32,
404 day: u32,
405 latitude: f64,
406 longitude: f64,
407 delta_t: f64,
408 horizon: crate::Horizon,
409) -> Result<crate::SunriseResult<crate::HoursUtc>> {
410 sunrise_sunset_utc(
411 year,
412 month,
413 day,
414 latitude,
415 longitude,
416 delta_t,
417 horizon.elevation_angle(),
418 )
419}
420
421#[cfg(feature = "chrono")]
458#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
459#[allow(clippy::needless_pass_by_value)]
460pub fn sunrise_sunset<Tz: TimeZone>(
461 date: DateTime<Tz>,
462 latitude: f64,
463 longitude: f64,
464 delta_t: f64,
465 elevation_angle: f64,
466) -> Result<crate::SunriseResult<DateTime<Tz>>> {
467 check_coordinates(latitude, longitude)?;
468
469 let day_start = truncate_to_day_start(&date);
470
471 calculate_sunrise_sunset_spa_algorithm(day_start, latitude, longitude, delta_t, elevation_angle)
473}
474
475#[allow(clippy::unnecessary_wraps)]
479fn calculate_sunrise_sunset_core(
480 jd_midnight: JulianDate,
481 latitude: f64,
482 longitude: f64,
483 delta_t: f64,
484 elevation_angle: f64,
485) -> Result<crate::SunriseResult<crate::HoursUtc>> {
486 let jce_day = jd_midnight.julian_ephemeris_century();
488 let x_terms = calculate_nutation_terms(jce_day);
489 let delta_psi_epsilon = calculate_delta_psi_epsilon(jce_day, &x_terms);
490 let epsilon_degrees =
491 calculate_true_obliquity_of_ecliptic(&jd_midnight, delta_psi_epsilon.delta_epsilon);
492 let nu_degrees = calculate_apparent_sidereal_time_at_greenwich(
493 &jd_midnight,
494 delta_psi_epsilon.delta_psi,
495 epsilon_degrees,
496 );
497
498 let mut alpha_deltas = [AlphaDelta {
500 alpha: 0.0,
501 delta: 0.0,
502 }; 3];
503 for (i, alpha_delta) in alpha_deltas.iter_mut().enumerate() {
504 let current_jd = jd_midnight.add_days((i as f64) - 1.0);
505 let current_jme = current_jd.julian_ephemeris_millennium();
506 let ad = calculate_alpha_delta(current_jme, delta_psi_epsilon.delta_psi, epsilon_degrees);
507 *alpha_delta = ad;
508 }
509
510 let m0 = (alpha_deltas[1].alpha - longitude - nu_degrees) / 360.0;
512
513 let polar_type = check_polar_conditions_type(latitude, elevation_angle, alpha_deltas[1].delta);
515
516 let (m_values, _h0_degrees) =
518 calculate_approximate_times(m0, latitude, elevation_angle, alpha_deltas[1].delta);
519
520 let (t_frac, r_frac, s_frac) = calculate_final_time_fractions(
522 m_values,
523 nu_degrees,
524 delta_t,
525 latitude,
526 longitude,
527 elevation_angle,
528 alpha_deltas,
529 );
530
531 let transit_hours = crate::HoursUtc::from_hours(t_frac * 24.0);
533 let sunrise_hours = crate::HoursUtc::from_hours(r_frac * 24.0);
534 let sunset_hours = crate::HoursUtc::from_hours(s_frac * 24.0);
535
536 match polar_type {
538 Some(PolarType::AllDay) => Ok(crate::SunriseResult::AllDay {
539 transit: transit_hours,
540 }),
541 Some(PolarType::AllNight) => Ok(crate::SunriseResult::AllNight {
542 transit: transit_hours,
543 }),
544 None => Ok(crate::SunriseResult::RegularDay {
545 sunrise: sunrise_hours,
546 transit: transit_hours,
547 sunset: sunset_hours,
548 }),
549 }
550}
551
552#[cfg(feature = "chrono")]
554fn calculate_sunrise_sunset_spa_algorithm<Tz: TimeZone>(
555 day: DateTime<Tz>,
556 latitude: f64,
557 longitude: f64,
558 delta_t: f64,
559 elevation_angle: f64,
560) -> Result<crate::SunriseResult<DateTime<Tz>>> {
561 let day_start = truncate_to_day_start(&day);
562
563 let (nu_degrees, delta_psi_epsilon, epsilon_degrees) =
565 calculate_sidereal_time_and_nutation(day_start.clone(), delta_t)?;
566
567 let alpha_deltas = calculate_alpha_deltas_for_three_days(
569 day_start,
570 delta_t,
571 delta_psi_epsilon,
572 epsilon_degrees,
573 )?;
574
575 let m0 = (alpha_deltas[1].alpha - longitude - nu_degrees) / 360.0;
577
578 let polar_type = check_polar_conditions_type(latitude, elevation_angle, alpha_deltas[1].delta);
580
581 let (m_values, _h0_degrees) =
583 calculate_approximate_times(m0, latitude, elevation_angle, alpha_deltas[1].delta);
584
585 let final_result = calculate_final_times(FinalTimeParams {
587 day,
588 m_values,
589 nu_degrees,
590 delta_t,
591 latitude,
592 longitude,
593 elevation_angle,
594 alpha_deltas,
595 });
596
597 match polar_type {
599 Some(PolarType::AllDay) => {
600 if let crate::SunriseResult::RegularDay { transit, .. } = final_result {
601 Ok(crate::SunriseResult::AllDay { transit })
602 } else {
603 Ok(final_result)
604 }
605 }
606 Some(PolarType::AllNight) => {
607 if let crate::SunriseResult::RegularDay { transit, .. } = final_result {
608 Ok(crate::SunriseResult::AllNight { transit })
609 } else {
610 Ok(final_result)
611 }
612 }
613 None => Ok(final_result),
614 }
615}
616#[cfg(feature = "chrono")]
617#[allow(clippy::needless_pass_by_value)]
618fn calculate_sidereal_time_and_nutation<Tz: TimeZone>(
619 day_start: DateTime<Tz>,
620 delta_t: f64,
621) -> Result<(f64, DeltaPsiEpsilon, f64)> {
622 let jd = JulianDate::from_datetime(&day_start, delta_t)?;
623 let jce_day = jd.julian_ephemeris_century();
624 let x_terms = calculate_nutation_terms(jce_day);
625 let delta_psi_epsilon = calculate_delta_psi_epsilon(jce_day, &x_terms);
626 let epsilon_degrees =
627 calculate_true_obliquity_of_ecliptic(&jd, delta_psi_epsilon.delta_epsilon);
628 let nu_degrees = calculate_apparent_sidereal_time_at_greenwich(
629 &jd,
630 delta_psi_epsilon.delta_psi,
631 epsilon_degrees,
632 );
633 Ok((nu_degrees, delta_psi_epsilon, epsilon_degrees))
634}
635
636#[cfg(feature = "chrono")]
637#[allow(clippy::needless_pass_by_value)]
638fn calculate_alpha_deltas_for_three_days<Tz: TimeZone>(
639 day_start: DateTime<Tz>,
640 delta_t: f64,
641 delta_psi_epsilon: DeltaPsiEpsilon,
642 epsilon_degrees: f64,
643) -> Result<[AlphaDelta; 3]> {
644 let base_jd = JulianDate::from_datetime(&day_start, delta_t)?;
645
646 let mut alpha_deltas = [AlphaDelta {
647 alpha: 0.0,
648 delta: 0.0,
649 }; 3];
650 for (i, alpha_delta) in alpha_deltas.iter_mut().enumerate() {
651 let current_jd = base_jd.add_days((i as f64) - 1.0);
652 let current_jme = current_jd.julian_ephemeris_millennium();
653 let ad = calculate_alpha_delta(current_jme, delta_psi_epsilon.delta_psi, epsilon_degrees);
654 *alpha_delta = ad;
655 }
656 Ok(alpha_deltas)
657}
658
659#[derive(Debug, Clone, Copy)]
666enum PolarType {
667 AllDay,
668 AllNight,
669}
670
671fn check_polar_conditions_type(
673 latitude: f64,
674 elevation_angle: f64,
675 delta1: f64,
676) -> Option<PolarType> {
677 let phi = degrees_to_radians(latitude);
678 let elevation_rad = degrees_to_radians(elevation_angle);
679 let delta1_rad = degrees_to_radians(delta1);
680
681 let acos_arg =
682 mul_add(sin(phi), -sin(delta1_rad), sin(elevation_rad)) / (cos(phi) * cos(delta1_rad));
683
684 if acos_arg < -1.0 {
685 Some(PolarType::AllDay)
686 } else if acos_arg > 1.0 {
687 Some(PolarType::AllNight)
688 } else {
689 None
690 }
691}
692
693fn calculate_approximate_times(
695 m0: f64,
696 latitude: f64,
697 elevation_angle: f64,
698 delta1: f64,
699) -> ([f64; 3], f64) {
700 let phi = degrees_to_radians(latitude);
701 let delta1_rad = degrees_to_radians(delta1);
702 let elevation_rad = degrees_to_radians(elevation_angle);
703
704 let acos_arg =
705 mul_add(sin(phi), -sin(delta1_rad), sin(elevation_rad)) / (cos(phi) * cos(delta1_rad));
706 let h0 = acos(acos_arg);
707 let h0_degrees = radians_to_degrees(h0).min(180.0);
708
709 let mut m = [0.0; 3];
710 m[0] = normalize_to_unit_range(m0);
711 m[1] = normalize_to_unit_range(m0 - h0_degrees / 360.0);
712 m[2] = normalize_to_unit_range(m0 + h0_degrees / 360.0);
713
714 (m, h0_degrees)
715}
716
717fn calculate_final_time_fractions(
720 m_values: [f64; 3],
721 nu_degrees: f64,
722 delta_t: f64,
723 latitude: f64,
724 longitude: f64,
725 elevation_angle: f64,
726 alpha_deltas: [AlphaDelta; 3],
727) -> (f64, f64, f64) {
728 let mut nu = [0.0; 3];
730 for (i, nu_item) in nu.iter_mut().enumerate() {
731 *nu_item = mul_add(360.985647f64, m_values[i], nu_degrees);
732 }
733
734 let mut n = [0.0; 3];
736 for (i, n_item) in n.iter_mut().enumerate() {
737 *n_item = m_values[i] + delta_t / 86400.0;
738 }
739
740 let alpha_delta_primes = calculate_interpolated_alpha_deltas(&alpha_deltas, &n);
742
743 let mut h_prime = [0.0; 3];
745 for i in 0..3 {
746 let h_prime_i = nu[i] + longitude - alpha_delta_primes[i].alpha;
747 h_prime[i] = limit_h_prime(h_prime_i);
748 }
749
750 let phi = degrees_to_radians(latitude);
752 let mut h = [0.0; 3];
753 for i in 0..3 {
754 let delta_prime_rad = degrees_to_radians(alpha_delta_primes[i].delta);
755 h[i] = radians_to_degrees(asin(mul_add(
756 sin(phi),
757 sin(delta_prime_rad),
758 cos(phi) * cos(delta_prime_rad) * cos(degrees_to_radians(h_prime[i])),
759 )));
760 }
761
762 let t = m_values[0] - h_prime[0] / 360.0;
764 let r = m_values[1]
765 + (h[1] - elevation_angle)
766 / (360.0
767 * cos(degrees_to_radians(alpha_delta_primes[1].delta))
768 * cos(phi)
769 * sin(degrees_to_radians(h_prime[1])));
770 let s = m_values[2]
771 + (h[2] - elevation_angle)
772 / (360.0
773 * cos(degrees_to_radians(alpha_delta_primes[2].delta))
774 * cos(phi)
775 * sin(degrees_to_radians(h_prime[2])));
776
777 (t, r, s)
778}
779
780#[cfg(feature = "chrono")]
782fn calculate_final_times<Tz: TimeZone>(
783 params: FinalTimeParams<Tz>,
784) -> crate::SunriseResult<DateTime<Tz>> {
785 let mut nu = [0.0; 3];
787 for (i, nu_item) in nu.iter_mut().enumerate() {
788 *nu_item = mul_add(360.985647f64, params.m_values[i], params.nu_degrees);
789 }
790
791 let mut n = [0.0; 3];
793 for (i, n_item) in n.iter_mut().enumerate() {
794 *n_item = params.m_values[i] + params.delta_t / 86400.0;
795 }
796
797 let alpha_delta_primes = calculate_interpolated_alpha_deltas(¶ms.alpha_deltas, &n);
799
800 let mut h_prime = [0.0; 3];
802 for i in 0..3 {
803 let h_prime_i = nu[i] + params.longitude - alpha_delta_primes[i].alpha;
804 h_prime[i] = limit_h_prime(h_prime_i);
805 }
806
807 let phi = degrees_to_radians(params.latitude);
809 let mut h = [0.0; 3];
810 for i in 0..3 {
811 let delta_prime_rad = degrees_to_radians(alpha_delta_primes[i].delta);
812 h[i] = radians_to_degrees(asin(mul_add(
813 sin(phi),
814 sin(delta_prime_rad),
815 cos(phi) * cos(delta_prime_rad) * cos(degrees_to_radians(h_prime[i])),
816 )));
817 }
818
819 let t = params.m_values[0] - h_prime[0] / 360.0;
821 let r = params.m_values[1]
822 + (h[1] - params.elevation_angle)
823 / (360.0
824 * cos(degrees_to_radians(alpha_delta_primes[1].delta))
825 * cos(phi)
826 * sin(degrees_to_radians(h_prime[1])));
827 let s = params.m_values[2]
828 + (h[2] - params.elevation_angle)
829 / (360.0
830 * cos(degrees_to_radians(alpha_delta_primes[2].delta))
831 * cos(phi)
832 * sin(degrees_to_radians(h_prime[2])));
833
834 let sunrise = add_fraction_of_day(params.day.clone(), r);
835 let transit = add_fraction_of_day(params.day.clone(), t);
836 let sunset = add_fraction_of_day(params.day, s);
837
838 crate::SunriseResult::RegularDay {
839 sunrise,
840 transit,
841 sunset,
842 }
843}
844
845fn calculate_interpolated_alpha_deltas(
847 alpha_deltas: &[AlphaDelta; 3],
848 n: &[f64; 3],
849) -> [AlphaDelta; 3] {
850 let a = limit_if_necessary(alpha_deltas[1].alpha - alpha_deltas[0].alpha);
851 let a_prime = limit_if_necessary(alpha_deltas[1].delta - alpha_deltas[0].delta);
852
853 let b = limit_if_necessary(alpha_deltas[2].alpha - alpha_deltas[1].alpha);
854 let b_prime = limit_if_necessary(alpha_deltas[2].delta - alpha_deltas[1].delta);
855
856 let c = b - a;
857 let c_prime = b_prime - a_prime;
858
859 let mut alpha_delta_primes = [AlphaDelta {
860 alpha: 0.0,
861 delta: 0.0,
862 }; 3];
863 for i in 0..3 {
864 alpha_delta_primes[i].alpha =
865 alpha_deltas[1].alpha + (n[i] * (mul_add(c, n[i], a + b))) / 2.0;
866 alpha_delta_primes[i].delta =
867 alpha_deltas[1].delta + (n[i] * (mul_add(c_prime, n[i], a_prime + b_prime))) / 2.0;
868 }
869 alpha_delta_primes
870}
871
872#[derive(Debug, Clone, Copy)]
873struct AlphaDelta {
874 alpha: f64,
875 delta: f64,
876}
877
878#[derive(Debug, Clone)]
880#[cfg(feature = "chrono")]
881struct FinalTimeParams<Tz: TimeZone> {
882 day: DateTime<Tz>,
883 m_values: [f64; 3],
884 nu_degrees: f64,
885 delta_t: f64,
886 latitude: f64,
887 longitude: f64,
888 elevation_angle: f64,
889 alpha_deltas: [AlphaDelta; 3],
890}
891
892#[cfg(feature = "chrono")]
930#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
931pub fn sunrise_sunset_for_horizon<Tz: TimeZone>(
932 date: DateTime<Tz>,
933 latitude: f64,
934 longitude: f64,
935 delta_t: f64,
936 horizon: Horizon,
937) -> Result<crate::SunriseResult<DateTime<Tz>>> {
938 sunrise_sunset(
939 date,
940 latitude,
941 longitude,
942 delta_t,
943 horizon.elevation_angle(),
944 )
945}
946
947fn calculate_alpha_delta(jme: f64, delta_psi: f64, epsilon_degrees: f64) -> AlphaDelta {
950 let b_terms = calculate_lbr_terms(jme, TERMS_B);
954 let b_degrees = lbr_to_normalized_degrees(jme, &b_terms, TERMS_B.len());
955
956 let r_terms = calculate_lbr_terms(jme, TERMS_R);
958 let r = calculate_lbr_polynomial(jme, &r_terms, TERMS_R.len());
959 assert!(
960 r != 0.0,
961 "Earth radius vector is zero - astronomical impossibility"
962 );
963
964 let l_terms = calculate_lbr_terms(jme, TERMS_L);
966 let l_degrees = lbr_to_normalized_degrees(jme, &l_terms, TERMS_L.len());
967
968 let theta_degrees = normalize_degrees_0_to_360(l_degrees + 180.0);
970
971 let beta_degrees = -b_degrees;
973 let beta = degrees_to_radians(beta_degrees);
974 let epsilon = degrees_to_radians(epsilon_degrees);
975
976 let delta_tau = ABERRATION_CONSTANT / (SECONDS_PER_HOUR * r);
978
979 let lambda_degrees = theta_degrees + delta_psi + delta_tau;
981 let lambda = degrees_to_radians(lambda_degrees);
982
983 let alpha_degrees = calculate_geocentric_sun_right_ascension(beta, epsilon, lambda);
985 let delta_degrees =
986 radians_to_degrees(calculate_geocentric_sun_declination(beta, epsilon, lambda));
987
988 AlphaDelta {
989 alpha: alpha_degrees,
990 delta: delta_degrees,
991 }
992}
993
994fn normalize_to_unit_range(val: f64) -> f64 {
996 let divided = val;
997 let limited = divided - floor(divided);
998 if limited < 0.0 {
999 limited + 1.0
1000 } else {
1001 limited
1002 }
1003}
1004
1005#[cfg(feature = "chrono")]
1006fn truncate_to_day_start<Tz: TimeZone>(datetime: &DateTime<Tz>) -> DateTime<Tz> {
1007 let tz = datetime.timezone();
1008 let local_midnight = datetime
1009 .date_naive()
1010 .and_hms_opt(0, 0, 0)
1011 .expect("midnight is always valid");
1012
1013 match tz.from_local_datetime(&local_midnight) {
1014 LocalResult::Single(dt) => dt,
1015 LocalResult::Ambiguous(earliest, _) => earliest,
1016 LocalResult::None => {
1017 let utc_midnight = datetime
1019 .with_timezone(&Utc)
1020 .date_naive()
1021 .and_hms_opt(0, 0, 0)
1022 .expect("midnight is always valid")
1023 .and_utc();
1024
1025 tz.from_utc_datetime(&utc_midnight.naive_utc())
1026 }
1027 }
1028}
1029
1030#[cfg(feature = "chrono")]
1031#[allow(clippy::needless_pass_by_value)]
1032fn add_fraction_of_day<Tz: TimeZone>(day: DateTime<Tz>, fraction: f64) -> DateTime<Tz> {
1033 const MS_PER_DAY: i32 = 24 * 60 * 60 * 1000; let millis_plus = (f64::from(MS_PER_DAY) * fraction) as i32; let day_start = truncate_to_day_start(&day);
1041
1042 day_start + chrono::Duration::milliseconds(i64::from(millis_plus))
1043}
1044
1045fn limit_if_necessary(val: f64) -> f64 {
1047 if val.abs() > 2.0 {
1048 normalize_to_unit_range(val)
1049 } else {
1050 val
1051 }
1052}
1053
1054fn limit_h_prime(h_prime: f64) -> f64 {
1056 let normalized = h_prime / 360.0;
1057 let limited = 360.0 * (normalized - floor(normalized));
1058
1059 if limited < -180.0 {
1060 limited + 360.0
1061 } else if limited > 180.0 {
1062 limited - 360.0
1063 } else {
1064 limited
1065 }
1066}
1067
1068#[cfg(feature = "chrono")]
1117#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
1118pub fn sunrise_sunset_multiple<Tz, H>(
1119 date: DateTime<Tz>,
1120 latitude: f64,
1121 longitude: f64,
1122 delta_t: f64,
1123 horizons: H,
1124) -> impl Iterator<Item = Result<(Horizon, crate::SunriseResult<DateTime<Tz>>)>>
1125where
1126 Tz: TimeZone,
1127 H: IntoIterator<Item = Horizon>,
1128{
1129 let precomputed = (|| -> Result<_> {
1131 check_coordinates(latitude, longitude)?;
1132
1133 let day_start = truncate_to_day_start(&date);
1134 let (nu_degrees, delta_psi_epsilon, epsilon_degrees) =
1135 calculate_sidereal_time_and_nutation(day_start.clone(), delta_t)?;
1136 let alpha_deltas = calculate_alpha_deltas_for_three_days(
1137 day_start,
1138 delta_t,
1139 delta_psi_epsilon,
1140 epsilon_degrees,
1141 )?;
1142
1143 Ok((nu_degrees, alpha_deltas))
1144 })();
1145
1146 horizons.into_iter().map(move |horizon| match &precomputed {
1147 Ok((nu_degrees, alpha_deltas)) => {
1148 let result = calculate_sunrise_sunset_spa_algorithm_with_precomputed(
1149 date.clone(),
1150 latitude,
1151 longitude,
1152 delta_t,
1153 horizon.elevation_angle(),
1154 *nu_degrees,
1155 *alpha_deltas,
1156 );
1157 Ok((horizon, result))
1158 }
1159 Err(e) => Err(e.clone()),
1160 })
1161}
1162#[cfg(feature = "chrono")]
1163#[allow(clippy::needless_pass_by_value)]
1164fn calculate_sunrise_sunset_spa_algorithm_with_precomputed<Tz: TimeZone>(
1165 day: DateTime<Tz>,
1166 latitude: f64,
1167 longitude: f64,
1168 delta_t: f64,
1169 elevation_angle: f64,
1170 nu_degrees: f64,
1171 alpha_deltas: [AlphaDelta; 3],
1172) -> crate::SunriseResult<DateTime<Tz>> {
1173 let day_start = truncate_to_day_start(&day);
1174
1175 let m0 = (alpha_deltas[1].alpha - longitude - nu_degrees) / 360.0;
1177 let polar_type = check_polar_conditions_type(latitude, elevation_angle, alpha_deltas[1].delta);
1178
1179 let (m_values, _h0_degrees) =
1181 calculate_approximate_times(m0, latitude, elevation_angle, alpha_deltas[1].delta);
1182
1183 let final_result = calculate_final_times(FinalTimeParams {
1185 day: day_start,
1186 m_values,
1187 nu_degrees,
1188 delta_t,
1189 latitude,
1190 longitude,
1191 elevation_angle,
1192 alpha_deltas,
1193 });
1194
1195 match polar_type {
1196 Some(PolarType::AllDay) => {
1197 if let crate::SunriseResult::RegularDay { transit, .. } = final_result {
1198 crate::SunriseResult::AllDay { transit }
1199 } else {
1200 final_result
1201 }
1202 }
1203 Some(PolarType::AllNight) => {
1204 if let crate::SunriseResult::RegularDay { transit, .. } = final_result {
1205 crate::SunriseResult::AllNight { transit }
1206 } else {
1207 final_result
1208 }
1209 }
1210 None => final_result,
1211 }
1212}
1213
1214#[cfg(feature = "chrono")]
1252#[cfg_attr(docsrs, doc(cfg(feature = "chrono")))]
1253#[allow(clippy::needless_pass_by_value)]
1254pub fn spa_time_dependent_parts<Tz: TimeZone>(
1255 datetime: DateTime<Tz>,
1256 delta_t: f64,
1257) -> Result<SpaTimeDependent> {
1258 let jd = JulianDate::from_datetime(&datetime, delta_t)?;
1259 spa_time_dependent_from_julian(jd)
1260}
1261
1262pub fn spa_time_dependent_from_julian(jd: JulianDate) -> Result<SpaTimeDependent> {
1272 let jme = jd.julian_ephemeris_millennium();
1273 let jce = jd.julian_ephemeris_century();
1274
1275 let l_terms = calculate_lbr_terms(jme, TERMS_L);
1277 let l_degrees = lbr_to_normalized_degrees(jme, &l_terms, TERMS_L.len());
1278
1279 let b_terms = calculate_lbr_terms(jme, TERMS_B);
1281 let b_degrees = lbr_to_normalized_degrees(jme, &b_terms, TERMS_B.len());
1282
1283 let r_terms = calculate_lbr_terms(jme, TERMS_R);
1285 let r = calculate_lbr_polynomial(jme, &r_terms, TERMS_R.len());
1286
1287 assert!(
1289 r != 0.0,
1290 "Earth radius vector is zero - astronomical impossibility"
1291 );
1292
1293 let theta_degrees = normalize_degrees_0_to_360(l_degrees + 180.0);
1295 let beta_degrees = -b_degrees;
1297
1298 let x_terms = calculate_nutation_terms(jce);
1300 let delta_psi_epsilon = calculate_delta_psi_epsilon(jce, &x_terms);
1301
1302 let epsilon_degrees =
1304 calculate_true_obliquity_of_ecliptic(&jd, delta_psi_epsilon.delta_epsilon);
1305
1306 let delta_tau = ABERRATION_CONSTANT / (SECONDS_PER_HOUR * r);
1308
1309 let lambda_degrees = theta_degrees + delta_psi_epsilon.delta_psi + delta_tau;
1311
1312 let nu_degrees = calculate_apparent_sidereal_time_at_greenwich(
1314 &jd,
1315 delta_psi_epsilon.delta_psi,
1316 epsilon_degrees,
1317 );
1318
1319 let beta = degrees_to_radians(beta_degrees);
1321 let epsilon = degrees_to_radians(epsilon_degrees);
1322 let lambda = degrees_to_radians(lambda_degrees);
1323 let alpha_degrees = calculate_geocentric_sun_right_ascension(beta, epsilon, lambda);
1324
1325 let delta_degrees =
1327 radians_to_degrees(calculate_geocentric_sun_declination(beta, epsilon, lambda));
1328
1329 Ok(SpaTimeDependent {
1330 r,
1331 nu_degrees,
1332 alpha_degrees,
1333 delta_degrees,
1334 })
1335}
1336
1337pub fn spa_with_time_dependent_parts(
1375 latitude: f64,
1376 longitude: f64,
1377 elevation: f64,
1378 refraction: Option<RefractionCorrection>,
1379 time_dependent: &SpaTimeDependent,
1380) -> Result<SolarPosition> {
1381 check_coordinates(latitude, longitude)?;
1382
1383 let nu_degrees = time_dependent.nu_degrees;
1386
1387 let h_degrees =
1389 normalize_degrees_0_to_360(nu_degrees + longitude - time_dependent.alpha_degrees);
1390 let h = degrees_to_radians(h_degrees);
1391
1392 let xi_degrees = 8.794 / (3600.0 * time_dependent.r);
1394 let xi = degrees_to_radians(xi_degrees);
1395 let phi = degrees_to_radians(latitude);
1396 let delta = degrees_to_radians(time_dependent.delta_degrees);
1397
1398 let u = atan(EARTH_FLATTENING_FACTOR * tan(phi));
1399 let y = mul_add(
1400 EARTH_FLATTENING_FACTOR,
1401 sin(u),
1402 (elevation / EARTH_RADIUS_METERS) * sin(phi),
1403 );
1404 let x = mul_add(elevation / EARTH_RADIUS_METERS, cos(phi), cos(u));
1405
1406 let delta_alpha_prime_degrees = radians_to_degrees(atan2(
1407 -x * sin(xi) * sin(h),
1408 mul_add(x * sin(xi), -cos(h), cos(delta)),
1409 ));
1410
1411 let delta_prime = radians_to_degrees(atan2(
1412 mul_add(y, -sin(xi), sin(delta)) * cos(degrees_to_radians(delta_alpha_prime_degrees)),
1413 mul_add(x * sin(xi), -cos(h), cos(delta)),
1414 ));
1415
1416 let h_prime_degrees = h_degrees - delta_alpha_prime_degrees;
1418
1419 let zenith_angle = radians_to_degrees(acos(mul_add(
1421 sin(degrees_to_radians(latitude)),
1422 sin(degrees_to_radians(delta_prime)),
1423 cos(degrees_to_radians(latitude))
1424 * cos(degrees_to_radians(delta_prime))
1425 * cos(degrees_to_radians(h_prime_degrees)),
1426 )));
1427
1428 let azimuth = normalize_degrees_0_to_360(
1430 180.0
1431 + radians_to_degrees(atan2(
1432 sin(degrees_to_radians(h_prime_degrees)),
1433 cos(degrees_to_radians(h_prime_degrees)) * sin(degrees_to_radians(latitude))
1434 - tan(degrees_to_radians(delta_prime)) * cos(degrees_to_radians(latitude)),
1435 )),
1436 );
1437
1438 let elevation_angle = 90.0 - zenith_angle;
1440 let final_zenith = refraction.map_or(zenith_angle, |correction| {
1441 if elevation_angle > SUNRISE_SUNSET_ANGLE {
1442 let pressure = correction.pressure();
1443 let temperature = correction.temperature();
1444 zenith_angle
1446 - (pressure / 1010.0) * (283.0 / (273.0 + temperature)) * 1.02
1447 / (60.0
1448 * tan(degrees_to_radians(
1449 elevation_angle + 10.3 / (elevation_angle + 5.11),
1450 )))
1451 } else {
1452 zenith_angle
1453 }
1454 });
1455
1456 SolarPosition::new(azimuth, final_zenith)
1457}
1458
1459#[cfg(all(test, feature = "chrono", feature = "std"))]
1460mod tests {
1461 use super::*;
1462 use chrono::{DateTime, FixedOffset, NaiveTime};
1463 use std::collections::HashSet;
1464
1465 #[test]
1466 fn test_spa_basic_functionality() {
1467 let datetime = "2023-06-21T12:00:00Z"
1468 .parse::<DateTime<FixedOffset>>()
1469 .unwrap();
1470
1471 let result = solar_position(
1472 datetime,
1473 37.7749, -122.4194,
1475 0.0,
1476 69.0,
1477 Some(RefractionCorrection::new(1013.25, 15.0).unwrap()),
1478 );
1479
1480 assert!(result.is_ok());
1481 let position = result.unwrap();
1482 assert!(position.azimuth() >= 0.0 && position.azimuth() <= 360.0);
1483 assert!(position.zenith_angle() >= 0.0 && position.zenith_angle() <= 180.0);
1484 }
1485
1486 #[test]
1487 fn test_sunrise_sunset_multiple() {
1488 let datetime = "2023-06-21T12:00:00Z"
1489 .parse::<DateTime<FixedOffset>>()
1490 .unwrap();
1491 let horizons = [
1492 Horizon::SunriseSunset,
1493 Horizon::CivilTwilight,
1494 Horizon::NauticalTwilight,
1495 ];
1496
1497 let results: Result<Vec<_>> = sunrise_sunset_multiple(
1498 datetime,
1499 37.7749, -122.4194, 69.0, horizons.iter().copied(),
1503 )
1504 .collect();
1505
1506 let results = results.unwrap();
1507
1508 assert_eq!(results.len(), 3);
1510
1511 let returned_horizons: HashSet<_> = results.iter().map(|(h, _)| *h).collect();
1513 for expected_horizon in horizons {
1514 assert!(returned_horizons.contains(&expected_horizon));
1515 }
1516
1517 for (horizon, bulk_result) in &results {
1519 let individual_result =
1520 sunrise_sunset_for_horizon(datetime, 37.7749, -122.4194, 69.0, *horizon).unwrap();
1521
1522 match (&individual_result, bulk_result) {
1524 (
1525 crate::SunriseResult::RegularDay {
1526 sunrise: s1,
1527 transit: t1,
1528 sunset: ss1,
1529 },
1530 crate::SunriseResult::RegularDay {
1531 sunrise: s2,
1532 transit: t2,
1533 sunset: ss2,
1534 },
1535 ) => {
1536 assert_eq!(s1, s2);
1537 assert_eq!(t1, t2);
1538 assert_eq!(ss1, ss2);
1539 }
1540 (
1541 crate::SunriseResult::AllDay { transit: t1 },
1542 crate::SunriseResult::AllDay { transit: t2 },
1543 )
1544 | (
1545 crate::SunriseResult::AllNight { transit: t1 },
1546 crate::SunriseResult::AllNight { transit: t2 },
1547 ) => {
1548 assert_eq!(t1, t2);
1549 }
1550 _ => panic!("Bulk and individual results differ in type for {horizon:?}"),
1551 }
1552 }
1553 }
1554
1555 #[test]
1556 fn test_sunrise_sunset_multiple_polar_consistency() {
1557 let datetime = "2023-06-21T12:00:00Z"
1558 .parse::<DateTime<FixedOffset>>()
1559 .unwrap();
1560
1561 let individual = sunrise_sunset_for_horizon(
1562 datetime,
1563 80.0, 0.0,
1565 69.0,
1566 Horizon::SunriseSunset,
1567 )
1568 .unwrap();
1569
1570 let bulk_results: Result<Vec<_>> =
1571 sunrise_sunset_multiple(datetime, 80.0, 0.0, 69.0, [Horizon::SunriseSunset]).collect();
1572
1573 let (_, bulk) = bulk_results.unwrap().into_iter().next().unwrap();
1574
1575 match (bulk, individual) {
1576 (
1577 crate::SunriseResult::AllDay { transit: t1 },
1578 crate::SunriseResult::AllDay { transit: t2 },
1579 )
1580 | (
1581 crate::SunriseResult::AllNight { transit: t1 },
1582 crate::SunriseResult::AllNight { transit: t2 },
1583 ) => assert_eq!(t1, t2),
1584 _ => panic!("expected matching polar-day/night results between bulk and individual"),
1585 }
1586 }
1587
1588 #[test]
1589 fn test_spa_no_refraction() {
1590 let datetime = "2023-06-21T12:00:00Z"
1591 .parse::<DateTime<FixedOffset>>()
1592 .unwrap();
1593
1594 let result = solar_position(datetime, 37.7749, -122.4194, 0.0, 69.0, None);
1595
1596 assert!(result.is_ok());
1597 let position = result.unwrap();
1598 assert!(position.azimuth() >= 0.0 && position.azimuth() <= 360.0);
1599 assert!(position.zenith_angle() >= 0.0 && position.zenith_angle() <= 180.0);
1600 }
1601
1602 #[test]
1603 fn test_spa_coordinate_validation() {
1604 let datetime = "2023-06-21T12:00:00Z"
1605 .parse::<DateTime<FixedOffset>>()
1606 .unwrap();
1607
1608 assert!(solar_position(
1610 datetime,
1611 95.0,
1612 0.0,
1613 0.0,
1614 0.0,
1615 Some(RefractionCorrection::new(1013.25, 15.0).unwrap())
1616 )
1617 .is_err());
1618
1619 assert!(solar_position(
1621 datetime,
1622 0.0,
1623 185.0,
1624 0.0,
1625 0.0,
1626 Some(RefractionCorrection::new(1013.25, 15.0).unwrap())
1627 )
1628 .is_err());
1629 }
1630
1631 #[test]
1632 fn test_sunrise_sunset_basic() {
1633 let date = "2023-06-21T00:00:00Z"
1634 .parse::<DateTime<FixedOffset>>()
1635 .unwrap();
1636
1637 let result = sunrise_sunset(date, 37.7749, -122.4194, 69.0, -0.833);
1638 assert!(result.is_ok());
1639
1640 let result =
1641 sunrise_sunset_for_horizon(date, 37.7749, -122.4194, 69.0, Horizon::SunriseSunset);
1642 assert!(result.is_ok());
1643 }
1644
1645 #[test]
1646 fn test_horizon_enum() {
1647 assert_eq!(Horizon::SunriseSunset.elevation_angle(), -0.83337);
1648 assert_eq!(Horizon::CivilTwilight.elevation_angle(), -6.0);
1649 assert_eq!(Horizon::Custom(-10.5).elevation_angle(), -10.5);
1650 }
1651
1652 #[test]
1653 fn truncate_to_day_start_keeps_local_midnight() {
1654 let tz = FixedOffset::east_opt(2 * 3600).unwrap();
1655 let midnight = tz.with_ymd_and_hms(2024, 3, 30, 0, 0, 0).unwrap();
1656
1657 let truncated = truncate_to_day_start(&midnight);
1658
1659 assert_eq!(truncated.date_naive(), midnight.date_naive());
1660 assert_eq!(truncated.time(), NaiveTime::from_hms_opt(0, 0, 0).unwrap());
1661 }
1662}