API Reference

AstroStates.AlternateEquinoctialStateType
AlternateEquinoctialState(a, h, k, altp, altq, λ)

Alternate equinoctial elements representation.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • a: Semi-major axis. Defines orbit size and energy.
    • Range: a ≠ 0. If a > 0: elliptic orbit. If a < 0: hyperbolic orbit.
  • h: Eccentricity vector h-component. h = e⋅sin(ω + Ω).
    • Range: Any real value. Related to eccentricity and orientation.
  • k: Eccentricity vector k-component. k = e⋅cos(ω + Ω).
    • Range: Any real value. Related to eccentricity and orientation.
  • altp: Alternate inclination vector p-component. altp = sin(i/2)⋅cos(Ω).
    • Range: Any real value. Alternative inclination parameterization.
  • altq: Alternate inclination vector q-component. altq = sin(i/2)⋅sin(Ω).
    • Range: Any real value. Alternative inclination parameterization.
  • mlong: Mean longitude (rad). Combined angle measure.
    • Range: [0, 2π). Normalized mean longitude.

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Alternative to standard equinoctial elements using sin(i/2) instead of tan(i/2).
  • Singularities reduced compared to classical elements.

Examples

alteq = AlternateEquinoctialState(7000.0, 0.01, 0.0, 0.05, 0.0, π/4)
source
AstroStates.CartesianStateType
CartesianState(posvel)
CartesianState(pos, vel)

Cartesian position and velocity state representation.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.

Fields

  • position::SVector{3,T}: Position vector [x, y, z]
    • Range: Any real values.
  • velocity::SVector{3,T}: Velocity vector [vx, vy, vz]
    • Range: Any real values.

Properties

  • posvel: Backward-compatible property returning combined state as Vector{Float64}

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Internal storage uses StaticArrays for efficient coordinate transformations.
  • The posvel property provides backward compatibility with code expecting a 6-element vector.
  • See docstrings for other constructor options for flexibility.

Examples

# From 6-element vector
cart = CartesianState([6778.0, 0.0, 0.0, 0.0, 7.66, 0.0])

# From separate position and velocity
using StaticArrays
pos = SVector(6778.0, 0.0, 0.0)
vel = SVector(0.0, 7.66, 0.0)
cart = CartesianState(pos, vel)

# Access fields
x, y, z = cart.position
vx, vy, vz = cart.velocity
posvel_vec = cart.posvel  # Returns Vector{Float64}
source
AstroStates.CartesianStateMethod
CartesianState(pos::AbstractVector, vel::AbstractVector)

Construct CartesianState from separate 3-element position and velocity vectors.

source
AstroStates.CartesianStateMethod
CartesianState(pos::SVector{3}, vel::SVector{3})

Construct CartesianState from separate position and velocity static vectors.

source
AstroStates.EquinoctialStateType
EquinoctialState(a, h, k, p, q, λ)

Equinoctial orbital elements representation.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • a: Semi-major axis. Defines orbit size and energy.
    • Range: a ≠ 0. If a > 0: elliptic orbit. If a < 0: hyperbolic orbit.
  • h: Eccentricity vector h-component. h = e⋅sin(ω + Ω).
    • Range: Any real value. Related to eccentricity and orientation.
  • k: Eccentricity vector k-component. k = e⋅cos(ω + Ω).
    • Range: Any real value. Related to eccentricity and orientation.
  • p: Inclination vector p-component. p = tan(i/2)⋅cos(Ω).
    • Range: Any real value. Related to inclination and node orientation.
  • q: Inclination vector q-component. q = tan(i/2)⋅sin(Ω).
    • Range: Any real value. Related to inclination and node orientation.
  • mlong: Mean longitude (rad). Ω + ω + ν combined angle measure.
    • Range: [0, 2π). Normalized mean longitude.

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Singularities are reduced compared to classical elements but still exist for retrograde equatorial orbits.

Examples

eq = EquinoctialState(7000.0, 0.01, 0.0, 0.1, 0.0, π/4)
source
AstroStates.IncomingAsymptoteStateType
IncomingAsymptoteState(rp, c3, rla, dla, bpa, ta)

