SuperNOVAS v1.5
The NOVAS C library, made better
Loading...
Searching...
No Matches
earth.c File Reference

Functions

double era (double jd_ut1_high, double jd_ut1_low)
 
short geo_posvel (double jd_tt, double ut1_to_tt, enum novas_accuracy accuracy, const observer *restrict obs, double *restrict pos, double *restrict vel)
 
int limb_angle (const double *pos_src, const double *pos_obs, double *restrict limb_ang, double *restrict nadir_ang)
 
int novas_diurnal_eop (double gmst, const novas_delaunay_args *restrict delaunay, double *restrict dxp, double *restrict dyp, double *restrict dut1)
 
int novas_diurnal_eop_at_time (const novas_timespec *restrict time, double *restrict dxp, double *restrict dyp, double *restrict dut1)
 
int novas_diurnal_libration (double gmst, const novas_delaunay_args *restrict delaunay, double *restrict dxp, double *restrict dyp, double *restrict dut1)
 
int novas_diurnal_ocean_tides (double gmst, const novas_delaunay_args *restrict delaunay, double *restrict dxp, double *restrict dyp, double *restrict dut1)
 
double novas_gast (double jd_ut1, double ut1_to_tt, enum novas_accuracy accuracy)
 
double novas_gmst (double jd_ut1, double ut1_to_tt)
 
short sidereal_time (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_equinox_type gst_type, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, double *restrict gst)
 
int terra (const on_surface *restrict location, double lst, double *restrict pos, double *restrict vel)
 
int wobble (double jd_tt, enum novas_wobble_direction direction, double xp, double yp, const double *in, double *out)
 

Detailed Description

Date
Created on Mar 6, 2025
Author
G. Kaplan and Attila Kovacs

Various finctions relating to Earth position and orientation

Function Documentation

◆ era()

double era ( double  jd_ut1_high,
double  jd_ut1_low 
)

Returns the value of the Earth Rotation Angle (θ) for a given UT1 Julian date. The expression used is taken from the note to IAU Resolution B1.8 of 2000. The input Julian date cane be split into an into high and low order parts (e.g. integer and fractional parts) for improved accuracy, or else one of the components (e.g. the low part) can be set to zero if no split is desired.

The algorithm used here is equivalent to the canonical θ = 0.7790572732640 + 1.00273781191135448 * t, where t is the time in days from J2000 (t = jd_high + jd_low - JD_J2000), but it avoids many two-PI 'wraps' that decrease precision (adopted from SOFA Fortran routine iau_era00; see also expression at top of page 35 of IERS Conventions (1996)).

REFERENCES:

  1. IERS Conventions, Chapter 5, Eq. 5.15
  2. IAU Resolution B1.8, adopted at the 2000 IAU General Assembly, Manchester, UK.
  3. Kaplan, G. (2005), US Naval Observatory Circular 179.
Parameters
jd_ut1_high[day] High-order part of UT1 Julian date.
jd_ut1_low[day] Low-order part of UT1 Julian date.
Returns
[deg] The Earth Rotation Angle (theta) in degrees [0:360].
See also
sidereal_time()
cirs_to_itrs()
itrs_to_cirs()

◆ geo_posvel()

short geo_posvel ( double  jd_tt,
double  ut1_to_tt,
enum novas_accuracy  accuracy,
const observer *restrict  obs,
double *restrict  pos,
double *restrict  vel 
)

Computes the geocentric GCRS position and velocity of an observer.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date.
ut1_to_tt[s] TT - UT1 time difference in seconds
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
obsObserver location
[out]pos[AU] Position 3-vector of observer, with respect to origin at geocenter, referred to GCRS axes, components in AU. (It may be NULL if not required.)
[out]vel[AU/day] Velocity 3-vector of observer, with respect to origin at geocenter, referred to GCRS axes, components in AU/day. (It must be distinct from the pos output vector, and may be NULL if not required)
Returns
0 if successful, -1 if the 'obs' is NULL or the two output vectors are the same, or else 1 if 'accuracy' is invalid, or 2 if 'obserrver->where' is invalid.
See also
place()
make_observer()
get_ut1_to_tt()

