Coordinate Systems and Frames

In SatMAD, frames and coordinate systems as well as conversions between them are handled through Astropy. This section introduces the additional frames defined by SatMAD.

Local Frames of Celestial Bodies (CelestialBody, CelestialBodyJ2000Equatorial and CelestialBodyTODEquatorial)

In addition to Earth, it is possible to define a Celestial Body in space through the CelestialBody class. For each such Celestial Body, it is then possible to realise local reference frames.

The first of them is a local inertial or “Celestial Reference System” (CRS) (CelestialBodyCRS). This class simply translates the ICRS coordinate system to the centre of the celestial body. As such, the XY plane does not correspond to the equator plane. This enables the user to define coordinates in this local inertial coordinate system and then (as an example) run an orbit propagation around it or run an analysis. For example, the following would define a “Sun CRS” (equivalent to Heliocentric Celestial Reference System), and a “Moon CRS” by simply subclassing CelestialBodyCRS.

from satmad.coordinates.frames import CelestialBodyCRS

class SunCRS(CelestialBodyCRS):
    body = "Sun"


class MoonCRS(CelestialBodyCRS):
    body = "Moon"
    ephemeris_type = "jpl"

In the latter, the optional ephemeris_type parameter determines which ephemeris to use when computing the location of the Celestial Bodies. Note that, while some frames like MoonCRS and MarsCRS are predefined, and other Planets can be easily created as shown.

Similarly, one other inertial and one nutating-precessing frame can be defined for each Celestial Body: CelestialBodyJ2000Equatorial and CelestialBodyTODEquatorial. The former is an inertial equatorial frame with the alignment fixed at J2000 Epoch. The latter is an equatorial frame but its orientation is computed with the instantaneous orientation due to planetary nutation and precession:

\[ \vec{r}_{CRS} = R_x(90+ \alpha) R_z(90- \delta)\times \vec{r}_{TOD} \]

New frames can be defined, inspecting the preset Mars frames:

from satmad.coordinates.frames import CelestialBodyJ2000Equatorial, CelestialBodyTODEquatorial, MarsCRS

class MarsTODEquatorial(CelestialBodyTODEquatorial):
    body_name = "Mars"
    cb_crs = MarsCRS


class MarsJ2000Equatorial(CelestialBodyJ2000Equatorial):
    body_name = "Mars"
    cb_crs = MarsCRS

Similar to the CelestialBody class, body_name defines the name of the body, such that its coordinates can be computed within the ICRS using builtin or JPL DE4XX ephemerides. The cb_crs parameter links this frame to CRS inertial frame of this body, as SatMAD (or Astropy) would not know how to chain the coordinate conversions from this frame to (for example) HCRS.

The definition of these coordinate systems require the modelling of how the planet is oriented with respect to ICRS. This is done via the Celestial Body Rotation parameters: right ascension and declination of the celestial body North Pole, and prime meridian angle. As can be imagined, these are different for each planet - the values used in this model are taken from Report on IAU Working Group on Cartographic Coordinates and Rotational Elements: 2015 [TCF2], except for the Moon, which is taken from the 2009 version of the same report. As such, these values are incompatible with, yet much more up-to-date than GMAT and NASA HORIZONS Web Interface, which use the much simpler 2000 version of these models. The TOD and J2000 Equatorial frame definitions are given in [OM1], [OM3] and [OM4].

The CelestialBodyCRS class corresponds to setting the frame as “ICRF” around a central body in GMAT. The CelestialBodyJ2000Equatorial is equivalent to “Body Inertial” and CelestialBodyTODEquatorial to “Body Equatorial” in GMAT.

Earth Based Additional Frames (J2000)

The built-in frames offered by Astropy do not include some frames that are used in satellite applications. To bridge this gap, this package offers Mean Pole and Equinox at J2000.0 Reference System (J2000).

J2000 coordinate frame is similar to GCRS but rotated by a constant frame bias [TCF1]:

