Library

Documentation for AstroCoords.jl.

AstroCoords.AstroCoordTransformationType
abstract type AstroCoordTransformation <: AstrodynamicsTransformation

An abstract type representing a Transformation Between Astrodynamics Coordinates with no Time Regularization

source
AstroCoords.CartesianType
Cartesian{T} <: AstroCoord

Cartesian Orbital Elements. 6D parameterziation of the orbit. x - X-position y - Y-position z - Z-position ẋ - X-velocity ẏ - Y-velocity ż - Z-velocity

Constructors Cartesian(x, y, z, ẋ, ẏ, ż) Cartesian(X::AbstractArray) Cartesian(X::AstroCoord, μ::Number)

source
AstroCoords.CoordinateType
abstract type AstroCoord{N, T} <: StaticMatrix{N, 1, T}

An abstract type representing a N-Dimensional Coordinate Set

source
AstroCoords.CylindricalType
Cylindrical{T} <: AstroCoord

Cylindrical Orbital Elements. 6D parameterziation of the orbit ρ - in-plane radius θ - in-plane angle z - out-of-plane distance ρdot - instantaneous rate of change of in-plsne radius θdot - instantaneous rate of change of θ ż - instantaneous rate of change of z

Constructors Cylindrical(r, θ, z, ṙ, θdot, ż) Cylindrical(X::AbstractArray) Cylindrical(X::AstroCoord, μ::Number)

source
AstroCoords.DelaunayType
Delaunay{T} <: AstroCoord

Delaunay Orbital Elements. 6D parameterziation of the orbit. L - Canonical Keplerian Energy G - Canonical Total Angular Momentum H - Canonical Normal Angular Momentum (Relative to Equator) M - Mean Anomaly ω - Argument of Periapsis Ω - Right Ascention of the Ascending Node

Constructors Delaunay(L, G, H, M, ω, Ω) Delaunay(X::AbstractVector{<:Number}) Delaunay(X::AstroCoord, μ::Number)

source
AstroCoords.J2EqOEType
J2EqOE{T} <: AstroCoord

Modified Equinoctial Orbital Elements. 6D parameterziation of the orbit. n - mean motion h - eccetricity projection onto longitude of perigee k - eccetricity projection onto ⟂ longitude of perigee p - projection of half inclination onto RAAN q - projection of half inclination onto ⟂ RAAN L - true longitude

Constructors J2EqOE(n, h, k, p, q, L) J2EqOE(X::AbstractArray) J2EqOE(X::AstroCoord, μ::Number)

source
AstroCoords.KeplerianType
Keplerian{T} <: AstroCoord

Keplerian Orbital Elements. 6D parameterziation of the orbit. a - semi-major axis e - eccetricity i - inclination Ω - Right Ascension of Ascending Node ω - Argument of Perigee f - True Anomaly

Constructors Keplerian(a, e, i, Ω, ω, f) Keplerian(X::AbstractArray) Keplerian(X::AstroCoord, μ::Number)

source
AstroCoords.MilankovichType
Milankovich{T} <: AstroCoord

Milankovich Orbital Elements. 7D parameterziation of the orbit. hx - X-component of Angular Momentum Vector hy - Y-component of Angular Momentum Vector hz - Z-component of Angular Momentum Vector ex - X-component of Eccentricity Vector ey - Y-component of Eccentricity Vector ez - Z-component of Eccentricity Vector L - True Longitude

Constructors Milankovich(hx, hy, hz, ex, ey, ez, L) Milankovich(X::AbstractArray) Milankovich(X::AstroCoord, μ::Number)

source
AstroCoords.ModEqType
ModEq{T} <: AstroCoord

Modified Equinoctial Orbital Elements. 6D parameterziation of the orbit. p - semi-parameter f - eccetricity projection onto longitude of perigee g - eccetricity projection onto ⟂ longitude of perigee h - projection of half inclination onto RAAN k - projection of half inclination onto ⟂ RAAN L - true longitude

Constructors ModEq(p, f, g, h, k, l) ModEq(X::AbstractArray) ModEq(X::AstroCoord, μ::Number)