Incoming asymptote parameters for hyperbolic trajectories.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • rp: Periapsis radius.
    • Range: rp > 0.
  • c3: Characteristic energy. Specific energy at infinity (v∞²).
    • -∞ < c3 < ∞, c3 > 0 for hyperbolic trajectories.
  • rla: Right ascension of incoming asymptote (rad).
    • Range: [0, 2π).
  • dla: Declination of incoming asymptote (rad).
    • Range: [-π/2, π/2].
  • bpa: B-plane angle (rad).
    • Range: [0, 2π).
  • ta: True anomaly at asymptote (rad).
    • Range: [0, 2π).

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Singularities exist when trajectory becomes parabolic (c3 ≈ 0).
  • Singularities exist when e ≈ 0. Note for elliptic orbits the apsides vector is used in place of the asymptote direction.
  • Singular when asymptote direction is aligned with the reference z-axis (dla = ±π/2).

Examples

inasym = IncomingAsymptoteState(6778.0, 5.0, 0.0, π/4, π/2, π/2)
source
AstroStates.KeplerianStateType
KeplerianState(a, e, i, raan, aop, ta)

Keplerian orbital elements representation using classical osculating elements.

Units

  • Distance/time units must be consistent with the gravitational parameter μ used in conversions.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • sma: Semi-major axis (must be nonzero). Defines orbit size and energy.
    • If sma > 0: elliptic orbit (bound). If sma < 0: hyperbolic orbit (unbound).
    • For elliptic orbits, sma is half the major axis length.
  • ecc: Eccentricity. Defines orbit shape.
    • ecc = 0: circular orbit. 0 < ecc < 1: elliptical orbit.
    • ecc = 1: parabolic orbit. ecc > 1: hyperbolic orbit.
    • Valid range: [0, ∞), but ecc ≈ 1 results in infinite sma.
  • inc: Inclination (rad). Angle between orbit plane and reference xy-plane.
    • Range: [0, π]. If inc < π/2: prograde orbit. If inc > π/2: retrograde orbit.
  • raan: Right ascension of ascending node (rad). Orients the orbit plane.
    • Range: [0, 2π). Angle from +x axis to ascending node, measured in xy-plane.
    • Defines where orbit plane intersects reference plane (ascending crossing).
    • Undefined for equatorial orbits (inc ≈ 0 or π).
  • aop: Argument of periapsis (rad). Orients the orbit within its plane.
    • Range: [0, 2π). Angle from ascending node to periapsis point.
    • Defines orientation of orbit ellipse within the orbital plane.
    • Undefined for circular orbits (ecc ≈ 0).
  • ta: True anomaly (rad). Spacecraft position within the orbit.
    • Range: [0, 2π). Angle from periapsis to current spacecraft position.
    • ta = 0: at periapsis. ta = π: at apoapsis (for elliptic orbits).
    • Undefined for circular orbits (ecc ≈ 0).

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Classical Keplerian elements have well-known singularities for special cases.

Examples

k = KeplerianState(7000.0, 0.01, π/4, 0.0, 0.0, π/3)
source
AstroStates.ModifiedEquinoctialStateType
ModifiedEquinoctialState(p, f, g, h, k, L)

Modified equinoctial elements representation.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • p: Semi-latus rectum. Parameter defining orbit size and shape.
    • Range: p > 0. Related to semi-major axis and eccentricity: p = a(1-e²).
  • f: Eccentricity vector f-component. f = e⋅cos(ω + Ω).
    • Range: Any real value. Eccentricity component in perifocal frame.
  • g: Eccentricity vector g-component. g = e⋅sin(ω + Ω).
    • Range: Any real value. Eccentricity component in perifocal frame.
  • h: Inclination vector h-component. h = tan(i/2)⋅cos(Ω).
    • Range: Any real value. Inclination component related to ascending node.
  • k: Inclination vector k-component. k = tan(i/2)⋅sin(Ω).
    • Range: Any real value. Inclination component related to ascending node.
  • L: True longitude (rad). Ω + ω + ν combined angle measure.
    • Range: [0, 2π). True longitude from reference direction.

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Singularities reduced compared to classical elements, but still exist for retrograde equatorial orbits.

Examples

mee = ModifiedEquinoctialState(7000.0, 0.01, 0.0, 0.1, 0.0, π/4)
source
AstroStates.ModifiedKeplerianStateType
ModifiedKeplerianState(rp, ra, inc, raan, aop, ta)