\[ \vec{r}_{J2000} = B \times \vec{r}_{GCRS} \]

This rotation is applicable only to the equinox based approach, and is only an approximation. The difference between GCRS and J2000 is less than 1 metre for the Low Earth Orbit, therefore these two can be used interchangeably with a small error.

The J2000 class is similar to (and compatible with) the Astropy Built-in Frames.

from astropy import units as u
from astropy.coordinates import CartesianRepresentation, CartesianDifferential
from satmad.coordinates.frames import J2000
from astropy.time import Time

time = Time("2020-01-11T11:00:00.000", scale="utc")
v_j2000 = CartesianDifferential([-4.7432196000, 0.7905366000, 5.5337561900], unit=u.km / u.s)
r_j2000 = CartesianRepresentation([5102.50960000, 6123.01152000, 6378.13630000], unit=u.km)
rv_j2000 = J2000(r_j2000.with_differentials(v_j2000), obstime=time, representation_type="cartesian", differential_type="cartesian")

Reference / API

Coordinate systems and frames defined by satmad.

class satmad.coordinates.frames.CelestialBodyCRS(*args, **kwargs)

A coordinate frame in the generic Celestial Reference System (CRS). This CRS is derived from ICRS by simply carrying the origin to the origin of the celestial body.

Uses the builtin ephemeris to compute planetary positions in ICRS, unless the coordinate definition explicitly specifies an ephemeris type.

property default_differential
property default_representation
ephemeris_type = 'builtin'
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'celestialbodycrs'
obstime = <Time object: scale='tt' format='jyear_str' value=J2000.000>
class satmad.coordinates.frames.CelestialBodyJ2000Equatorial(*args, **kwargs)

A coordinate frame representing the Equatorial System of the Celestial Body at J2000 Epoch. This is not defined for the Earth as it has its own J2000 Equatorial frame.

The axis orientations are:

  • x-axis: Along the line formed by the intersection of the body equator and the x-y plane of the FK5 system, at the J2000 epoch

  • y-axis: Completes the right-handed set.

  • z-axis: Along the instantaneous body spin axis direction at the J2000 epoch

This corresponds to “Body Inertial” in GMAT.

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'celestialbodyj2000equatorial'
obstime = <Time object: scale='tt' format='jyear_str' value=J2000.000>
class satmad.coordinates.frames.CelestialBodyTODEquatorial(*args, **kwargs)

A coordinate frame representing the True-of-Date Equatorial System of the Celestial Body. This is not defined for the Earth as it has its own TOD

The axis orientations are:

  • x-axis: Along the line formed by the intersection of the body equator and the ecliptic planes.

  • y-axis: Completes the right-handed set.

  • z-axis: Normal to the equatorial plane.

This corresponds to “Equator System” in GMAT.

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'celestialbodytodequatorial'
obstime = <Time object: scale='tt' format='jyear_str' value=J2000.000>
class satmad.coordinates.frames.J2000(*args, copy=True, representation_type=None, differential_type=None, **kwargs)

A coordinate or frame in the Mean Pole and Equinox at J2000.0 Reference System (J2000).

This coordinate frame is similar to GCRS but rotated by the frame bias. This rotation is applicable only to the equinox based approach, and is only an approximation. The difference betweenGCRS and J2000 is less than 1m for the Low Earth Orbit, therefore these two can be used interchangeably with a small error.

This frame is for legacy applications that store data in J2000. GCRS should be used wherever possible.

References

The definitions and conversions are from IERS Conventions 2010, Chapter 5 [TCF1] in References.

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'j2000'
obstime = <Time object: scale='tt' format='jyear_str' value=J2000.000>
class satmad.coordinates.frames.MarsCRS(*args, **kwargs)

Mars Celestial Reference System. This is simply the ICRS shifted to the centre of the Mars with the velocity adjusted with respect to the Mars.