source
AstroCoords.SphericalType
Spherical{T} <: AstroCoord

Spherical Orbital Elements. 6D parameterziation of the orbit r - radius θ - in-plane angle ϕ - out-of-plane angle ṙ - instantaneous rate of change of radius θdot - instantaneous rate of change of θ ϕdot - instantaneous rate of change of ϕ

Constructors Spherical(r, θ, ϕ, ṙ, ω, Ω) Spherical(X::AbstractArray) Spherical(X::AstroCoord, μ::Number)

source
AstroCoords.TransformationType

The Transformation supertype defines a simple interface for performing transformations. Subtypes should be able to apply a coordinate system transformation on the correct data types by overloading the call method, and usually would have the corresponding inverse transformation defined by Base.inv(). Efficient compositions can optionally be defined by compose() (equivalently ).

source
AstroCoords.USM6Type
USM6{T} <: AstroCoord

Unified State Model Orbital Elements. 6D parameterziation of the orbit using Velocity Hodograph and MRP's C - Velocity Hodograph Component Normal to the Radial Vector Laying in the Orbital Plane Rf1 - Velocity Hodograph Component 90 degrees ahead of the Eccentricity Vector - Along the Intermediate Rotating Frame X-Axis Rf2 - Velocity Hodograph Component 90 degrees ahead of the Eccentricity Vector - Along the Intermediate Rotating Frame Y-Axis σ1 - First Modified Rodriguez Parameter σ2 - Second Modified Rodriguez Parameter σ3 - Third Modified Rodriguez Parameter

Constructors USM6(C, Rf1, Rf2, σ1, σ2, σ3) USM6(X::AbstractArray) USM6(X::AstroCoord, μ::Number)

source
AstroCoords.USM7Type
USM7{T} <: AstroCoord

Unified State Model Orbital Elements. 7D parameterziation of the orbit using Velocity Hodograph and Quaternions C - Velocity Hodograph Component Normal to the Radial Vector Laying in the Orbital Plane Rf1 - Velocity Hodograph Component 90 degrees ahead of the Eccentricity Vector - Along the Intermediate Rotating Frame X-Axis Rf2 - Velocity Hodograph Component 90 degrees ahead of the Eccentricity Vector - Along the Intermediate Rotating Frame Y-Axis ϵO1 - First Imaginary Quaternion Component ϵO2 - Second Imaginary Quaternion Component ϵO3 - Third Imaginary Quaternion Component η0 - Real Quaternion Component

Constructors USM7(C, Rf1, Rf2, ϵO1, ϵO2, ϵO3, η0) USM7(X::AbstractArray) USM7(X::AstroCoord, μ::Number)

source
AstroCoords.USMEMType
USMEM{T} <: AstroCoord

Unified State Model Orbital Elements. 6D parameterziation of the orbit using Velocity Hodograph and Exponential Mapping C - Velocity Hodograph Component Normal to the Radial Vector Laying in the Orbital Plane Rf1 - Velocity Hodograph Component 90 degrees ahead of the Eccentricity Vector - Along the Intermediate Rotating Frame X-Axis Rf2 - Velocity Hodograph Component 90 degrees ahead of the Eccentricity Vector - Along the Intermediate Rotating Frame Y-Axis a1 - First Exponential Mapping Component a2 - Second Exponential Mapping Component a3 - Third Exponential Mapping Component

Constructors USMEM(C, Rf1, Rf2, a1, a2, a3) USMEM(X::AbstractArray) USMEM(X::AstroCoord, μ::Number)

source
AstroCoords.EP2MRPMethod
EP2MRP(β::AbstractVector{<:Number})

Converts Euler Parameter rotation description into Modified Rodriguez Parameters.

Arguments

  • β::AbstractVector{<:Number}: The Euler Parameter description of a rotation.

Returns

  • σ::AbstractVector{<:Number}: The Modified Rodriguez Parameter description of a rotation.
source
AstroCoords.IOE2J2IOEMethod
function IOE2J2IOE(u::IOE{T}, μ::V) where {T<:Number, V<:Number}