Modified Keplerian orbital elements using periapsis and apoapsis radii instead of semi-major axis and eccentricity.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • rp: Radius of periapsis. Closest approach distance to central body.
    • Range: rp > 0. Typically rp > central body radius for physical orbits.
  • ra: Radius of apoapsis. Farthest distance from central body (elliptic orbits only).
    • Range: ra ≥ rp for elliptic orbits. For hyperbolic orbits: ra = ∞ (not used).
  • inc: Inclination (rad). Angle between orbit plane and reference xy-plane.
    • Range: [0, π]. If inc < π/2: prograde orbit. If inc > π/2: retrograde orbit.
    • inc = 0: equatorial orbit in xy-plane. inc = π/2: polar orbit.
  • raan: Right ascension of ascending node (rad). Orients the orbit plane.
    • Range: [0, 2π). Angle from +x axis to ascending node, measured in xy-plane.
    • Defines where orbit plane intersects reference plane (ascending crossing).
    • Undefined for equatorial orbits (inc ≈ 0 or π).
  • aop: Argument of periapsis (rad). Orients the orbit within its plane.
    • Range: [0, 2π). Angle from ascending node to periapsis point.
    • Defines orientation of orbit ellipse within the orbital plane.
    • Undefined for circular orbits (rp ≈ ra).
  • ta: True anomaly (rad). Spacecraft position within the orbit.
    • Range: [0, 2π). Angle from periapsis to current spacecraft position.
    • ta = 0: at periapsis. ta = π: at apoapsis (for elliptic orbits).
    • Undefined for circular orbits (rp ≈ ra).

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Alternative to classical Keplerian elements using radius parameters instead of sma/ecc.
  • Shares same singularities as classical Keplerian elements for circular and equatorial orbits.

Examples

# 400 km x 35,786 km orbit (GTO-like)
modkep = ModifiedKeplerianState(6778.0, 42164.0, π/6, 0.0, 0.0, 0.0)
source
AstroStates.OrbitStateType
OrbitState(state::AbstractVector, statetype::AbstractOrbitStateType)

A generic wrapper for orbital state vectors and their associated metadata.

Arguments

  • state: The numerical state vector (e.g., position/velocity, orbital elements). Must be a vector of real numbers.
  • statetype: Marker instance indicating the state representation (e.g., Cartesian(), Keplerian(), etc.).
  •          - to see all state types, use `subtypes(AbstractOrbitStateType)`.

Example

state_vec = [7000.0, 0.0, 0.0, 0.0, 7.5, 0.0]  # Cartesian position and velocity
os = OrbitState(state_vec, Cartesian())

Notes

  • Allows easy switching between different state representations without type instability.
  • The struct is internally parameterized for performance, and type safety, and differentiation, but users should construct it using the outer constructor as shown above.
source
AstroStates.OutGoingAsymptoteStateType
OutGoingAsymptoteState(rp, c3, rla, dla, bpa, ta)

Outgoing asymptote parameters for hyperbolic trajectories.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • rp: Periapsis radius.
    • Range: rp > 0.
  • c3: Characteristic energy. Specific energy at infinity (v∞²).
    • -∞ < c3 < ∞, c3 > 0 for hyperbolic trajectories.
  • rla: Right ascension of asymptote (rad).
    • Range: [0, 2π).
  • dla: Declination of asymptote (rad).
    • Range: [-π/2, π/2].
  • bpa: B-plane angle (rad).
    • Range: [0, 2π).
  • ta: True anomaly at asymptote (rad).
    • Range: [0, 2π).

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Singularities exist when trajectory becomes parabolic (c3 ≈ 0).
  • Singularities exist when e ≈ 0. Note for elliptic orbits the apsides vector is used in place of the asymptote direction.
  • Singular when asymptote direction is aligned with the reference z-axis (dla = ±π/2).

Examples

outasym = OutGoingAsymptoteState(6778.0, 5.0, 0.0, π/4, π/2, π/2)
source
AstroStates.SphericalAZIFPAStateType
SphericalAZIFPAState(r, ra, dec, v, vazi, fpa)