body_name = 'Mars'
property default_differential
property default_representation
ephemeris_type = 'jpl'
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'marscrs'
class satmad.coordinates.frames.MarsJ2000Equatorial(*args, **kwargs)

A coordinate frame representing the Equatorial System of Mars at J2000 Epoch.

body_name = 'Mars'
cb_crs

alias of satmad.coordinates.frames.MarsCRS

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'marsj2000equatorial'
class satmad.coordinates.frames.MarsTODEquatorial(*args, **kwargs)

A coordinate frame representing the True-of-Date Equatorial System of Mars.

body_name = 'Mars'
cb_crs

alias of satmad.coordinates.frames.MarsCRS

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'marstodequatorial'
class satmad.coordinates.frames.MoonCRS(*args, **kwargs)

Moon Celestial Reference System. This is simply the ICRS shifted to the centre of the Moon with the velocity adjusted with respect to the Moon.

This uses the ephemeris type jpl to ensure that Moon velocity is computed (builtin cannot compute velocity). This is critical for coordinate transformations that involve velocity.

body_name = 'Moon'
property default_differential
property default_representation
ephemeris_type = 'jpl'
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'mooncrs'
class satmad.coordinates.frames.MoonJ2000Equatorial(*args, **kwargs)

A coordinate frame representing the Equatorial System of Moon at J2000 Epoch.

This corresponds to the Moon Principal Axis (PA) system. This has a fixed and small rotational offset with respect to the Mean Earth/rotation (ME) system.

body_name = 'Moon'
cb_crs

alias of satmad.coordinates.frames.MoonCRS

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'moonj2000equatorial'
class satmad.coordinates.frames.TIRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)

A coordinate or frame in the Terrestrial Intermediate Reference System (TIRS).

References

The definitions and conversions are from IERS Conventions 2010, Chapter 5 [TCF1] in References.

property default_differential
property default_representation
frame_attributes = {'obstime': <astropy.coordinates.attributes.TimeAttribute object>}
property frame_specific_representation_info
name = 'tirs'
obstime = <Time object: scale='tt' format='jyear_str' value=J2000.000>
satmad.coordinates.frames.cb_crs_to_cb_j2000_eq(cb_crs_coord, _)

Conversion from Celestial Reference System of a Central Body to its Local Equatorial at J2000 Epoch Reference System.

satmad.coordinates.frames.cb_crs_to_cb_tod_eq(cb_crs_coord, _)

Conversion from Celestial Reference System of a Central Body to its Local Equatorial True-of-Date Reference System.

satmad.coordinates.frames.cb_crs_to_icrs(cb_crs_coord, _)

Conversion from Celestial Reference System of a Central Body to ICRS.

satmad.coordinates.frames.cb_j2000_eq_to_cb_crs(cb_eq_coord, _)

Conversion from Local Equatorial at J2000 Epoch Reference System of a Central Body to its Celestial Reference System.

satmad.coordinates.frames.cb_tod_eq_to_cb_crs(cb_eq_coord, _)

Conversion from Local True-of-Date Equatorial Reference System of a Central Body to its Celestial Reference System.

satmad.coordinates.frames.gcrs_to_j2000()

Constant conversion matrix (Frame Bias) from GCRS to J2000.

satmad.coordinates.frames.icrs_to_cb_crs(icrs_coord, cb_crs_frame)

Conversion from ICRS to Celestial Reference System of a Central Body.

satmad.coordinates.frames.init_pvt(frame, time, pos, vel=None, copy=False)

Convenience method to initialise a SkyCoord object using the inputs.

If velocity input is not given (or is None), then velocity input is not set. If a velocity value of (0,0,0) is desired, these values should be explicitly assigned and the corresponding velocity vector should be given as input.

Note that, the position vector input may be including the velocity inputs. If velocity input is not given, then these velocity values are not overwritten.

Also note that, Time, CartesianRepresentation and CartesianDifferential objects as well as the resulting SkyCoord object support multiple value inputs, as long as the input sizes match. As an example, it is possible to input 30 time, position and velocity values, representing discrete points on a trajectory.

