Core routines

Parameters & Interpolations

Main.HighFidelityEphemerisModel.HighFidelityEphemerisModelParametersType

Construct HighFidelityEphemerisModelParameters struct.

Arguments

  • et0::Float64: reference epoch in seconds past J2000
  • DU::Real: canonical distance unit
  • GMs::Vector{Float64}: gravitational constants of the bodies, in km^3/s^2
  • naif_ids::Vector{String}: NAIF IDs of the bodies
  • naif_frame::String: inertial frame in which dynamics is integrated
  • abcorr::String: aberration correction for querying ephemerides of third bodies
  • filepath_spherical_harmonics::Union{Nothing,String}: path to spherical harmonics data file
  • nmax::Int: maximum degree of spherical harmonics to be included
  • frame_PCPF::Union{Nothing,String}: NAIF frame of planet-centered planet-fixed frame
  • get_jacobian_func::Bool: whether to construct symbolic Jacobian function (only for Nbody dynamics)
  • interpolate_ephem_span::Union{Nothing,Vector{Float64}}: span of epochs to interpolate ephemerides
  • interpolation_time_step::Real: time step for interpolation
source
Base.showMethod

Overload method for showing InterpolatedEphemeris

source
Main.HighFidelityEphemerisModel.InterpolatedEphemerisType

InterpolatedEphemeris struct

Fields

  • naif_id::String: NAIF ID of the body
  • et_range::Tuple{Float64, Float64}: span of epochs to interpolate
  • splines::Array{Spline1D, 1}: splines for the interpolated ephemeris
  • rescale_epoch::Bool: whether to rescale the epoch to the canonical time unit
  • tstar::Float64: canonical time unit

Arguments

  • naif_id::String: NAIF ID of the body
  • ets::Vector{Float64}: epochs to interpolate
  • rvs::Array{Float64, 2}: position and velocity vectors of the body, in km and km/s
  • rescale_epoch::Bool: whether to rescale the epoch to the canonical time unit
  • tstar::Float64: canonical time unit
  • spline_order::Int: order of the spline
source
Base.showMethod

Overload method for showing InterpolatedTransformation

source
Main.HighFidelityEphemerisModel.pxformMethod

Interpolate transformation matrix at a given epoch

Arguments

  • transformation::InterpolatedTransformation: interpolated transformation struct
  • et::Float64: epoch to interpolate
source
Main.HighFidelityEphemerisModel.InterpolatedTransformationType

InterpolatedTransformation struct

Fields

  • et_range::Tuple{Float64, Float64}: span of epochs to interpolate
  • frame_from::String: frame from which the transformation is computed
  • frame_to::String: frame to which the transformation is computed
  • axis_sequence::Tuple{Int, Int, Int}: sequence of axes for the Euler angles
  • splines::Array{Spline1D, 1}: splines for the interpolated transformation
  • rescale_epoch::Bool: whether to rescale the epoch to the canonical time unit
  • TU::Float64: canonical time unit

Arguments

  • ets::Vector{Float64}: epochs to interpolate
  • frame_from::String: frame from which the transformation is computed
  • frame_to::String: frame to which the transformation is computed
  • rescale_epoch::Bool: whether to rescale the epoch to the canonical time unit
  • TU::Float64: canonical time unit
  • spline_order::Int: order of the spline
source

Equations of motion

Jacobians & Hessians

Main.HighFidelityEphemerisModel.eom_hessian_fdMethod
eom_hessian_fd(eom::Function, x, u, params, t)

Evaluate Hessian of equations of motion using ForwardDiff. The second argument u is a place-holder for control input.

Arguments

  • eom: "static" equations of motion, with signature eom(x, params, t) –> dx
  • x: State vector
  • u: Control input
  • params: Parameters
  • t: Time

Returns

  • hess: Hessian of the equations of motion
source
Main.HighFidelityEphemerisModel.eom_jacobian_fdMethod
eom_jacobian_fd(eom::Function, x, u, params, t)

Evaluate Jacobian of equations of motion using ForwardDiff. The second argument u is a place-holder for control input.

Arguments

  • eom: "static" equations of motion, with signature eom(x, params, t) –> dx
  • x: State vector
  • u: Control input
  • params: Parameters
  • t: Time

Returns

  • jac: Jacobian of the equations of motion
source

Perturbations

Main.HighFidelityEphemerisModel.get_VWnmMethod

Get multipliers into gravity potential c.f. Montenbruck & Gill pg.66 eqn (3.27)

Arguments

  • phi::Real: latitude
  • lambda::Real: longitude
  • R::Real: radius of Earth
  • r::Real: radius of spacecraft
  • n::Int: degree
  • m::Int: order

Returns

  • Vnm::Real: Vnm component
  • Wnm::Real: Wnm component
source
Main.HighFidelityEphemerisModel.legendreMethod

Compute Legendre function of degree n order m

Arguments

  • n::Int: degree (must be >= m)
  • m::Int: order
  • t::Real: argument (cosine of latitude)

Returns

  • Pnm::Real: Legendre function value
source
Main.HighFidelityEphemerisModel.spherical_harmonics_accel_PCPFMethod

Get acceleration due to potential up to degree nmax in planet-centered planet-fixed frame

Arguments

  • phi::Real: latitude
  • lambda::Real: longitude
  • CS::Matrix: coefficient matrix
  • GM::Real: gravitational parameter
  • R::Real: reference radius
  • r::Real: radius of spacecraft
  • nmax::Int: maximum degree

Returns

  • accel_PCPF::Vector{Real}: acceleration in planet-centered planet-fixed frame
source
Main.HighFidelityEphemerisModel.spherical_harmonics_nm_accel_PCPFMethod

Get acceleration due to (n,m) potential in planet-centered planet-fixed frame c.f. Montenbruck & Gill pg.68 eqn (3.33)

Arguments

  • phi::Real: latitude
  • lambda::Real: longitude
  • Cnm_dict::Dict: dictionary of Cnm coefficients (de-normalized)
  • Snm_dict::Dict: dictionary of Snm coefficients (de-normalized)
  • GM::Real: gravitational parameter
  • R::Real: reference radius
  • r::Real: radius of spacecraft
  • n::Int: degree
  • m::Int: order

Returns

  • accel_PCPF::Vector{Real}: acceleration in planet-centered planet-fixed frame
source
Main.HighFidelityEphemerisModel.third_body_accelMethod
third_body_accel(r_spacecraft, r_3body, mu_3body)

Compute third-body acceleration via Battin's formula (eqn (8.60) ~ (8.62) in [1])

[1] Battin, R. H. (1987). An introduction to the mathematics and methods of astrodynamics. American Institute of Aeronautics and Astronautics, Inc.

source

Events

Main.HighFidelityEphemerisModel.cart2trueanomalyMethod

Compute osculating true anomaly from state and mu

Arguments

  • state::Vector: state vector in Cartesian coordinates, in order [x, y, z, vx, vy, vz]
  • mu::Float64: gravitational parameter
  • to2pi::Bool: if true, return the true anomaly in the range [0, 2π]
source
Main.HighFidelityEphemerisModel.get_trueanomaly_eventMethod

Get event function to detect target osculating true anomaly

Arguments

  • θ_target::Real: target osculating true anomaly
  • t_bounds::Tuple{Real,Real}: time bounds
  • radius_bounds::Tuple{Real,Real}: radius bounds
  • θ_check_range::Float64: check range for true anomaly

Returns:

  • _condition: event function with signature (x, t, integrator) -> Union{Float64, NaN}
source