References ephemeris(), geo_posvel(), NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_FULL_ACCURACY, novas_gast(), NOVAS_OBSERVER_AT_GEOCENTER, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_REDUCED_ACCURACY, NOVAS_SOLAR_SYSTEM_OBSERVER, terra(), tod_to_gcrs(), tt2tdb(), and observer::where.

◆ limb_angle()

int limb_angle ( const double *  pos_src,
const double *  pos_obs,
double *restrict  limb_ang,
double *restrict  nadir_ang 
)

Determines the angle of an object above or below the Earth's limb (horizon). The geometric limb is computed, assuming the Earth to be an airless sphere (no refraction or oblateness is included). The observer can be on or above the Earth. For an observer on the surface of the Earth, this function returns the approximate unrefracted elevation.

Parameters
pos_src[AU] Position 3-vector of observed object, with respect to origin at geocenter, components in AU.
pos_obs[AU] Position 3-vector of observer, with respect to origin at geocenter, components in AU.
[out]limb_ang[deg] Angle of observed object above (+) or below (-) limb in degrees, or NAN if reurning with an error. It may be NULL if not required.
[out]nadir_angNadir angle of observed object as a fraction of apparent radius of limb: lt;1.0 if below the limb; 1.0 on the limb; or >1.0 if above the limb. Returns NAN in case of an error return. It may be NULL if not required.
Returns
0 if successful, or -1 if either of the input vectors is NULL or if either input position is a null vector (at the geocenter).
See also
place()

References M_PI, and novas_vlen().

◆ novas_diurnal_eop()

int novas_diurnal_eop ( double  gmst,
const novas_delaunay_args *restrict  delaunay,
double *restrict  dxp,
double *restrict  dyp,
double *restrict  dut1 
)

Calculate corrections to the Earth orientation parameters (EOP) due to short term (diurnal and semidiurnal) librations and the ocean tides. See Chapters 5 and 8 of the IERS Conventions. This one is faster than novas_diurnal_eop_at_time() if the GMST and/or fundamental arguments are readily avaialbe, e.g. from a previous calculation.

These corrections are typically at the sub-milliarcsecond level in polar offsets x and y and sub-millisecond level in UT1. Thus, these diurnal and semi-diurnal variations are important for the highest precision astrometry only, when converting between ITRS and TIRS frames.

REFERENCES:

  1. IERS Conventions Chapter 5, https://iers-conventions.obspm.fr/content/chapter5/icc5.pdf
  2. IERS Conventions Chapter 8, https://iers-conventions.obspm.fr/content/chapter8/icc8.pdf
Parameters
gmst[h] Greenwich Mean Sidereal Time of observation, e.g. as obtained by novas_gmst().
delaunay[rad] Delaunay arguments, e.g. as returned by fund_args().
[out]dxp[arcsec] x-pole correction for libration and ocean tides, or NULL if not required.
[out]dyp[arcsec] y-pole corrections for libration and ocean tides, or NULL if not required.
[out]dut1[s] UT1 correction for libration and ocean tides, or NULL if not required.
Returns
0 if successful, or else -1 if the delaunay arguments are NULL (errno set to EINVAL).
Since
1.5
Author
Attila Kovacs
See also
novas_diurnal_eop_at_time()
novas_diurnal_librations()
novas_diurnal_ocean_tides()

◆ novas_diurnal_eop_at_time()

int novas_diurnal_eop_at_time ( const novas_timespec *restrict  time,
double *restrict  dxp,
double *restrict  dyp,
double *restrict  dut1 
)

Calculate corrections to the Earth orientation parameters (EOP) due to short term (diurnal and semidiurnal) librations and the ocean tides at a given astromtric time. See Chapters 5 and 8 of the IERS Conventions.

REFERENCES:

  1. IERS Conventions Chapter 5, https://iers-conventions.obspm.fr/content/chapter5/icc5.pdf
  2. IERS Conventions Chapter 8, https://iers-conventions.obspm.fr/content/chapter8/icc8.pdf
Parameters
timeAstrometric time specification
dxp[arcsec] x-pole correction for libration and ocean tides, or NULL if not required.
dyp[arcsec] y-pole corrections for libration and ocean tides, or NULL if not required.
dut1[s] UT1 correction for libration and ocean tides, or NULL if not required.
Returns
0 if successful, or else -1 if the delaunay arguments are NULL (errno set to EINVAL).
Since
1.5
Author
Attila Kovacs
See also
novas_diurnal_eop()