Spherical coordinates with azimuth and flight path angle velocity representation.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • r: Radial distance from central body center. Magnitude of position vector.
    • Range: r > 0.
  • ra: Right ascension (rad). Azimuthal angle of position in xy-plane.
    • Range: [0, 2π). Measured counterclockwise from +x axis.
  • dec: Declination (rad). Elevation angle of position above/below xy-plane.
    • Range: [-π/2, π/2]. dec = 0: position in xy-plane. dec = ±π/2: at poles.
  • v: Magnitude of velocity vector.
    • Range: v ≥ 0.
  • vazi: Velocity azimuth (rad). Horizontal direction of velocity vector.
    • Range: [0, 2π). Angle east of north in local horizontal plane.
    • vazi = 0: northward velocity. vazi = π/2: eastward velocity.
  • fpa: Flight path angle (rad). Elevation of velocity above local horizontal.
    • Range: [-π/2, π/2]. fpa = 0: horizontal flight. fpa > 0: climbing.
    • fpa < 0: descending flight. fpa = ±π/2: purely radial motion.

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Singularities exist at r = 0 or v = 0 due to undefined angles.
  • Flight path angle and azimuth provide intuitive velocity direction description.

Examples

sph = SphericalAZIFPAState(6478.0, 0.0, π/4, 7.5, π/4, π/4)
source
AstroStates.SphericalRADECStateType
SphericalRADECState(r, dec, ra, v, decv, rav)

Spherical coordinates state with right ascension and declination components.

Units

  • Distance and time units must be consistent with the gravitational parameter μ used in the simulation.
  • All angular quantities are in radians.

Fields (all ::T where T<:Real)

  • r: Radial distance from central body center.
    • Range: r > 0.
  • dec: Declination (rad). Elevation angle of position above/below xy-plane.
    • Range: [-π/2, π/2]. dec = 0: equatorial. dec = ±π/2: polar.
  • ra: Right ascension (rad). Azimuthal angle of position in xy-plane.
    • Range: [0, 2π). Measured counterclockwise from +x axis.
  • v: Velocity magnitude.
    • Range: v ≥ 0.
  • decv: Declination of velocity (rad). Elevation angle of velocity vector.
    • Range: [-π/2, π/2]. Angle between velocity vector and xy-plane.
  • rav: Right ascension of velocity (rad). Azimuthal angle of velocity.
    • Range: [0, 2π). Direction of velocity in xy-plane.

Notes

  • Parametric so automatic differentiation and high-precision types are supported.
  • Singularities exist at r = 0 or v = 0 due to undefined angles.

Examples

radec = SphericalRADECState(6778.0, π/4, 0.0, 7.5, 0.0, π/2)
source
AstroStates.alt_equinoctial_to_equinoctialMethod
alt_equinoctial_to_equinoctial(alt::Vector{<:Real}; tol::Float64 = 1e-12)

Convert alternate equinoctial elements to standard equinoctial elements.

Arguments

  • alt::Vector{<:Real}: Alternate equinoctial state [a, h, k, altp, altq, λ]
    • a : semi-major axis [length]
    • h : e⋅g component of eccentricity vector
    • k : e⋅f component of eccentricity vector
    • altp : sin(i/2)⋅cos(Ω)
    • altq : sin(i/2)⋅sin(Ω)
    • λ : mean longitude [rad]
  • tol::Float64: tolerance for inclination singularity (default = 1e-12)

Returns

Standard equinoctial state [a, h, k, p, q, λ] - a : semi-major axis [length] - h : e⋅g component of eccentricity vector - k : e⋅f component of eccentricity vector - p : tan(i/2)⋅cos(Ω) - q : tan(i/2)⋅sin(Ω) - λ : mean longitude [rad]

Notes

  • Fails if inclination approaches 180°.
  • All angles in radians. Units consistent with μ.

Examples

alt_eq = [7000.0, 0.01, 0.0, 0.05, 0.0, π/4]
std_eq = alt_equinoctial_to_equinoctial(alt_eq)
source
AstroStates.cart_to_equinoctialMethod
cart_to_equinoctial(cart::Vector{<:Real}, μ::Real; tol::Float64 = 1e-12)

Convert a Cartesian state vector to equinoctial orbital elements.

Arguments

  • cart::Vector{<:Real}: Cartesian state [x, y, z, vx, vy, vz]
  • μ::Real: gravitational parameter [length³/time²]
  • tol::Float64: tolerance for singularity detection (default = 1e-12)

Returns

Equinoctial state [a, h, k, p, q, λ]:

  • a : semi-major axis [length]
  • h : e⋅g component of eccentricity vector
  • k : e⋅f component of eccentricity vector
  • p : tan(i/2)⋅cos(Ω)
  • q : tan(i/2)⋅sin(Ω)
  • λ : mean longitude [rad]

Notes

  • Fails gracefully with NaN if orbit is hyperbolic, singular, or edge-case unsupported.
  • All angles are in radians. Units consistent with μ.
  • Note that in most cases in states, h is the magnitude of angular momentum. But, not for equinoctial elements.

