use super::instant::TimeScale;
use qtty::Days;
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct JD;
impl TimeScale for JD {
const LABEL: &'static str = "Julian Day:";
#[inline(always)]
fn to_jd_tt(value: Days) -> Days {
value
}
#[inline(always)]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct JDE;
impl TimeScale for JDE {
const LABEL: &'static str = "JDE";
#[inline(always)]
fn to_jd_tt(value: Days) -> Days {
value
}
#[inline(always)]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct MJD;
const MJD_EPOCH: Days = Days::new(2_400_000.5);
impl TimeScale for MJD {
const LABEL: &'static str = "MJD";
#[inline(always)]
fn to_jd_tt(value: Days) -> Days {
value + MJD_EPOCH
}
#[inline(always)]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt - MJD_EPOCH
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct TDB;
#[inline]
pub(crate) fn tdb_minus_tt_days(jd_tt: Days) -> Days {
let t = (jd_tt.value() - 2_451_545.0) / 36_525.0;
let m_e = (357.5291092 + 35999.0502909 * t).to_radians();
let m_j = (246.4512 + 3035.2335 * t).to_radians();
let d = (297.8502042 + 445267.1115168 * t).to_radians();
let om = (125.0445550 - 1934.1362091 * t).to_radians();
let dt_sec = 0.001_657 * (m_e + 0.01671 * m_e.sin()).sin()
+ 0.000_022 * (d - m_e).sin()
+ 0.000_014 * (2.0 * d).sin()
+ 0.000_005 * m_j.sin()
+ 0.000_005 * om.sin();
Days::new(dt_sec / 86_400.0)
}
impl TimeScale for TDB {
const LABEL: &'static str = "TDB";
#[inline]
fn to_jd_tt(tdb_value: Days) -> Days {
tdb_value - tdb_minus_tt_days(tdb_value)
}
#[inline]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt + tdb_minus_tt_days(jd_tt)
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct TT;
impl TimeScale for TT {
const LABEL: &'static str = "TT";
#[inline(always)]
fn to_jd_tt(value: Days) -> Days {
value }
#[inline(always)]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct TAI;
const TT_MINUS_TAI: Days = Days::new(32.184 / 86_400.0);
impl TimeScale for TAI {
const LABEL: &'static str = "TAI";
#[inline(always)]
fn to_jd_tt(value: Days) -> Days {
value + TT_MINUS_TAI
}
#[inline(always)]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt - TT_MINUS_TAI
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct TCG;
const L_G: f64 = 6.969_290_134e-10;
const TCG_EPOCH_T0: f64 = 2_443_144.500_372_5;
impl TimeScale for TCG {
const LABEL: &'static str = "TCG";
#[inline]
fn to_jd_tt(tcg_value: Days) -> Days {
let jd_tcg = tcg_value.value();
Days::new(jd_tcg - L_G * (jd_tcg - TCG_EPOCH_T0))
}
#[inline]
fn from_jd_tt(jd_tt: Days) -> Days {
let tt = jd_tt.value();
Days::new(tt + L_G * (tt - TCG_EPOCH_T0) / (1.0 - L_G))
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct TCB;
const L_B: f64 = 1.550_519_768e-8;
impl TimeScale for TCB {
const LABEL: &'static str = "TCB";
#[inline]
fn to_jd_tt(tcb_value: Days) -> Days {
let jd_tcb = tcb_value.value();
Days::new(jd_tcb - L_B * (jd_tcb - TCG_EPOCH_T0))
}
#[inline]
fn from_jd_tt(jd_tt: Days) -> Days {
let tt = jd_tt.value();
Days::new(tt + L_B * (tt - TCG_EPOCH_T0))
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct GPS;
const GPS_EPOCH_JD_TT: Days = Days::new(2_444_244.5 + 51.184 / 86_400.0);
impl TimeScale for GPS {
const LABEL: &'static str = "GPS";
#[inline(always)]
fn to_jd_tt(value: Days) -> Days {
value + GPS_EPOCH_JD_TT
}
#[inline(always)]
fn from_jd_tt(jd_tt: Days) -> Days {
jd_tt - GPS_EPOCH_JD_TT
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct UnixTime;
const UNIX_EPOCH_JD: Days = Days::new(2_440_587.5);
const LEAP_SECONDS: [(f64, f64); 28] = [
(2_441_317.5, 10.0), (2_441_499.5, 11.0), (2_441_683.5, 12.0), (2_442_048.5, 13.0), (2_442_413.5, 14.0), (2_442_778.5, 15.0), (2_443_144.5, 16.0), (2_443_509.5, 17.0), (2_443_874.5, 18.0), (2_444_239.5, 19.0), (2_444_786.5, 20.0), (2_445_151.5, 21.0), (2_445_516.5, 22.0), (2_446_247.5, 23.0), (2_447_161.5, 24.0), (2_447_892.5, 25.0), (2_448_257.5, 26.0), (2_448_804.5, 27.0), (2_449_169.5, 28.0), (2_449_534.5, 29.0), (2_450_083.5, 30.0), (2_450_630.5, 31.0), (2_451_179.5, 32.0), (2_453_736.5, 33.0), (2_454_832.5, 34.0), (2_456_109.5, 35.0), (2_457_204.5, 36.0), (2_457_754.5, 37.0), ];
#[inline]
pub fn tai_minus_utc(jd_utc: f64) -> f64 {
let mut lo = 0usize;
let mut hi = LEAP_SECONDS.len();
while lo < hi {
let mid = lo + (hi - lo) / 2;
if LEAP_SECONDS[mid].0 <= jd_utc {
lo = mid + 1;
} else {
hi = mid;
}
}
if lo == 0 {
10.0
} else {
LEAP_SECONDS[lo - 1].1
}
}
const TT_MINUS_TAI_SECS: f64 = 32.184;
impl TimeScale for UnixTime {
const LABEL: &'static str = "Unix";
#[inline]
fn to_jd_tt(value: Days) -> Days {
let jd_utc = value.value() + UNIX_EPOCH_JD.value();
let ls = tai_minus_utc(jd_utc);
Days::new(jd_utc + (ls + TT_MINUS_TAI_SECS) / 86_400.0)
}
#[inline]
fn from_jd_tt(jd_tt: Days) -> Days {
let approx_utc = jd_tt.value() - (37.0 + TT_MINUS_TAI_SECS) / 86_400.0;
let ls = tai_minus_utc(approx_utc);
let jd_utc = jd_tt.value() - (ls + TT_MINUS_TAI_SECS) / 86_400.0;
Days::new(jd_utc - UNIX_EPOCH_JD.value())
}
}
#[derive(Debug, Copy, Clone, PartialEq, PartialOrd)]
pub struct UT;
impl TimeScale for UT {
const LABEL: &'static str = "UT";
#[inline]
fn to_jd_tt(ut_value: Days) -> Days {
let jd_ut = super::instant::Time::<JD>::from_days(ut_value);
let dt_secs = super::delta_t::delta_t_seconds_from_ut(jd_ut);
ut_value + dt_secs.to::<qtty::Day>()
}
#[inline]
fn from_jd_tt(jd_tt: Days) -> Days {
let mut ut = jd_tt;
for _ in 0..3 {
let jd_ut = super::instant::Time::<JD>::from_days(ut);
let dt_days = super::delta_t::delta_t_seconds_from_ut(jd_ut).to::<qtty::Day>();
ut = jd_tt - dt_days;
}
ut
}
}
macro_rules! impl_time_conversions {
($single:ty) => {};
($first:ty, $($rest:ty),+ $(,)?) => {
$(
impl From<super::instant::Time<$first>> for super::instant::Time<$rest> {
#[inline]
fn from(t: super::instant::Time<$first>) -> Self {
t.to::<$rest>()
}
}
impl From<super::instant::Time<$rest>> for super::instant::Time<$first> {
#[inline]
fn from(t: super::instant::Time<$rest>) -> Self {
t.to::<$first>()
}
}
)+
impl_time_conversions!($($rest),+);
};
}
impl_time_conversions!(JD, JDE, MJD, TDB, TT, TAI, TCG, TCB, GPS, UnixTime, UT);
#[cfg(test)]
mod tests {
use super::super::instant::Time;
use super::*;
use qtty::{Day, Second, Seconds};
#[test]
fn jd_mjd_roundtrip() {
let jd = Time::<JD>::new(2_451_545.0);
let mjd: Time<MJD> = jd.to::<MJD>();
assert!((mjd.quantity() - Days::new(51_544.5)).abs() < Days::new(1e-10));
let back: Time<JD> = mjd.to::<JD>();
assert!((back.quantity() - Days::new(2_451_545.0)).abs() < Days::new(1e-10));
}
#[test]
fn jd_mjd_from_into() {
let jd = Time::<JD>::new(2_451_545.0);
let mjd: Time<MJD> = jd.into();
assert!((mjd.quantity() - Days::new(51_544.5)).abs() < Days::new(1e-10));
let back: Time<JD> = Time::from(mjd);
assert!((back.quantity() - Days::new(2_451_545.0)).abs() < Days::new(1e-10));
}
#[test]
fn tai_tt_offset() {
let tai = Time::<TAI>::new(2_451_545.0);
let tt: Time<TT> = tai.to::<TT>();
let expected_offset = Seconds::new(32.184).to::<Day>();
assert!((tt.quantity() - (tai.quantity() + expected_offset)).abs() < Days::new(1e-15));
}
#[test]
fn gps_epoch_roundtrip() {
let gps_zero = Time::<GPS>::new(0.0);
let jd: Time<JD> = gps_zero.to::<JD>();
let expected = Days::new(2_444_244.5) + Seconds::new(51.184).to::<Day>();
assert!((jd.quantity() - expected).abs() < Days::new(1e-12));
}
#[test]
fn unix_epoch_roundtrip() {
let unix_zero = Time::<UnixTime>::new(0.0);
let jd: Time<JD> = unix_zero.to::<JD>();
let expected = Days::new(2_440_587.5) + Seconds::new(42.184).to::<Day>();
assert!(
(jd.quantity() - expected).abs() < Days::new(1e-10),
"Unix epoch JD(TT) = {}, expected {}",
jd.quantity(),
expected
);
}
#[test]
fn unix_2020_leap_seconds() {
let unix_2020 = Time::<UnixTime>::new(18262.0);
let jd: Time<JD> = unix_2020.to::<JD>();
let expected = Days::new(2_458_849.5) + Seconds::new(69.184).to::<Day>();
assert!(
(jd.quantity() - expected).abs() < Days::new(1e-10),
"2020 Unix JD(TT) = {}, expected {}",
jd.quantity(),
expected
);
}
#[test]
fn tdb_tt_offset() {
let jd = Time::<JD>::new(2_451_545.0);
let tdb: Time<TDB> = jd.to::<TDB>();
let offset_secs = (tdb.quantity() - jd.quantity()).to::<Second>();
assert!(
offset_secs.abs() < Seconds::new(0.002),
"TDB−TT offset at J2000 = {} s, expected < 2 ms",
offset_secs
);
}
#[test]
fn tdb_tt_roundtrip() {
let jd = Time::<JD>::new(2_451_545.0);
let tdb: Time<TDB> = jd.to::<TDB>();
let back: Time<JD> = tdb.to::<JD>();
assert!(
(back.quantity() - jd.quantity()).abs() < Days::new(1e-14),
"TDB→TT roundtrip error: {} days",
(back.quantity() - jd.quantity()).abs()
);
}
#[test]
fn tcg_tt_offset_at_j2000() {
let tt = Time::<TT>::new(2_451_545.0);
let tcg: Time<TCG> = tt.to::<TCG>();
let offset_days = tcg.quantity() - tt.quantity();
let offset_secs = offset_days.to::<Second>();
assert!(
(offset_secs - Seconds::new(0.506)).abs() < Seconds::new(0.01),
"TCG−TT offset at J2000 = {} s, expected ~0.506 s",
offset_secs
);
}
#[test]
fn tcg_tt_roundtrip() {
let tt = Time::<TT>::new(2_451_545.0);
let tcg: Time<TCG> = tt.to::<TCG>();
let back: Time<TT> = tcg.to::<TT>();
assert!(
(back.quantity() - tt.quantity()).abs() < Days::new(1e-12),
"TCG→TT roundtrip error: {} days",
(back.quantity() - tt.quantity()).abs()
);
}
#[test]
fn tcb_tdb_offset_at_j2000() {
let tt = Time::<TT>::new(2_451_545.0);
let tcb: Time<TCB> = tt.to::<TCB>();
let offset_days = tcb.quantity() - tt.quantity();
let offset_secs = offset_days.to::<Second>();
assert!(
(offset_secs - Seconds::new(11.25)).abs() < Seconds::new(0.5),
"TCB−TT offset at J2000 = {} s, expected ~11.25 s",
offset_secs
);
}
#[test]
fn tcb_tt_roundtrip() {
let tt = Time::<TT>::new(2_458_000.0);
let tcb: Time<TCB> = tt.to::<TCB>();
let back: Time<TT> = tcb.to::<TT>();
assert!(
(back.quantity() - tt.quantity()).abs() < Days::new(1e-10),
"TCB→TT roundtrip error: {} days",
(back.quantity() - tt.quantity()).abs()
);
}
#[test]
fn ut_to_jd_applies_delta_t() {
let ut = Time::<UT>::new(2_451_545.0);
let jd: Time<JD> = ut.to::<JD>();
let offset_secs = (jd.quantity() - ut.quantity()).to::<Second>();
assert!(
(offset_secs - Seconds::new(63.83)).abs() < Seconds::new(1.0),
"UT→JD offset = {} s, expected ~63.83 s",
offset_secs
);
}
#[test]
fn ut_jd_roundtrip() {
let jd = Time::<JD>::new(2_451_545.0);
let ut: Time<UT> = jd.to::<UT>();
let back: Time<JD> = ut.to::<JD>();
assert!(
(back.quantity() - jd.quantity()).abs() < Days::new(1e-12),
"roundtrip error: {} days",
(back.quantity() - jd.quantity()).abs()
);
}
#[test]
fn ut_from_into() {
let ut = Time::<UT>::new(2_451_545.0);
let jd: Time<JD> = ut.into();
let back: Time<UT> = jd.into();
assert!((back.quantity() - ut.quantity()).abs() < Days::new(1e-12));
}
#[test]
fn jde_roundtrip() {
let jde = Time::<JDE>::new(2_451_545.0);
let jd: Time<JD> = jde.to::<JD>();
assert!((jd.quantity() - Days::new(2_451_545.0)).abs() < Days::new(1e-10));
let back: Time<JDE> = jd.to::<JDE>();
assert!((back.quantity() - jde.quantity()).abs() < Days::new(1e-10));
}
#[test]
fn tt_to_tai() {
let tt = Time::<TT>::new(2_451_545.0);
let tai: Time<TAI> = tt.to::<TAI>();
let back: Time<TT> = tai.to::<TT>();
assert!(
(back.quantity() - tt.quantity()).abs() < Days::new(1e-15),
"TT → TAI → TT roundtrip error: {}",
(back.quantity() - tt.quantity()).abs()
);
}
#[test]
fn gps_from_jd() {
let gps_zero = Time::<GPS>::new(0.0);
let jd: Time<JD> = gps_zero.to::<JD>();
let back: Time<GPS> = jd.to::<GPS>();
assert!((back.quantity() - gps_zero.quantity()).abs() < Days::new(1e-12));
}
#[test]
fn unix_from_jd() {
let unix_2020 = Time::<UnixTime>::new(18262.0); let jd: Time<JD> = unix_2020.to::<JD>();
let back: Time<UnixTime> = jd.to::<UnixTime>();
assert!(
(back.quantity() - unix_2020.quantity()).abs() < Days::new(1e-10),
"UnixTime roundtrip error: {} days",
(back.quantity() - unix_2020.quantity()).abs()
);
}
#[test]
fn tai_minus_utc_pre_1972_returns_10() {
assert_eq!(tai_minus_utc(2_400_000.0), 10.0);
}
}