References fund_args(), novas_diurnal_eop(), NOVAS_JD_J2000, NOVAS_REDUCED_ACCURACY, and novas_time_gst().

◆ novas_diurnal_libration()

int novas_diurnal_libration ( double  gmst,
const novas_delaunay_args *restrict  delaunay,
double *restrict  dxp,
double *restrict  dyp,
double *restrict  dut1 
)

Calculate diurnal and semi-diurnal libration corrections to the Earth orientation parameters (EOP) for the non-rigid Earth. Only terms with periods less than 2 days are considered, since the longer period terms are included in the IAU2000A nutation model. See Chapter 5 of the IERS conventions. The typical amplitude of librations is tens of micro-arcseconds (μas) in polar offsets x and y and tens of micro-seconds in UT1.

Normally, you would need the combined effect from librations and ocean tides together to apply diurnal corrections to the Earth orientation parameters (EOP) published by IERS. For that, novas_diurnal_eop() or novas_diurnal_eop_at_time() is more directly suited.

REFERENCES:

  1. IERS Conventions Chapter 5, Tablea 5.1a and 5.1b, see https://iers-conventions.obspm.fr/content/chapter5/icc5.pdf
Parameters
gmst[h] Greenwich Mean Sidereal Time of observation, e.g. as obtained by novas_gmst().
delaunay[rad] Delaunay arguments, e.g. as returned by fund_args().
[out]dxp[arcsec] x-pole correction due to librations, or NULL if not required.
[out]dyp[arcsec] y-pole correction due to librations, or NULL if not required.
[out]dut1[s] UT1 correction due to librations, or NULL if not required.
Returns
0 if successful, or else -1 if the delaunay arguments are NULL (errno set to EINVAL).
Since
1.5
Author
Attila Kovacs
See also
novas_diurnal_eop()
novas_diurnal_ocean_tides()

◆ novas_diurnal_ocean_tides()

int novas_diurnal_ocean_tides ( double  gmst,
const novas_delaunay_args *restrict  delaunay,
double *restrict  dxp,
double *restrict  dyp,
double *restrict  dut1 
)

Calculate corrections to the Earth orientation parameters (EOP) due to the ocean tides. See Chapter 8 of the IERS Conventions. Ocean tides manifest at the sub-milliarcsecond level in polar offsets x and y and sub-millisecond level in UT1.

Normally, you would need the combined effect from librations and ocean tides together to apply diurnal corrections to the Earth orientation parameters (EOP) published by IERS. For that, novas_diurnal_eop() or novas_diurnal_eop_at_time() is more directly suited.

REFERENCES:

  1. IERS Conventions Chapter 8, Tables 8.2a, 8.2b, 8.3a, and 8.3b, see https://iers-conventions.obspm.fr/content/chapter8/icc8.pdf
Parameters
gmst[h] Greenwich Mean Sidereal Time of observation, e.g. as obtained by novas_gmst().
delaunay[rad] Delaunay arguments, e.g. as returned by fund_args().
[out]dxp[arcsec] x-pole correction due to ocean tides, or NULL if not required.
[out]dyp[arcsec] y-pole correction due to ocean tides, or NULL if not required.
[out]dut1[s] UT1 correction due to ocean tides, or NULL if not required.
Returns
0 if successful, or else -1 if the delaunay arguments are NULL (errno set to EINVAL).
Since
1.5
Author
Attila Kovacs
See also
novas_diurnal_eop()
novas_diurnal_librations()

◆ novas_gast()

double novas_gast ( double  jd_ut1,
double  ut1_to_tt,
enum novas_accuracy  accuracy 
)

Returns the Greenwich Apparent Sidereal Time (GAST) for a given UT1 date.

Parameters
jd_ut1[day] UT1-based Julian Date.
ut1_to_tt[s] UT1 - UTC time difference.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
Returns
[h] The Greenwich Apparent Sidereal Time (GAST) in the 0-24 range.