Examples

cart = [6778.0, 0.0, 0.0, 0.0, 7.66, 0.0]
equinoctial = cart_to_equinoctial(cart, 398600.4418)
source
AstroStates.cart_to_inasymptoteMethod
cart_to_inasymptote(cart::Vector{<:Real}, μ::Real; tol::Float64 = 1e-12)

Convert a Cartesian state to incoming asymptote parameters.

Arguments

  • cart::Vector{<:Real}: Cartesian state [x, y, z, vx, vy, vz]
  • μ::Real: gravitational parameter
  • tol::Real: tolerance for detecting singularities

Returns

A 6-element vector: [rₚ, C₃, λₐ, δₐ, θᵦ, ν]

  • rₚ : radius of periapsis
  • C₃ : characteristic energy
  • λₐ : right ascension of incoming asymptote
  • δₐ : declination of incoming asymptote
  • θᵦ : B-plane angle
  • ν : true anomaly

Notes

  • Angles are in radians.
  • Dimensional quantities are consistent units with μ.

Examples

cart = [10000.0, 0.0, 0.0, 0.0, 12.0, 0.0]  # Hyperbolic trajectory
inasym = cart_to_inasymptote(cart, 398600.4418)
source
AstroStates.cart_to_kepMethod
function cart_to_kep(cart::Vector{<:Real}, μ::Real; tol::Float64=1e-12)

Convert a Cartesian state vector to Keplerian orbital elements.

Arguments

  • cart::Vector{Float64}: Cartesian state [x, y, z, vx, vy, vz]
  • μ::Real: Gravitational parameter
  • tol::Float64: Tolerance for detecting circular or equatorial orbits (default: 1e-12)

Returns

A vector [a, e, i, Ω, ω, ν] where:

  • a: semi-major axis
  • e: eccentricity
  • i: inclination
  • Ω: right ascension of ascending node (RAAN)
  • ω: argument of periapsis
  • ν: true anomaly

Notes

  • Angles must be in radians.
  • Dimensional quantities must use consistent units with μ.

Examples

cart = [6778.0, 0.0, 0.0, 0.0, 7.66, 0.0]
kep = cart_to_kep(cart, 398600.4418)
source
AstroStates.cart_to_meeMethod
cart_to_mee(cart::Vector{<:Real}, μ::Real; j::Float64 = 1.0)

Convert Cartesian state to Modified Equinoctial Elements (MEE).

Arguments

  • cart::Vector{<:Real}: 6-element vector [x, y, z, vx, vy, vz]
  • μ::Real: Gravitational parameter
  • j::Float64=1.0: Optional constant (1 for prograde, -1 for retrograde), defaults to 1.0
  • tol::Float64: Optional tolerance for singularity checking

Returns

  • A 6-element vector [p, f, g, h, k, L] representing the modified equinoctial elements.

Examples

cart = [6778.0, 0.0, 0.0, 0.0, 7.66, 0.0]
mee = cart_to_mee(cart, 398600.4418)
source
AstroStates.cart_to_outasymptoteMethod
cart_to_outasymptote(cart::Vector{<:Real}, μ::Real; tol::Float64 = 1e-12) -> Vector{Float64}

Convert a Cartesian state to outgoing asymptote parameters.

Arguments

  • cart::Vector{<:Real}: Cartesian state [x, y, z, vx, vy, vz]
  • μ::Real: gravitational parameter
  • tol::Float64: tolerance for detecting singularities

Returns

A 6-element vector: [rₚ, C₃, λₐ, δₐ, θᵦ, ν]

  • rₚ : radius of periapsis
  • C₃ : characteristic energy
  • λₐ : right ascension of outgoing asymptote
  • δₐ : declination of outgoing asymptote
  • θᵦ : B-plane angle
  • ν : true anomaly

Notes

  • Angles are in radians.
  • Dimensional quantities are consistent with units of μ

Examples

cart = [10000.0, 0.0, 0.0, 0.0, 12.0, 0.0]  # Hyperbolic trajectory
outasym = cart_to_outasymptote(cart, 398600.4418)
source
AstroStates.cart_to_sphazfpaMethod
cart_to_sphazfpa(cart::Vector{<:Real}; tol::Float64 = 1e-12)

Convert a Cartesian state to Spherical AZ-FPA representation.

