Orbital States
Representations of orbital states in Cartesian and Keplerian forms.
Quick Example
import lox_space as lox
t = lox.Time("TAI", 2024, 1, 1)
# Array pattern — position in meters, velocity in m/s
state = lox.Cartesian(
time=t,
position=[6678e3, 0.0, 0.0],
velocity=[0.0, 7730.0, 0.0],
)
# Elementwise pattern — unitful keyword arguments
state = lox.Cartesian(
time=t,
x=6678.0 * lox.km, y=0.0 * lox.km, z=0.0 * lox.km,
vx=0.0 * lox.km_per_s, vy=7.73 * lox.km_per_s, vz=0.0 * lox.km_per_s,
)
# Access as numpy arrays (meters, m/s)
pos = state.position() # shape (3,)
vel = state.velocity() # shape (3,)
# Access as unitful scalars
print(f"x = {state.x.to_kilometers():.1f} km")
# Convert to Keplerian elements
kep = state.to_keplerian()
print(f"a = {kep.semi_major_axis().to_kilometers():.1f} km")
print(f"e = {kep.eccentricity():.6f}")
print(f"i = {kep.inclination().to_degrees():.2f} deg")
# Create from Keplerian elements (semi-major axis + eccentricity)
kep = lox.Keplerian(
time=t,
semi_major_axis=6678.0 * lox.km,
eccentricity=0.001,
inclination=51.6 * lox.deg,
)
state = kep.to_cartesian()
# Create from periapsis/apoapsis radii
kep = lox.Keplerian(
time=t,
periapsis_radius=7000.0 * lox.km,
apoapsis_radius=7400.0 * lox.km,
)
# Create from periapsis/apoapsis altitudes
kep = lox.Keplerian(
time=t,
periapsis_altitude=600.0 * lox.km,
apoapsis_altitude=1000.0 * lox.km,
)
# Create a circular orbit from altitude
circ = lox.Keplerian.circular(
t,
altitude=800 * lox.km,
inclination=51.6 * lox.deg,
)
# Create a Sun-synchronous orbit from altitude
eop = lox.EOPProvider("finals.csv")
sso = lox.Keplerian.sso(
t,
altitude=800 * lox.km,
ltan=(13, 30),
provider=eop,
)
# Work with trajectories
trajectory = lox.Trajectory([state1, state2, state3])
interpolated = trajectory.interpolate(t + lox.TimeDelta(100))
Cartesian
Represents an orbital state (position and velocity) at a specific time.
Parameters:
-
–timeThe epoch of this state.
-
–positionPosition vector as array-like [x, y, z] in meters.
-
–velocityVelocity vector as array-like [vx, vy, vz] in m/s.
-
–x, y, zIndividual position components as Distance (keyword-only, alternative to position).
-
–vx, vy, vzIndividual velocity components as Velocity (keyword-only, alternative to velocity).
-
–originCentral body (default: Earth).
-
–frameReference frame (default: ICRF).
Examples:
>>> t = lox.Time("TAI", 2024, 1, 1)
>>> # Array pattern (meters, m/s)
>>> state = lox.Cartesian(
... t,
... position=[6678e3, 0.0, 0.0],
... velocity=[0.0, 7730.0, 0.0],
... )
>>> # Elementwise pattern (unitful)
>>> state = lox.Cartesian(
... t,
... x=6678.0 * lox.km, y=0.0 * lox.km, z=0.0 * lox.km,
... vx=0.0 * lox.km_per_s, vy=7.73 * lox.km_per_s, vz=0.0 * lox.km_per_s,
... )
Methods:
-
origin–Return the central body (origin) of this state.
-
position–Return the position vector as a numpy array in meters, shape (3,).
-
reference_frame–Return the reference frame of this state.
-
rotation_lvlh–Compute the rotation matrix from inertial to LVLH frame.
-
time–Return the epoch of this state.
-
to_frame–Transform this state to a different reference frame.
-
to_ground_location–Convert this state to a ground location.
-
to_keplerian–Convert this Cartesian state to Keplerian orbital elements.
-
to_origin–Transform this state to a different central body.
-
velocity–Return the velocity vector as a numpy array in m/s, shape (3,).
Attributes:
-
vx(Velocity) –Return the x component of the velocity.
-
vy(Velocity) –Return the y component of the velocity.
-
vz(Velocity) –Return the z component of the velocity.
-
x(Distance) –Return the x component of the position.
-
y(Distance) –Return the y component of the position.
-
z(Distance) –Return the z component of the position.
position
position() -> ndarray
Return the position vector as a numpy array in meters, shape (3,).
rotation_lvlh
rotation_lvlh() -> ndarray
Compute the rotation matrix from inertial to LVLH frame.
to_frame
to_frame(frame: str | Frame, provider: EOPProvider | None = None) -> Self
Transform this state to a different reference frame.
to_keplerian
to_keplerian() -> Keplerian
Convert this Cartesian state to Keplerian orbital elements.
to_origin
Transform this state to a different central body.
velocity
velocity() -> ndarray
Return the velocity vector as a numpy array in m/s, shape (3,).
Keplerian
Represents an orbit using Keplerian (classical) orbital elements.
The orbital shape can be specified in three ways:
semi_major_axis+eccentricityperiapsis_radius+apoapsis_radius(keyword-only)periapsis_altitude+apoapsis_altitude(keyword-only)
Parameters:
-
–timeEpoch of the elements.
-
–semi_major_axisSemi-major axis as Distance.
-
–eccentricityOrbital eccentricity (0 = circular, <1 = elliptical).
-
–inclinationInclination as Angle (default 0).
-
–longitude_of_ascending_nodeRAAN as Angle (default 0).
-
–argument_of_periapsisArgument of periapsis as Angle (default 0).
-
–true_anomalyTrue anomaly as Angle (default 0).
-
–originCentral body (default: Earth).
-
–periapsis_radiusPeriapsis radius as Distance (keyword-only).
-
–apoapsis_radiusApoapsis radius as Distance (keyword-only).
-
–periapsis_altitudePeriapsis altitude as Distance (keyword-only).
-
–apoapsis_altitudeApoapsis altitude as Distance (keyword-only).
-
–mean_anomalyMean anomaly as Angle (keyword-only, mutually exclusive with true_anomaly).
Examples:
>>> t = lox.Time("TAI", 2024, 1, 1)
>>> orbit = lox.Keplerian(
... t,
... semi_major_axis=6678.0 * lox.km,
... eccentricity=0.001,
... inclination=51.6 * lox.deg,
... )
From radii:
>>> orbit = lox.Keplerian(
... t,
... periapsis_radius=7000.0 * lox.km,
... apoapsis_radius=7400.0 * lox.km,
... )
From altitudes:
>>> orbit = lox.Keplerian(
... t,
... periapsis_altitude=600.0 * lox.km,
... apoapsis_altitude=1000.0 * lox.km,
... )
Methods:
-
argument_of_periapsis–Return the argument of periapsis.
-
circular–Construct a circular orbit.
-
eccentricity–Return the orbital eccentricity.
-
inclination–Return the inclination.
-
longitude_of_ascending_node–Return the longitude of the ascending node (RAAN).
-
orbital_period–Return the orbital period.
-
origin–Return the central body (origin) of this orbit.
-
semi_major_axis–Return the semi-major axis.
-
sso–Construct a Sun-synchronous orbit.
-
time–Return the epoch of these elements.
-
to_cartesian–Convert these Keplerian elements to a Cartesian state.
-
true_anomaly–Return the true anomaly.
circular
classmethod
circular(
time: Time,
*,
semi_major_axis: Distance | None = None,
altitude: Distance | None = None,
inclination: Angle | None = None,
longitude_of_ascending_node: Angle | None = None,
true_anomaly: Angle | None = None,
origin: str | int | Origin | None = None
) -> Self
Construct a circular orbit.
Exactly one of semi_major_axis or altitude must be provided.
Eccentricity is always 0 and argument of periapsis is always 0.
Parameters:
-
(timeTime) –Epoch of the orbit.
-
(semi_major_axisDistance | None, default:None) –Semi-major axis (mutually exclusive with altitude).
-
(altitudeDistance | None, default:None) –Orbital altitude (mutually exclusive with semi_major_axis).
-
(inclinationAngle | None, default:None) –Inclination (default 0).
-
(longitude_of_ascending_nodeAngle | None, default:None) –RAAN (default 0).
-
(true_anomalyAngle | None, default:None) –True anomaly (default 0).
-
(originstr | int | Origin | None, default:None) –Central body (default: Earth).
Examples:
>>> t = lox.Time("TAI", 2024, 1, 1)
>>> orbit = lox.Keplerian.circular(
... t,
... altitude=800 * lox.km,
... inclination=51.6 * lox.deg,
... )
longitude_of_ascending_node
longitude_of_ascending_node() -> Angle
Return the longitude of the ascending node (RAAN).
sso
classmethod
sso(
time: Time,
*,
altitude: Distance | None = None,
semi_major_axis: Distance | None = None,
inclination: Angle | None = None,
eccentricity: float = 0.0,
ltan: tuple[int, int] | None = None,
ltdn: tuple[int, int] | None = None,
argument_of_periapsis: Angle | None = None,
true_anomaly: Angle | None = None,
provider: EOPProvider | None = None
) -> Self
Construct a Sun-synchronous orbit.
Exactly one of altitude, semi_major_axis, or inclination
must be provided. The remaining orbital elements are derived from the
SSO constraint.
Parameters:
-
(timeTime) –Epoch of the orbit.
-
(altitudeDistance | None, default:None) –Orbital altitude (mutually exclusive with semi_major_axis/inclination).
-
(semi_major_axisDistance | None, default:None) –Semi-major axis (mutually exclusive with altitude/inclination).
-
(inclinationAngle | None, default:None) –Inclination (mutually exclusive with altitude/semi_major_axis).
-
(eccentricityfloat, default:0.0) –Eccentricity (default 0.0).
-
(ltantuple[int, int] | None, default:None) –Local time of ascending node as
(hours, minutes)tuple. -
(ltdntuple[int, int] | None, default:None) –Local time of descending node as
(hours, minutes)tuple. -
(argument_of_periapsisAngle | None, default:None) –Argument of periapsis (default 0.0).
-
(true_anomalyAngle | None, default:None) –True anomaly (default 0.0).
-
(providerEOPProvider | None, default:None) –EOP provider for time scale conversions.
Examples:
>>> eop = lox.EOPProvider("finals.csv")
>>> t = lox.UTC(2020, 2, 18).to_scale("TDB")
>>> orbit = lox.Keplerian.sso(
... t,
... altitude=800 * lox.km,
... ltan=(13, 30),
... provider=eop,
... )
Trajectory
A time-series of orbital states with interpolation support.
Parameters:
-
–statesList of Cartesian objects in chronological order.
Examples:
>>> trajectory = propagator.propagate(times)
>>> state = trajectory.interpolate(t)
>>> arr = trajectory.to_numpy()
Methods:
-
find_events–Find events where a function crosses zero.
-
find_windows–Find time windows where a function is positive.
-
from_numpy–Create a Trajectory from a numpy array.
-
interpolate–Interpolate the trajectory at a specific time.
-
origin–Return the central body (origin) of this trajectory.
-
reference_frame–Return the reference frame of this trajectory.
-
states–Return the list of states in this trajectory.
-
to_frame–Transform all states to a different reference frame.
-
to_numpy–Export trajectory to a numpy array with columns [t(s), x(m), y(m), z(m), vx(m/s), vy(m/s), vz(m/s)].
-
to_origin–Transform all states to a different central body.
find_events
Find events where a function crosses zero.
find_windows
Find time windows where a function is positive.
from_numpy
classmethod
from_numpy(
start_time: Time,
states: ndarray,
origin: str | int | Origin | None = None,
frame: str | Frame | None = None,
) -> Self
Create a Trajectory from a numpy array.
Parameters:
-
(start_timeTime) –Reference epoch for the trajectory.
-
(statesndarray) –2D array with columns [t(s), x(m), y(m), z(m), vx(m/s), vy(m/s), vz(m/s)].
-
(originstr | int | Origin | None, default:None) –Central body (default: Earth).
-
(framestr | Frame | None, default:None) –Reference frame (default: ICRF).
interpolate
Interpolate the trajectory at a specific time.
to_frame
to_frame(frame: str | Frame, provider: EOPProvider | None = None) -> Self
Transform all states to a different reference frame.
to_numpy
to_numpy() -> ndarray
Export trajectory to a numpy array with columns [t(s), x(m), y(m), z(m), vx(m/s), vy(m/s), vz(m/s)].