REFERENCES:

  1. Capitaine, N. et al. (2003), Astronomy and Astrophysics 412, 567-586, eq. (42).
  2. Kaplan, G. (2005), US Naval Observatory Circular 179.
  3. IERS Conventions Chapter 5, Eq. 5.30 and Table 5.2e, see at https://iers-conventions.obspm.fr/content/chapter5/icc5.pdf and https://iers-conventions.obspm.fr/content/chapter5/additional_info/tab5.2e.txt
Since
1.5
Author
Attila Kovacs
See also
novas_gmst()

References e_tilt(), and novas_gmst().

◆ novas_gmst()

double novas_gmst ( double  jd_ut1,
double  ut1_to_tt 
)

Returns the Greenwich Mean Sidereal Time (GMST) for a given UT1 date, using eq. 42 from Capitaine et al. (2003).

REFERENCES:

  1. Capitaine, N. et al. (2003), Astronomy and Astrophysics 412, 567-586, eq. (42).
  2. IERS Conventions Chapter 5, Eq. 5.32, see at https://iers-conventions.obspm.fr/content/chapter5/icc5.pdf
Parameters
jd_ut1[day] UT1-based Julian Date.
ut1_to_tt[s] UT1 - UTC time difference.
Returns
[h] The Greenwich Mean Sidereal Time (GMST) in the 0-24 range.
Since
1.5
Author
Attila Kovacs
See also
novas_gast()

References era(), and novas_gmst_prec().

◆ sidereal_time()

short sidereal_time ( double  jd_ut1_high,
double  jd_ut1_low,
double  ut1_to_tt,
enum novas_equinox_type  gst_type,
enum novas_earth_rotation_measure  erot,
enum novas_accuracy  accuracy,
double *restrict  gst 
)
Deprecated:
Use novas_gmst() or novas_gast() instead to get the same results simpler.

Computes the Greenwich sidereal time, either mean or apparent, at the specified Julian date. The Julian date can be broken into two parts if convenient, but for the highest precision, set 'jd_high' to be the integral part of the Julian date, and set 'jd_low' to be the fractional part.

NOTES:

  1. Contains fix for known sidereal time units bug.
  2. As of version 1.5.0, the function uses Eq 42 of Capitaine et al. (2003) exclusively. In prior versions, including NOVAS C the same method was provided twice, with the only difference being an unnecessary change of coordinate basis to GCRS. This latter, mor roundabout way of getting the same answer has been eliminated.<.li>

REFERENCES:

  1. Capitaine, N. et al. (2003), Astronomy and Astrophysics 412, 567-586, eq. (42).
  2. Kaplan, G. (2005), US Naval Observatory Circular 179.
Parameters
jd_ut1_high[day] High-order part of UT1 Julian date.
jd_ut1_low[day] Low-order part of UT1 Julian date. (You can leave it at zero if 'jd_high' specified the date with sufficient precision)
ut1_to_tt[s] TT - UT1 Time difference in seconds
gst_typeNOVAS_MEAN_EQUINOX (0) or NOVAS_TRUE_EQUINOX (1), depending on whether wanting mean or apparent GST, respectively.
erotUnused as of 1.5.0.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]gst[h] Greenwich (mean or apparent) sidereal time, in hours [0:24]. (In case the returned error code is >1 the gst value will be set to NAN.)
Returns
0 if successful, or -1 if the 'gst' argument is NULL, 1 if 'accuracy' is
See also
novas_time_gst()
era()
tod_to_itrs()
itrs_to_tod()
cel_pole()
get_ut1_to_tt()

References NOVAS_FULL_ACCURACY, novas_gast(), novas_gmst(), NOVAS_REDUCED_ACCURACY, and NOVAS_TRUE_EQUINOX.

◆ terra()

int terra ( const on_surface *restrict  location,
double  lst,
double *restrict  pos,
double *restrict  vel 
)

Computes the position and velocity vectors of a terrestrial observer with respect to the center of the Earth, based on the GRS80 reference ellipsoid, consistent with the IERS conventions.

This function ignores polar motion, unless the observer's longitude and latitude have been corrected for it, and variation in the length of day (angular velocity of earth).