Computes the J2 Perturbed Intermediate Orbit Elements from a Intermediate Orbit Element set.

Arguments

-u::AbstractVector{<:Number}: The Intermediate Orbit Element vector [I1; I2; I3; I4; I5; I6]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_J2IOR::SVector{6, <:Number}: The J2 Perturbed Intermediate Orbit Element vector [a; e; i; Ω(RAAN); ω(AOP); f(True Anomaly)].

source
AstroCoords.IOE2koeMethod
function IOE2koeM(u::AbstractVector{T}, μ::V) where {T<:Number, V<:Number}

Computes the Keplerian orbital elements from a Intermediate Orbit Element set.

Arguments

-u::AbstractVector{<:Number}: The Intermediate Orbit Element vector [I1; I2; I3; I4; I5; I6]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_koeM::SVector{6, <:Number}: The Keplerian state vector [a; e; i; Ω(RAAN); ω(AOP); M(Mean Anomaly)].

source
AstroCoords.IOE2modEqNMethod
function IOE2modEq(u::AbstractVector{T}, μ::V) where {T<:Number, V<:Number}

Computes the Modified Equinoctial orbital elements from a Intermediate Orbit Element set.

Arguments

-u::AbstractVector{<:Number}: The Intermediate Orbit Element vector [I1; I2; I3; I4; I5; I6]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_modEq::SVector{6, <:Number}: The Modified Equinoctial state vector [n; f; g; h; k; L].

source
AstroCoords.J2EqOE2cartMethod
function J2EqOE2cart(u::AbstractVector{<:Number}, μ::Number)

Computes the Cartesian state vector from a J2 Perturbed Equinoctial Orbit Elements.

Arguments

-u::AbstractVector{<:Number}: The J2 Perturbed Equinoctial Orbit Element vector [n; h; k; p; q; L]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_cart::SVector{6, <:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż].

source
AstroCoords.J2IOE2IOEMethod
function J2IOE2IOE(J2::J2IOE{T}, μ::V) where {T<:Number, V<:Number}

Computes the Intermediate Orbit Elements from a J2 Perturbed Intermediate Orbit Element set.

Arguments

-J2::J2IOE{<:Number}: The J2 Perturbed Intermediate Orbit Element vector [a; e; i; Ω(RAAN); ω(AOP); f(True Anomaly)]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_IOE::SVector{6, <:Number}: The Intermediate Orbit Element vector [I1; I2; I3; I4; I5; I6].

source
AstroCoords.KeplerSolverMethod
KeplerSolver(M::T, e::Number; tol::Float64=10 * eps(T)) where {T<:Number}

Solves for true anomaly given the mean anomaly and eccentricity of an orbit.

Arguments

-M::Number: Mean Anomaly of the orbit [radians]. -e::Number: Eccentricity of the orbit.

Keyword Arguments

-tol::Float64: Convergence tolerance of Kepler solver. [Default=10*eps(T)]

Returns