Parameters
  • frame (str or Type[BaseCoordinateFrame]) – Frame name (string) or frame object itself (e.g.GCRS)

  • time (Time) – time(s) associated with the coordinates

  • pos (CartesianRepresentation) – Position vector(s)

  • vel (CartesianDifferential or None) – Velocity vector(s)

  • copy (bool, optional) – If True (default), a copy of any coordinate data is made.

Returns

SkyCoord object representing the cartesian input

Return type

SkyCoord

satmad.coordinates.frames.itrs_to_tirs(itrs_coord, _)

Dynamic conversion matrix (Polar Motion) from ITRS to TIRS.

satmad.coordinates.frames.j2000_to_gcrs()

Constant conversion matrix (Frame Bias) from J2000 to GCRS.

satmad.coordinates.frames.tirs_to_itrs(tirs_coord, _)

Dynamic conversion matrix (Polar Motion) from TIRS to ITRS.

Planetary rotational parameters, defined relative to their mean axis of rotation and body-dependent definitions of longitude.

All formulas (with the exception of the Moon) are from: Archinal, B.A., Acton, C.H., A’Hearn, M.F. et al. Report of the IAU Working Group on Cartographic Coordinates and Rotational Elements: 2015. Celest Mech Dyn Astr 130, 22 (2018). https://doi.org/10.1007/s10569-017-9805-5

Moon formulas are from the 2009 edition as the 2015 version no longer lists the explicit lunar formulae - they refer to DE430 or the older DE421 values. Archinal, B.A., A’Hearn, M.F., Bowell, E. et al. Report of the IAU Working Group on Cartographic Coordinates and Rotational Elements: 2009. Celest Mech Dyn Astr 109, 101–135 (2011). https://doi.org/10.1007/s10569-010-9320-4

satmad.coordinates.planetary_rot_elems.celestial_body_rot_params(epoch, body_name)

Retrieves the rotational parameters for planets and satellites for the requested epoch.

These parameters link the Celestial Body ICRS to the equatorial inertial or body-fixed coordinates.

Parameters
  • epoch (time) – Epoch at which rotational parameters are evaluated, defaults to J2000 epoch

  • body_name (str) – Name of the celestial body (planet or satellite)

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[Quantity]

Raises

NotImplementedError – Rotational parameters not available for the requested body name

satmad.coordinates.planetary_rot_elems.jupiter_rot_params(t, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Jupiter.

Parameters
  • t (float) – Julian Centuries since J2000 epoch

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.mars_rot_params(t, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Mars.

Parameters
  • t (float) – Julian Centuries since J2000 epoch

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.mercury_rot_params(t, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Mercury.

Parameters
  • t (float) – Julian Centuries since J2000 epoch

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.moon_rot_params(t, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for the Moon.

These parameters are useful to convert local ICRS frame to the Principal Axis (PA) frame of the Moon. The Mean Earth/rotation axis (ME) has an additional fixed offset with respect to PA.

Parameters
  • t (float) – Julian Centuries since J2000 epoch

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.neptune_rot_params(t, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Neptune.

Parameters
  • t (float) – Julian Centuries since J2000 epoch

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.saturn_rot_params(t, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Saturn.

Parameters
  • t (float) – Julian Centuries since J2000 epoch

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.sun_rot_params(_, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for the Sun.

Parameters
  • _ – Empty parameter

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.uranus_rot_params(_, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Uranus.

Parameters
  • _ – Empty parameter

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]

satmad.coordinates.planetary_rot_elems.venus_rot_params(_, d)

Computes the North Pole rotational parameters (ra, dec, prime meridian angle) for Venus.

Parameters
  • _ – Empty parameter

  • d (float) – Days since J2000 epoch

Returns

ra, dec, w – Right ascension and declination of celestial body north pole, and prime meridian angle

Return type

Tuple[float]