Arguments

  • cart::Vector{<:Real}: Cartesian state [x, y, z, vx, vy, vz]
  • tol::Float64: Numerical tolerance for singularity checks (default: 1e-12)

Returns

A 6-element Spherical AZ-FPA state [r, λ, δ, v, αₚ, ψ]:

  • r : radial distance [length]
  • λ : right ascension [rad]
  • δ : declination [rad]
  • v : velocity magnitude [length/time]
  • αₚ : flight path azimuth [rad]
  • ψ : flight path angle [rad]

Notes

  • Returns fill(NaN, 6) if r or v are near zero or orbit is singular
  • All angles are in radians.

Examples

cart = [6778.0, 0.0, 0.0, 0.0, 7.66, 0.0]
sphazfpa = cart_to_sphazfpa(cart)
source
AstroStates.cart_to_sphradecMethod
cart_to_sphradec(state::Vector{<:Real}; tol::Float64=1e-12)

Convert a Cartesian state vector to a spherical RA/DEC state vector.

Arguments

  • state::Vector{<:Real}: A 6-element vector [x, y, z, vx, vy, vz] representing Cartesian position and velocity.

Returns

A 6-element vector [r, ra, dec, v, vra, vdec] where:

  • r : magnitude of position vector
  • λᵣ : right ascension (radians)
  • δᵣ : declination (radians)
  • v : magnitude of velocity
  • λᵥ : azimuthal direction of velocity (radians)
  • δᵥ : elevation direction of velocity (radians)

Notes

  • Assumes all angles are in radians.
  • Units must be consistent between position and velocity components.

Examples

cart = [6778.0, 0.0, 0.0, 0.0, 7.66, 0.0]
sphradec = cart_to_sphradec(cart)
source
AstroStates.equinoctial_to_alt_equinoctialMethod
equinoctial_to_alt_equinoctial(eq::Vector{<:Real}; tol::Float64 = 1e-12)

Convert equinoctial elements to alternate equinoctial elements.

Arguments

  • eq::Vector{<:Real}: Equinoctial state [a, h, k, p, q, λ]
    • a : semi-major axis [length]
    • h : e⋅g component of eccentricity vector
    • k : e⋅f component of eccentricity vector
    • p : tan(i/2)⋅cos(Ω)
    • q : tan(i/2)⋅sin(Ω)
    • λ : mean longitude [rad]
  • tol::Float64: tolerance for inclination singularity detection (default: 1e-12)

Returns

Alternate equinoctial state [a, h, k, altp, altq, λ] - a : semi-major axis [length] - h : e⋅g component of eccentricity vector - k : e⋅f component of eccentricity vector - altp : sin(i/2)⋅cos(Ω) - altq : sin(i/2)⋅sin(Ω) - λ : mean longitude [rad]

Notes

  • Singular when inclination is near 180°, returns fill(NaN, 6) in that case.
  • Units and angles must be consistent; λ in radians.

Examples

eq = [7000.0, 0.01, 0.0, 0.1, 0.0, π/4]
alt_eq = equinoctial_to_alt_equinoctial(eq)
source
AstroStates.equinoctial_to_cartMethod
equinoctial_to_cart(eq::Vector{<:Real}, μ::Real; tol::Float64 = 1e-12)

Convert equinoctial elements to Cartesian state.

Arguments

  • eq::Vector{<:Real}: Equinoctial state [a, h, k, p, q, λ]
    • a : semi-major axis [length]
    • h : e⋅g component of eccentricity vector
    • k : e⋅f component of eccentricity vector
    • p : tan(i/2)⋅cos(Ω)
    • q : tan(i/2)⋅sin(Ω)
    • λ : mean longitude [rad]
  • μ::Real: gravitational parameter [length³/time²]
  • tol::Float64: numerical tolerance (default: 1e-12)

Returns

Cartesian state [x, y, z, vx, vy, vz]

Notes

  • Returns fill(NaN, 6) if eccentricity is too high or radius becomes non-physical.
  • Assumes all angles are in radians and other units are consistent with μ.

Examples

eq = [7000.0, 0.01, 0.0, 0.1, 0.0, π/4]
cart = equinoctial_to_cart(eq, 398600.4418)
source
AstroStates.inasymptote_to_kepMethod
inasymptote_to_kep(outasym::Vector{<:Real}, μ::Real; tol::Float64=1e-12)

Convert incoming asymptote elements to Keplerian elements.