-f::Number`: True Anomaly of the orbit [radians]

source
AstroCoords.MRP2EPMethod
MRP2EP(σ::AbstractVector{<:Number})

Converts Modified Rodriguez Parameters rotation description into Euler Parameter.

Arguments

  • σ::AbstractVector{<:Number}: The Modified Rodriguez Parameter description of a rotation.

Returns

  • β::AbstractVector{<:Number}: The Euler Parameter description of a rotation.
source
AstroCoords.Mil2cartMethod

Mil2cart(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Converts Milankovich state vector into the Cartesian state vector.

Arguments

-u::AbstractVector{<:Number}: The Milankovich state vector [H; e; L]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_cart::SVector{6, <:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż].

source
AstroCoords.ModEq2koeMethod
ModEq2koe(u::AbstractVector{T}, μ::Number) where {T<:Number}

Converts Modified Equinoctial elements into the Keplerian elements.

Note

All angles are in radians.

Arguments

-u:AbstractVector{<:Number}: The Modified Equinoctial state vector [p; f; g; h; k; l]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_koe::SVector{6, <:Number}: The Keplerian state vector [a; e; i; Ω(RAAN); ω(AOP); ν(True Anomaly)].

source
AstroCoords.USM62USM7Method
USM62USM7(u::AbstractVector{T}, μ::Number) where {T<:Number}

Converts USM with Modified Rodrigue Parameters to USM with quaternions.

Arguments

-u::AbstractVector{<:Number}: The USM6 vector [C; Rf1; Rf2; σ1; σ2; σ3]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_USM::SVector{6, <:Number}: The Unified State Model vector [C; Rf1; Rf2; ϵO1; ϵO2; ϵO3; η0].

source
AstroCoords.USM72USM6Method
USM72USM6(u::AbstractVector{T}, μ::Number) where {T<:Number}

Converts USM with quaternions to USM with Modified Rodrigue Parameters.

Arguments

-u::AbstractVector{<:Number}: The Unified State Model vector [C; Rf1; Rf2; ϵO1; ϵO2; ϵO3; η0]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_USM6::SVector{6, <:Number}: The USM6 State vector [C; Rf1; Rf2; σ1; σ2; σ3].

source
AstroCoords.USM72USMEMMethod
USM72USMEM(u::AbstractVector{T}, μ::Number) where {T<:Number}

Converts USM with quaternions to USM with exponential mapping.

Arguments

-u::AbstractVector{<:Number}: The Unified State Model vector [C; Rf1; Rf2; ϵO1; ϵO2; ϵO3; η0]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_USMEM::SVector{6, <:Number}: The USMEM State Vector [C; Rf1; Rf2; a1; a2; a3, Φ].

source
AstroCoords.USM72koeMethod
USM72koe(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Converts Unified State Model elements into the Keplerian orbital element set. Van den Broeck, Michael. "An Approach to Generalizing Taylor Series Integration for Low-Thrust Trajectories." (2017). https://repository.tudelft.nl/islandora/object/uuid%3A2567c152-ab56-4323-bcfa-b076343664f9

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The Unified State Model vector [C; Rf1; Rf2; ϵO1; ϵO2; ϵO3; η0]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_koe:SVector{6, <:Number}: The Keplerian State vector [a; e; i; Ω(RAAN); ω(AOP); f(True Anomaly)].

source
AstroCoords.USMEM2USM7Method
USMEM2USM7(u::AbstractVector{T}, μ::Number) where {T<:Number}

Converts USM with exponential mapping to USM with quaternions.

Arguments

-u::AbstractVector{<:Number}: The USMEM vector [C; Rf1; Rf2; a1; a2; a3, Φ]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_USM::SVector{7, <:Number}: The Unified State Model vector [C; Rf1; Rf2; ϵO1; ϵO2; ϵO3; η0].

source
AstroCoords.angle_between_vectorsMethod
angle_between_vectors(
    v1::AbstractVector{T1}, v2::AbstractVector{T2}
) where {T1<:Number,T2<:Number}

Computes the angle between two vectors in a more numerically stable way than dot product.

Arguments

-v1::AbstractVector{<:Number}: The first vector of the computation -v2::AbstractVector{<:Number}: The second vector of the computation

Returns

-angle::Number: The angle between the two vectors

source
AstroCoords.angularMomentumQuantityMethod
angularMomentumQuantity(u::AbstractVector{<:Number})

Computes the instantaneous angular momentum.

Arguments

-u::AbstractVector{<:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż].

Returns

-angular_momentum::Number: Angular momentum of the body.

source
AstroCoords.angularMomentumQuantityMethod
angularMomentumQuantity(X::AstroCoord, μ::Number)

Computes the instantaneous angular momentum.

Arguments

-X::AstroCoord: An coordinate set describing the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

-angular_momentum::Number: Angular momentum of the body.

source
AstroCoords.angularMomentumVectorMethod
angularMomentumVector(u::AbstractVector{<:Number})

Computes the instantaneous angular momentum vector from a Cartesian state vector.

Arguments

-u::AbstractVector{<:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż].

Returns

-'angular_momentum::Vector{<:Number}': 3-Dimensional angular momemtum vector.

source
AstroCoords.angularMomentumVectorMethod
angularMomentumVector(X::AstroCoord, μ::Number)

Computes the instantaneous angular momentum vector from a Cartesian state vector.

Arguments

-X::AstroCoord: An coordinate set describing the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

-angular_momentum::Vector{<:Number}: 3-Dimensional angular momemtum vector.

source
AstroCoords.cart2J2EqOEMethod
function cart2J2EqOE(u::AbstractVector{<:Number}, μ::Number)

Computes the J2 Perturbed Equinoctial Orbit Elements from a Cartesian state vector.

Arguments

-u::AbstractVector{<:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_J2EqOE::SVector{6, <:Number}: The J2 Perturbed Equinoctial Orbit Element vector [n; h; k; p; q; L].

source
AstroCoords.cart2MilMethod
cart2Mil(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Converts Cartesian state vector into the Milankovich state vector.

Arguments

-u::AbstractVector{<:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż]. -μ::Number: Standard graviational parameter of central body.

Keyword Arguments

-equatorial_tol::Float64: The tolerance on what is considered an equatorial orbit (no inclination). [Default=1e-15] -circular_tol::Float64: The tolerance on what is considered a circular orbit (no eccentricity). [Default=1e-15]

Returns

-u_Mil::SVector{7, <:Number}: The Milankovich state vector [H; e; L].

source
AstroCoords.cart2cylindMethod

cart2cylind(u::AbstractVector{T}, μ::Number) where {T<:Number}

Computes the cylindrical orbital elements from a Cartesian set.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: Cartesian state vector [x; y; z; ẋ; ẏ; ż]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_cylind::SVector{6, <:Number}`: The cylindrical orbital element vector [r; θ; z; ṙ; θdot; ż].

