ECI & Sidereal Time#
This is a description of the sktimeutils.eci module. This module sprovides the capability to convert from Earth Centered Inertial, ECI, to geocentric coordinates, GEO. The ECI coordinates are assumed to be specified with the true equator and mean equinox. The conversion from ECI to GEO is a rotation around the Z axis based upon the Greenwich Mean Sidereal Time, see function gmst.
eciteme_to_itrf#
- eciteme_to_itrf(eciv: ndarray, utc: datetime64 | datetime | float, do_polar_motion: bool = False) ndarray[source]#
Converts ECI TEME vectors to ITRF/GEO geocentric vectors. Calculates greenwich mean sidereal time using IAU1982 and rotates the vector around the Earth’s Z axis to get Psuedo Earth Fixed (PEF). This is usually good enough for most application and by default does not account for Polar Motion, which is generally quite small, around a few meters. Applications that require higher precision can request that polar motion be included but it does slow the code down.
This code has been tested against the skyfield package using method skyfield.sgp4lib.TEME_to_ITRF. Agreement is generally very good, typically around the centimeter level or better for Low earth satellites. Much of the diffrence seems to be due to slight differences in the implementation details of GMST IAU1982. I think the differences are primarily numerical roundoff differences rather than algorithm differences. The polar motion has also been checked against skyfield and gives the same answer.
- Parameters:
- ecivnp.ndarray
A 1-D array of size[3,] or a 2-D array of size [3, N] where N is the number of vectors. The first column must be size 3 and stores the X,Y,Z components of each ECI vector.
- utcscalar or array
A scalar or array of coordinated universal times. The number of elements must match the number, N, of eci vectors The values must be convertible to astropy.time.Time values which includes, strings, datetime, numpy.datetime64, astropy.time.Time and floats. Float and strings should represent Coordinated Universal Time and floating point values must be expressed as a modified julian date.
- do_polar_motionbool
If true then apply a polar motion correction
- Returns:
- np.ndarray
Returns an array of geocentric vectors. The array is the same size, shape and units as parameter eciv
References
See Revisiting Spacetrack Report #3: Rev 1, D. Vallado, August 2006, AIAA 2006-6753 Appendix C. Also see function TEME_to_ITRF in package skyfield.sgp4lib.py and the matlab function teme2ecef on Matlab Central. The latter was written in 2007 by David Vallado, who seems to be the civilian expert on all matters relating to SGP4 and TEME ECI.
itrf_to_eciteme#
- itrf_to_eciteme(geov: ndarray, utc: datetime64 | datetime | float, do_polar_motion: bool = False) ndarray[source]#
Converts ITRF/GEO geocentric vectors to ECI TEME vectors. See the sister function
eciteme_to_itrffor specific details- Parameters:
- geovnp.ndarray
A 1-D array of size[3,] or a 2-D array of size [3, N] where N is the number of vectors. The first column must be size 3 and stores the X,Y,Z components of each ITRF vector.
- utcscalar or array
A scalar or array of coordinated universal times. The number of elements must match the number, N, of itrf vectors The values must be convertible to astropy.time.Time values which includes, strings, datetime, numpy.datetime64, astropy.time.Time and floats. Float and strings should represent Coordinated Universal Time and floating point values must be expressed as a modified julian date.
- do_polar_motionbool
If true then apply a polar motion correction
- Returns:
- np.ndarray
Returns an array of geocentric vectors. The array is the same size, shape and units as parameter geov.
gmst#
- gmst(utc: datetime64 | datetime | float | ndarray, **kwargs) float | ndarray[source]#
Calculates the Greenwich Mean Sidereal Time for the given Universal time(s). The Universal time is UTC. We use the astropy.time.Time.sidereal_time function to calculate sidereal time; this function uses IERS corrections for timescales and SOFA algorithms. It is is very accurate. The algorithm has given the same results as the novas.sidereal_time function to 10 or 11 decimal places.
- Parameters:
- utc: float, list or array of times.
A scalar,array or list of coordinated universal times. These Sideral time will be calculated for each of these elements. The UTC can be given in any of the formats supported by juliandate.ut_to_jd but are typically datetime.datetime or numpy.datetime.
- ** kwargsOptional keyword arguments
User can pass in optional keyword arguments for method astropy.time.Time.sideral_time. Typically used for applications that need to specify the IAU algorithm, eg IAU1982.
- Returns:
- float or numpy.ndarray
The greenwich mean sidereal time. This is expressed as a angle in radians between 0 of \(2\pi\).
polar_motion_rotation_matrix#
- polar_motion_rotation_matrix(utc, transpose=False) Tuple[float, float][source]#
Returns the polar motion rotation matrix for the given times. The rotation matrix is written as the matrix product of \(ROT1(y_p)*ROT2(x_p)\) where ROT1 is rotation around the global x axis and ROT2 is rotation around the global y axis.
\[\begin{split}\begin{equation} \mathbf{ROT1} = \left( \begin{array}{ccc} \cos x & 0 & -\sin x \\ 0 & 1 & 0 \\ \sin x & 0 & \cos x \end{array} \right) \end{equation}\end{split}\]\[\begin{split}\begin{equation} \mathbf{ROT2} = \left( \begin{array}{ccc} 1 & 0 & 0 \\ 0 & \cos y & \sin y \\ 0 & -\sin y & \cos y \end{array} \right) \end{equation}\end{split}\]- Parameters:
- utcscalar or array
A scalar or array [N] of coordinated universal times. The number of elements must match the number, N, of eci vectors The values must be convertible to astropy.time.Time values which includes, strings, datetime, numpy.datetime64, astropy.time.Time and floats. Float and strings should represent Coordinated Universal Time and floating point values must be expressed as a modified julian date.
- transposebool
If true then generate the transpose of the rotation matrix. Typically True is used when rotating from Psuedo Earth Fixed to ITRF. False will map ITRF to Psuedo Fixed Earth.
- Returns:
- np.ndarray
Returns an array of stacked rotation matrices. It will be of shape [3,3] if parameter utc is scalar. It will be of shape [N,3,3] if parameter utc is an array of N values.