Arguments

  • outasym: Vector{<:Real} of outgoing asymptote elements:

    • rₚ : periapsis radius [length]
    • C₃ : characteristic energy [length²/time²]
    • λₐ : right ascension of the asymptote [rad]
    • δₐ : declination of the asymptote [rad]
    • θᵦ : B-plane angle [rad]
    • ν : true anomaly [rad]
  • μ: Gravitational parameter [length³/time²]

  • tol: Singularity tolerance (default = 1e-12)

Returns

  • Keplerian state vector [a, e, i, Ω, ω, ν]

Notes

  • Returns fill(NaN, 6) if singularity is detected.
  • Angles in radians. Units consistent with μ.

Examples

inasym = [6778.0, 5.0, 0.0, π/4, π/2, π/2]
kep = inasymptote_to_kep(inasym, 398600.4418)
source
AstroStates.kep_to_cartMethod
kep_to_cart(state::Vector{<:Real}, μ::Real; tol::Float64=1e-12)

Convert a Keplerian state vector to a Cartesian state vector.

Arguments

  • state::Vector{<:Real}: Keplerian elements [a, e, i, Ω, ω, ν]
  • μ: Gravitational parameter
  • tol: Tolerance for singularities like p ≈ 0 (default: 1e-12)
  • a: semi-major axis
  • e: eccentricity
  • i: inclination
  • Ω: right ascension of ascending node
  • ω: argument of periapsis
  • ν: true anomaly

Returns

A 6-element vector [x, y, z, vx, vy, vz] representing Cartesian position and velocity.

Examples

kep = [7000.0, 0.01, π/4, 0.0, 0.0, π/3]
cart = kep_to_cart(kep, 398600.4418)

Notes

  • Angles must be in radians.
  • Dimensional quantities must be consistent units with μ.
  • Returns a vector of NaNs if conversion is undefined.
source
AstroStates.kep_to_modkepMethod
kep_to_modkep(kep::Vector{<:Real}; tol::Float64 = 1e-12)

Convert a classical Keplerian state to a Modified Keplerian state.

Arguments

  • kep::Vector{<:Real}: Keplerian state vector [a, e, i, Ω, ω, ν]

    • a : semi-major axis
    • e : eccentricity
    • i : inclination
    • Ω : right ascension of ascending node
    • ω : argument of periapsis
    • ν : true anomaly
  • tol::Float64: tolerance for singularity and consistency checks (default = 1e-12)

Returns

Modified Keplerian state [rₚ, rₐ, i, Ω, ω, ν] or fill(NaN, 6) if invalid.

Notes

  • Parabolic orbits (e ≈ 1) and singular conics are not supported.
  • Units must be consistent. Angles in radians.

Examples

kep = [7000.0, 0.01, π/4, 0.0, 0.0, π/3]
modkep = kep_to_modkep(kep)
source
AstroStates.mee_to_cartMethod
mee_to_cart(mod_equinoct::Vector{<:Real}, μ::Real; j::Real = 1.0)

Convert Modified Equinoctial Elements to Cartesian state.

Arguments

  • mod_equinoct::Vector{<:Real}: Vector containing the Modified Equinoctial Elements [p, f, g, h, k, L]
  • μ::Real: Gravitational parameter
  • j::Real=1.0: Optional constant (1 for prograde, -1 for retrograde), defaults to 1.0

Returns

  • A 6-element vector [x, y, z, vx, vy, vz] representing the Cartesian position and velocity.

Examples

mee = [7000.0, 0.01, 0.0, 0.1, 0.0, π/4]
cart = mee_to_cart(mee, 398600.4418)
source
AstroStates.modkep_to_kepMethod
modkep_to_kep(modkep::Vector{<:Real}; tol::Float64 = 1e-12)

Convert a Modified Keplerian state to a classical Keplerian state.

Arguments

  • modkep::Vector{<:Real}: Modified Keplerian state vector [rₚ, rₐ, i, Ω, ω, ν]

    • rₚ : radius of periapsis
    • rₐ : radius of apoapsis
    • i : inclination
    • Ω : right ascension of ascending node
    • ω : argument of periapsis
    • ν : true anomaly
  • tol::Float64: tolerance for singularity and consistency checks (default = 1e-12)

Returns

A classical Keplerian state [a, e, i, Ω, ω, ν] or fill(NaN, 6) if invalid.

Notes

  • Returns NaN vector if input describes a parabolic, undefined, or inconsistent orbit.
  • Assumes input angles are in radians and distances use consistent units.