source
AstroCoords.cart2delaunayMethod
cart2delaunay(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Computes the Delaunay orbital elements from a Cartesian set. Laskar, Jacques. "Andoyer construction for Hill and Delaunay variables." Celestial Mechanics and Dynamical Astronomy 128.4 (2017): 475-482.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The Cartesian orbital element vector [x; y; z; ẋ; ẏ; ż]. -μ::Number: Standard graviational parameter of central body.

Returns

-'u_cart::SVector{6, <:Number}': Delaunay Orbital Element Vector [L; G; H; M; ω; Ω]

source
AstroCoords.cart2koeMethod
function cart2koe(
    u::AbstractVector{T}, μ::V; equatorial_tol::Float64=1E-15, circular_tol::Float64=1E-15
) where {T<:Number,V<:Number}

Computes the Keplerian orbital elements from a cartesian set.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż]. -μ::Number: Standard graviational parameter of central body.

Keyword Arguments

-equatorial_tol::Float64: The tolerance on what is considered an equatorial orbit (no inclination). [Default=1e-15] -circular_tol::Float64: The tolerance on what is considered a circular orbit (no eccentricity). [Default=1e-15]

Returns

-u_koe::SVector{6, <:Number}`: Keplerian orbital element vector [a; e; i; Ω(RAAN); ω(AOP); f(True Anomaly)].

source
AstroCoords.cart2sphereMethod
cart2sphere(u::AbstractVector{T}, μ::Number) where {T<:Number}

Computes the spherical orbital elements from a spherical set.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The Cartesian state vector [x; y; z; ẋ; ẏ; ż]. -μ::Number: Standard graviational parameter of central body.

Returns

-'u_sphere::SVector{6, <:Number}': Spherical Orbital Element Vector [r; θ; ϕ; ṙ; θdot; ϕdot]

source
AstroCoords.cylind2cartMethod
cylind2cart(u::AbstractVector{T}, μ::Number) where {T<:Number}

Computes the Cartesian orbital elements from a cylindrical set.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The cylindrical state vector [r; θ; z; ṙ; θdot; ż]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_cart::SVector{6, <:Number}: The Cartesian orbital element vector [x; y; z; ẋ; ẏ; ż].

source
AstroCoords.delaunay2cartMethod
delaunay2cart(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Computes the Cartesian orbital elements from a Delaunay set. Laskar, Jacques. "Andoyer construction for Hill and Delaunay variables." Celestial Mechanics and Dynamical Astronomy 128.4 (2017): 475-482.

Note

All angles are in radians.

Argument

-u::AbstractVector{<:Number}: The Delaunay orbital element vector [L; G; H; M; ω; Ω]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_cart::SVector{6, <:Number}`: The cartesian orbital element vector [x; y; z; ẋ; ẏ; ż].

