adam_core.coordinates.transform module

adam_core.coordinates.transform.clear_translation_cache() None[source]

Clear the in-process translation cache used by cartesian_to_origin.

Primarily intended for testing and benchmarking.

adam_core.coordinates.transform.cartesian_to_geodetic(coords_cartesian: ndarray | Array, a: float, f: float, max_iter: int = 100, tol: float = 1e-15) Array[source]

Convert Cartesian coordinates to a geodetic coordinate.

Parameters:
  • coords_cartesian ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) – 3D Cartesian coordinate including time derivatives. x : x-position in units of distance. y : y-position in units of distance. z : z-position in units of distance. vx : x-velocity in the same units of x per arbitrary unit of time. vy : y-velocity in the same units of y per arbitrary unit of time. vz : z-velocity in the same units of z per arbitrary unit of time.

  • a (float (1)) – Semi-major axis of the Earth in units of distance.

  • f (float (1)) – Flattening of the Earth.

  • max_iter (int (1)) – Maximum number of iterations to perform.

  • tol (float (1)) – Tolerance for the iteration.

Returns:

coords_geodetic – 3D geodetic coordinate including time derivatives. alt : Altitude in units of distance. lon : Longitude in degrees. lat : Latitude in degrees. vup : Up velocity in the same units of x per arbitrary unit of time. veast : East velocity in degrees per arbitrary unit of time. vnorth : North velocity in degrees per arbitrary unit of time.

Return type:

~jax.numpy.ndarray (N, 6)

adam_core.coordinates.transform.cartesian_to_spherical(coords_cartesian: ndarray) ndarray[source]

Convert Cartesian coordinates to a spherical coordinates.

Parameters:

coords_cartesian ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) – 3D Cartesian coordinates including time derivatives. x : x-position in units of distance. y : y-position in units of distance. z : z-position in units of distance. vx : x-velocity in the same units of x per arbitrary unit of time. vy : y-velocity in the same units of y per arbitrary unit of time. vz : z-velocity in the same units of z per arbitrary unit of time.

Returns:

coords_spherical – 3D Spherical coordinates including time derivatives. rho : Radial distance in the same units of x, y, and z. lon : Longitude ranging from 0.0 to 360.0 degrees. lat : Latitude ranging from -90.0 to 90.0 degrees with 0 at the equator. vrho : Radial velocity in the same units as rho per arbitrary unit of time

(same unit of time as the x, y, and z velocities).

vlonLongitudinal velocity in degrees per arbitrary unit of time

(same unit of time as the x, y, and z velocities).

vlatLatitudinal velocity in degrees per arbitrary unit of time.

(same unit of time as the x, y, and z velocities).

Return type:

ndarray` (N, 6)

adam_core.coordinates.transform.spherical_to_cartesian(coords_spherical: ndarray | Array) Array[source]

Convert spherical coordinates to Cartesian coordinates.

Parameters:

coords_spherical ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) –

3D Spherical coordinates including time derivatives. rho : Radial distance in the same units of x, y, and z. lon : Longitude ranging from 0.0 to 360.0 degrees. lat : Latitude ranging from -90.0 to 90.0 degrees with 0 at the equator. vrho : Radial velocity in the same units as rho per arbitrary unit of time

(same unit of time as the x, y, and z velocities).

vlonLongitudinal velocity in degrees per arbitrary unit of time

(same unit of time as the x, y, and z velocities).

vlat :Latitudinal velocity in degrees per arbitrary unit of time.

(same unit of time as the x, y, and z velocities).

Returns:

coords_cartesian – 3D Cartesian coordinates including time derivatives. x : x-position in units of distance. y : y-position in units of distance. z : z-position in units of distance. vx : x-velocity in the same units of x per arbitrary unit of time. vy : y-velocity in the same units of y per arbitrary unit of time. vz : z-velocity in the same units of z per arbitrary unit of time.

Return type:

{~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)

adam_core.coordinates.transform.cartesian_to_keplerian(coords_cartesian: ndarray | Array, t0: ndarray | Array, mu: ndarray | Array) Array[source]

Convert Cartesian coordinates to Keplerian coordinates.

Parameters:
  • coords_cartesian ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) – 3D Cartesian coordinates including time derivatives. x : x-position in units of au. y : y-position in units of au. z : z-position in units of au. vx : x-velocity in units of au per day. vy : y-velocity in units of au per day. vz : z-velocity in units of au per day.

  • t0 ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Epoch at which cometary elements are defined in MJD TDB.

  • mu ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Gravitational parameter (GM) of the attracting body in units of au**3 / d**2.

Returns:

coords_keplerian – 13D Keplerian coordinates. a : semi-major axis in au. p : semi-latus rectum in au. q : periapsis distance in au. Q : apoapsis distance in au. e : eccentricity. i : inclination in degrees. raan : Right ascension (longitude) of the ascending node in degrees. ap : argument of periapsis in degrees. M : mean anomaly in degrees. nu : true anomaly in degrees. n : mean motion in degrees per day. P : period in days. tp : time of periapsis passage in days.

Return type:

~jax.numpy.ndarray (N, 13)

adam_core.coordinates.transform.keplerian_to_cartesian(coords_keplerian: ndarray | Array, mu: ndarray | Array, max_iter: int = 100, tol: float = 1e-15) Array[source]

Convert Keplerian coordinates to Cartesian coordinates.

Parabolic orbits (e = 1.0 +- 1e-15) with elements (a, e, i, raan, ap, M) cannot be converted to Cartesian orbits since their semi-major axes are by definition undefined. Please consider representing these orbits with Cometary elements and using those to convert to Cartesian. See ~adam_core.coordinates.cometary.cometary_to_cartesian.

Parameters:
  • coords_keplerian ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) – 6D Keplerian coordinate. a : semi-major axis in au. e : eccentricity. i : inclination in degrees. raan : Right ascension (longitude) of the ascending node in degrees. ap : argument of periapsis in degrees. M : mean anomaly in degrees.

  • mu ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Gravitational parameter (GM) of the attracting body in units of au**3 / d**2.

  • max_iter (int, optional) – Maximum number of iterations over which to converge. If number of iterations is exceeded, will use the value of the relevant anomaly at the last iteration.

  • tol (float, optional) – Numerical tolerance to which to compute anomalies using the Newtown-Raphson method.

Returns:

coords_cartesian – 3D Cartesian coordinates including time derivatives. x : x-position in units of au. y : y-position in units of au. z : z-position in units of au. vx : x-velocity in units of au per day. vy : y-velocity in units of au per day. vz : z-velocity in units of au per day.

Return type:

~jax.numpy.ndarray (N, 6)

Raises:

ValueError – When semi-major axis is less than 0 for elliptical orbits or when: semi-major axis is greater than 0 for hyperbolic orbits.

adam_core.coordinates.transform.cartesian_to_cometary(coords_cartesian: ndarray | Array, t0: ndarray | Array, mu: ndarray | Array) Array[source]

Convert Cartesian coordinates to Keplerian coordinates.

Parameters:
  • coords_cartesian ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) – 3D Cartesian coordinates including time derivatives. x : x-position in units of au. y : y-position in units of au. z : z-position in units of au. vx : x-velocity in units of au per day. vy : y-velocity in units of au per day. vz : z-velocity in units of au per day.

  • t0 ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Epoch at which cometary elements are defined in MJD TDB.

  • mu ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Gravitational parameter (GM) of the attracting body in units of au**3 / d**2.

Returns:

coords_cometary – 6D Cometary coordinates. q : periapsis distance in au. e : eccentricity. i : inclination in degrees. raan : Right ascension (longitude) of the ascending node in degrees. ap : argument of periapsis in degrees. tp : time of periapse passage in days.

Return type:

~jax.numpy.ndarray (N, 6)

adam_core.coordinates.transform.cometary_to_cartesian(coords_cometary: ndarray | Array, t0: ndarray | Array, mu: ndarray | Array, max_iter: int = 100, tol: float = 1e-15) Array[source]

Convert Cometary coordinates to Cartesian coordinates.

Parameters:
  • coords_cometary ({~numpy.ndarray, ~jax.numpy.ndarray} (N, 6)) – 6D Cometary coordinate. q : periapsis distance in au. e : eccentricity. i : inclination in degrees. raan : Right ascension (longitude) of the ascending node in degrees. ap : argument of periapsis in degrees. tp : time of periapse passage in days.

  • t0 ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Epoch at which cometary elements are defined in MJD TDB.

  • mu ({~numpy.ndarray, ~jax.numpy.ndarray} (N)) – Gravitational parameter (GM) of the attracting body in units of au**3 / d**2.

  • max_iter (int, optional) – Maximum number of iterations over which to converge. If number of iterations is exceeded, will use the value of the relevant anomaly at the last iteration.

  • tol (float, optional) – Numerical tolerance to which to compute anomalies using the Newtown-Raphson method.

Returns:

coords_cartesian – 3D Cartesian coordinates including time derivatives. x : x-position in units of au. y : y-position in units of au. z : z-position in units of au. vx : x-velocity in units of au per day. vy : y-velocity in units of au per day. vz : z-velocity in units of au per day.

Return type:

~jax.numpy.ndarray (N, 6)

adam_core.coordinates.transform.cartesian_to_origin(coords: CartesianCoordinates, origin: OriginCodes) CartesianCoordinates[source]

Translate coordinates to a different origin.

Parameters:
  • coords (~adam_core.coordinates.cartesian.CartesianCoordinates) – Cartesian coordinates and optionally their covariances.

  • origin (~adam_core.coordinates.origin.OriginCodes) – Desired origin. Input origins may be either ~adam_core.coordinates.origin.OriginCodes or a str of an observatory code, but the output origin (this kwarg) should always be an ~adam_core.coordinates.origin.OriginCodes. If you are looking to generate ephemerides for an observatory, please use a ~adam_core.propagator.propagator.Propagator instead.

Returns:

CartesianCoordinates – Translated Cartesian coordinates and their covariances.

Return type:

~adam_core.coordinates.cartesian.CartesianCoordinates

Raises:

ValueError – If origin is not a ~adam_core.coordinates.origin.OriginCodes object or a str of an observatory code.

adam_core.coordinates.transform.apply_time_varying_rotation(coords: CartesianCoordinates, frame_out: Literal['ecliptic', 'equatorial', 'itrf93']) CartesianCoordinates[source]

Apply a time-varying rotation to the Cartesian coordinates. At the moment only rotations between itrf93 and ecliptic/equatorial are supported.

Parameters:
  • coords (CartesianCoordinates) – The coordinates to rotate and translate to the desired frame.

  • frame_out (Literal["ecliptic", "equatorial", "itrf93"]) – The desired output frame.

Returns:

The rotated coordinates.

Return type:

CartesianCoordinates

adam_core.coordinates.transform.cartesian_to_frame(coords: CartesianCoordinates, frame: Literal['ecliptic', 'equatorial', 'itrf93']) CartesianCoordinates[source]

Rotate Cartesian coordinates and their covariances to the given frame.

Parameters:
  • coords (~adam_core.coordinates.cartesian.CartesianCoordinates) – Cartesian coordinates and optionally their covariances.

  • frame ({'ecliptic', 'equatorial', 'itrf93'}) – Desired reference frame of the output coordinates.

Returns:

CartesianCoordinates – Rotated Cartesian coordinates and their covariances.

Return type:

~adam_core.coordinates.cartesian.CartesianCoordinates

adam_core.coordinates.transform.transform_coordinates(coords: CoordinateType, representation_out: type[CoordinateType] | None = None, frame_out: Literal['ecliptic', 'equatorial', 'itrf93'] | None = None, origin_out: OriginCodes | None = None) CoordinateType[source]

Transform coordinates between frames (‘ecliptic’, ‘equatorial’, ‘itrf93’), origins, and/or representations (‘cartesian’, ‘spherical’, ‘keplerian’, ‘cometary’).

Input coordinates may be defined from multiple origins but if origin_out is specified, all coordinates will be transformed to that origin.

Parameters:
  • coords – Coordinates to transform between representations and frames.

  • representation_out – Desired coordinate type or representation of the output coordinates. If None, the output coordinates will be the same type as the input coordinates.

  • frame_out ({'ecliptic', 'equatorial'}) – Desired reference frame of the output coordinates.

  • origin_out (~adam_core.coordinates.origin.OriginCodes) – Desired origin. Input origins may be either ~adam_core.coordinates.origin.OriginCodes or a str of an observatory code, but the output origin (this kwarg) should always be an ~adam_core.coordinates.origin.OriginCodes. If you are looking to generate ephemerides for an observatory, please use a ~adam_core.propagator.propagator.Propagator instead.

Returns:

coords_out – Coordinates in desired output representation and frame.

Return type:

~adam_core.coordinates.Coordinates

Raises:
  • ValueError – If frame_in, frame_out are not one of ‘equatorial’, ‘ecliptic’. If representation_in, representation_out are not one of ‘cartesian’, ‘spherical’, ‘keplerian’, ‘cometary’, ‘geodetic’.

  • TypeError – If coords is not a ~adam_core.coordinates.Coordinates object. If origin_out is not a ~adam_core.coordinates.OriginCodes object.