Examples

modkep = [6778.0, 42164.0, π/6, 0.0, 0.0, 0.0]
kep = modkep_to_kep(modkep)
source
AstroStates.outasymptote_to_kepMethod
outasymptote_to_kep(outasym::Vector{<:Real}, μ::Real; tol::Float64=1e-12)

Convert outgoing asymptote elements to Keplerian elements.

Arguments

  • outasym: Vector{<:Real} of outgoing asymptote elements:

    • rₚ : periapsis radius [length]
    • C₃ : characteristic energy [length²/time²]
    • λₐ : right ascension of the asymptote [rad]
    • δₐ : declination of the asymptote [rad]
    • θᵦ : B-plane angle [rad]
    • ν : true anomaly [rad]
  • μ: Gravitational parameter [length³/time²]

  • tol: Singularity tolerance (default = 1e-12)

Returns

  • Keplerian state vector [a, e, i, Ω, ω, ν]

Notes

  • Returns fill(NaN, 6) if singularity is detected.
  • Angles in radians. Units consistent with μ.

Examples

outasym = [6778.0, 5.0, 0.0, π/4, π/2, π/2]
kep = outasymptote_to_kep(outasym, 398600.4418)
source
AstroStates.sphazfpa_to_cartMethod
sphazfpa_to_cart(spherical::Vector{<:Real}) -> Vector{<:Real}

Convert a Spherical AZ-FPA state to a Cartesian state vector.

Arguments

  • spherical::Vector{<:Real}: Spherical AZ-FPA state vector [r, λ, δ, v, αₚ, ψ]
    • r : radial distance [length]
    • λ : right ascension [rad]
    • δ : declination [rad]
    • v : velocity magnitude [length/time]
    • αₚ : flight path azimuth (angle east of north in local horizon) [rad]
    • ψ : flight path angle (angle above local horizon) [rad]

Returns

A 6-element Cartesian state vector [x, y, z, vx, vy, vz].

Notes

  • All angles must be in radians.
  • Velocity frame uses local vertical/local horizontal.

Examples

sphazfpa = [6478.0, 0.0, π/4, 7.5, π/4, π/4]
cart = sphazfpa_to_cart(sphazfpa)
source
AstroStates.sphradec_to_cartMethod
sphradec_to_cart(state::Vector{<:Real}) -> Vector{<:Real}

Convert a spherical RA/Dec state [r, λᵣ, δᵣ, v, λᵥ, δᵥ] to a Cartesian state [x, y, z, vx, vy, vz].

Arguments

  • state::Vector{<:Real}: length-6 vector where:
    • r = position magnitude
    • λᵣ = right ascension of position (radians)
    • δᵣ = declination of position (radians)
    • v = velocity magnitude
    • λᵥ = azimuth of velocity direction (radians)
    • δᵥ = elevation of velocity direction (radians)

Returns

  • Vector{<:Real}: length-6 Cartesian state [x, y, z, vx, vy, vz].

Notes

  • All angles are in radians.
  • Units must be consistent between position and velocity components.
  • Right ascension typically in [0, 2π); declination in [-π/2, π/2].

Examples

sphradec = [7000.0, 0.0, 0.0, 7.5, π/2, 0.0]
cart = sphradec_to_cart(sphradec)
source
AstroStates.to_stateMethod
to_state(os::OrbitState) -> AbstractState

Convert an OrbitState to its corresponding concrete state type (subtype of AbstractState), using the stored state vector and marker type.

Arguments

  • os::OrbitState: The OrbitState instance to convert.

Returns

  • The concrete state type (e.g., CartesianState, KeplerianState, etc.) corresponding to the statetype marker and the stored state vector.

Example

state_vec = [7000.0, 0.0, 0.0, 0.0, 7.5, 0.0]
os = OrbitState(state_vec, Cartesian())
cs = to_state(os)  # cs is a CartesianState

# For other state types:
os_kep = OrbitState([8000.0, 0.05, 0.1, 0.2, 0.3, 0.4], Keplerian())
ks = to_state(os_kep)  # ks is a KeplerianState

Notes

  • This function does not perform any physical conversions; it simply reconstructs the concrete state type from the stored vector and marker.
  • For conversions between different state types (e.g., CartesianState to KeplerianState), use the appropriate constructor or conversion function, providing additional parameters (such as μ) if required.
source

Index