source
AstroCoords.eccentricAnomaly2MeanAnomalyMethod
eccentricAnomaly2MeanAnomaly(E::Number, e::Number)

Converts the true anomaly into the mean anomaly.

Arguments

-E::Number: Eccentric anomaly of the orbit [radians]. -e::Number: Eccentricity of the orbit.

Returns

-'M::Number': Mean anomaly of the orbit [radians].

source
AstroCoords.eccentricAnomaly2TrueAnomalyMethod
eccentricAnomaly2TrueAnomaly(E::Number, e::Number)

Converts the eccentric anomaly into the true anomaly.

Arguments

-E::Number: Eccentric anomaly of the orbit [radians]. -e::Number: Eccentricity of the orbit.

Returns

-f::Number: True anomaly of the orbit [radians].

source
AstroCoords.koe2IOEMethod
function koeM2IOE(u::AbstractVector{T}, μ::V) where {T<:Number, V<:Number}

Computes the Intermediate Orbit Elements from a Keplerian set.

Arguments

-u::AbstractVector{<:Number}: The Keplerian state vector [a; e; i; Ω(RAAN); ω(AOP); M(Mean Anomaly)]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_IOE::SVector{6, <:Number}: The Intermediate Orbit Element vector [I1; I2; I3; I4; I5; I6].

source
AstroCoords.koe2ModEqMethod
koe2ModEq(u::AbstractVector{T}, μ::Number) where {T<:Number}

Converts Keplerian elements into the Modified Equinoctial elements.

Note

All angles are in radians.

Arguments

-u:AbstractVector{<:Number}: The Keplerian State vector [a; e; i; Ω(RAAN); ω(AOP); ν(True Anomaly)]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_ModEq::SVector{6, <:Number}: The Modified Equinoctial state vector [p; f; g; h; k; l].