The true equator and equinox of date do not form an inertial system. Therefore, with respect to an inertial system, the very small velocity component (several meters/day) due to the precession and nutation of the Earth's axis is not accounted for here.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
Parameters
locationLocation of observer in Earth's rotating frame
lst[h] Local apparent sidereal time at reference meridian in hours.
[out]pos[AU] Position vector of observer with respect to center of Earth, equatorial rectangular coordinates, referred to true equator and equinox of date, components in AU. If reference meridian is Greenwich and 'lst' = 0, 'pos' is effectively referred to equator and Greenwich. (It may be NULL if no position data is required).
[out]vel[AU/day] Velocity vector of observer with respect to center of Earth, equatorial rectangular coordinates, referred to true equator and equinox of date, components in AU/day. (It must be distinct from the pos output vector, and may be NULL if no velocity data is required).
Returns
0 if successful, or -1 if location is NULL or if the pos and vel output arguments are identical pointers.
See also
make_on_surface()
geo_posvel()
sidereal_time()

References NOVAS_KM.

◆ wobble()

int wobble ( double  jd_tt,
enum novas_wobble_direction  direction,
double  xp,
double  yp,
const double *  in,
double *  out 
)

Corrects a vector in the ITRS (rotating Earth-fixed system) for polar motion, and also corrects the longitude origin (by a tiny amount) to the Terrestrial Intermediate Origin (TIO). The ITRS vector is thereby transformed to the Terrestrial Intermediate Reference System (TIRS), or equivalently the Pseudo Earth Fixed (PEF), based on the true (rotational) equator and TIO; or vice versa. Because the true equator is the plane orthogonal to the direction of the Celestial Intermediate Pole (CIP), the components of the output vector are referred to z and x axes toward the CIP and TIO, respectively.

NOTES:

  1. You may obtain Earth orientation parameters from the IERS site and Bulletins. For sub-mas accuracy you should augment the (interpolated) published values for diurnal and semi-diurnal librations and the effect of ocean tides, e.g. by using novas_diurnal_eop() or novas_diurnal_eop_at_time().
  2. Generally, this function should not be called if global pole offsets were set via cel_pole() and then used via place() or one of its variants to calculate Earth orientation corrected (TOD or CIRS) apparent coordinates. In such cases, calling wobble() would apply duplicate corrections. It is generally best to forgo using cel_pole() going forward, and instead apply Earth orinetation corrections with wobble() only when converting vectors between the Earth-fixed ITRS and TIRS frames.

REFERENCES:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Lambert & Bizouard (2002), Astronomy and Astrophysics 394, 317-321.
  3. IERS Conventions Chapter 5, https://iers-conventions.obspm.fr/content/chapter5/icc5.pdf
  4. IERS Conventions Chapter 8, https://iers-conventions.obspm.fr/content/chapter8/icc8.pdf
Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date.
directionWOBBLE_ITRS_TO_TIRS (0) or WOBBLE_TIRS_TO_ITRS (1) to include corrections for the TIO's longitude (IAU 2006 method); or else WOBBLE_ITRS_TO_PEF (2) or WOBBLE_PEF_TO_ITRS (3) to correct for dx, dy but not for the TIO's longitude (old, pre IAU 2006 method). Negative values default to WOBBLE_TIRS_TO_ITRS.
xp[arcsec] Conventionally-defined X coordinate of Celestial Intermediate Pole with respect to ITRS pole, in arcseconds. As measured or else the interpolated published value from IERS, possibly augmented for diurnal variations caused by librations ocean tides if precision below the milliarcsecond level is required.
yp[arcsec] Conventionally-defined Y coordinate of Celestial Intermediate Pole with respect to ITRS pole, in arcseconds. As measured or else the interpolated published value from IERS, possibly augmented for diurnal variations caused by librations ocean tides if precision below the milliarcsecond level is required.
inInput position vector, geocentric equatorial rectangular coordinates, in the original system defined by 'direction'
[out]outOutput Position vector, geocentric equatorial rectangular coordinates, in the final system defined by 'direction'. It can be the same vector as the input.
Returns
0 if successful, or -1 if the direction is invalid output vector argument is NULL.
See also
novas_diurnal_eop()
novas_diurnal_eop_at_time()

References NOVAS_WOBBLE_DIRECTIONS, WOBBLE_ITRS_TO_PEF, WOBBLE_ITRS_TO_TIRS, and WOBBLE_TIRS_TO_ITRS.