source
AstroCoords.koe2USM7Method
koe2USM7(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Converts Keplerian orbital elements into the Unified State Model set. Van den Broeck, Michael. "An Approach to Generalizing Taylor Series Integration for Low-Thrust Trajectories." (2017). https://repository.tudelft.nl/islandora/object/uuid%3A2567c152-ab56-4323-bcfa-b076343664f9

Note

All angles are in radians.

Arguments

-u:AbstractVector{<:Number}: The Keplerian state vector [a; e; i; Ω(RAAN); ω(AOP); f(True Anomaly)]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_USM::SVector{7, <:Number}: The Unified State Model vector [C; Rf1; Rf2; ϵO1; ϵO2; ϵO3; η0].

source
AstroCoords.koe2cartMethod
koe2cart(u::AbstractVector{T}, μ::V) where {T<:Number,V<:Number}

Computes the Cartesian orbital elements from a Keplerian set.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The Keplerian state vector [a; e; i; Ω(RAAN); ω(AOP); f(True Anomaly)]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_cart::Vector{6, <:Number}: The Keplerian orbital element vector [x; y; z; ẋ; ẏ; ż].

source
AstroCoords.meanAnomaly2EccentricAnomalyMethod
meanAnomaly2EccentricAnomaly(M::T, e::Number; tol::Float64=10 * eps(T)) where {T<:Number}

Converts the Mean Anomaly into the Eccentric Anomaly

Arguments

-M::Number: Mean Anomaly of the orbit [radians] -e::Number: Eccentricity of the orbit

Keyword Arguments

-tol::Float64: Convergence tolerance of Kepler solver. [Default=10*eps(T)]

Returns

-E::Number: Eccentric Anomaly of the orbit [radians]

source
AstroCoords.meanAnomaly2TrueAnomalyMethod
meanAnomaly2TrueAnomaly(M::T, e::Number; tol::Float64=10 * eps(T)) where {T<:Number}

Converts the mean anomaly into the true anomaly.

Arguments

-M::Number: Mean anomaly of the orbit [radians]. -e::Number: Eccentricity of the orbit.

Keyword Arguments

-tol::Float64: Convergence tolerance of Kepler solver. [Default=10*eps(T)]

Returns

-f::Number: Mean anomaly of the orbit [radians].

source
AstroCoords.meanMotionMethod
meanMotion(X::AstroCoord, μ::Number)

Computes the Keplerian mean motion about a central body.

Arguments

-X::AstroCoord: An coordinate set describing the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

  • n::Number: The orbital mean motion.
source
AstroCoords.meanMotionMethod
meanMotion(a::Number, μ::Number)

Computes the Keplerian mean motion about a central body.

Arguments

-a::Number: The semi-major axis of the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

  • n::Number: The orbital mean motion.
source
AstroCoords.modEqN2IOEMethod
function modEq2IOE(u::AbstractVector{T}, μ::V) where {T<:Number, V<:Number}

Computes the Intermediate Orbit Elements from a Modified Equinoctial set.

Arguments

-u::AbstractVector{<:Number}: The Modified Equinoctial state vector [p; f; g; h; k; L]. -μ::Number: Standard graviational parameter of central body.

Returns

-u_IOE::SVector{6, <:Number}: The Intermediate Orbit Element vector [I1; I2; I3; I4; I5; I6].

source
AstroCoords.orbitalNRGMethod
orbitalNRG(X::AstroCoord, μ::Number)

Computes the keplerian orbital energy.

Arguments

-X::AstroCoord: An coordinate set describing the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

-NRG::Number: The orbital energy.

source
AstroCoords.orbitalNRGMethod
orbitalNRG(a::Number, μ::Number)

Computes the keplerian orbital energy.

Arguments

-a::Number: The semi-major axis of the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

-NRG::Number: The orbital energy.

source
AstroCoords.orbitalPeriodMethod
orbitalPeriod(X::AstroCoord, μ::Number)

Computes the Keplerian orbital period about a central body.

Arguments

-X::AstroCoord: An coordinate set describing the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

-T::Number: The orbital period.

source
AstroCoords.orbitalPeriodMethod
orbitalPeriod(a::Number, μ::Number)

Computes the Keplerian orbital period about a central body.

Arguments

-a::Number: The semi-major axis of the orbit. -μ::Number: Standard graviational parameter of central body.

Returns

-T::Number: The orbital period.

source
AstroCoords.sphere2cartMethod
sphere2cart(u::AbstractVector{T}, μ::Number) where {T<:Number}

Computes the Cartesian orbital elements from a spherical set.

Note

All angles are in radians.

Arguments

-u::AbstractVector{<:Number}: The spherical orbital element vector [r; θ; ϕ; ṙ; θdot; ϕdot]. -μ::Number: Standard graviational parameter of central body.

Returns

-'u_cart::SVector{6, <:Number}': The Cartesian orbital element vector [x; y; z; ẋ; ẏ; ż].

source
AstroCoords.trueAnomaly2EccentricAnomalyMethod
trueAnomaly2EccentricAnomaly(f::Number, e::Number)

Converts the true anomaly into the mean anomaly.

Arguments

-f::Number: True anomaly of the orbit [radians]. -e::Number: Eccentricity of the orbit.

Returns

-E::Number: Eccentric anomaly of the orbit [radians].

source
AstroCoords.trueAnomaly2MeanAnomalyMethod
trueAnomaly2MeanAnomaly(f::Number, e::Number)

Converts the true anomaly into the mean anomaly.

Arguments

-f::Number: True anomaly of the orbit [radians]. -e::Number: Eccentricity of the orbit.

Returns

-M::Number: Mean anomaly of the orbit [radians].

source
Base.invMethod
inv(trans::Transformation)

Returns the inverse (or reverse) of the transformation trans

source