Core routines
Parameters & Interpolations
Main.HighFidelityEphemerisModel.HighFidelityEphemerisModelParameters
— TypeConstruct HighFidelityEphemerisModelParameters struct.
Arguments
et0::Float64
: reference epoch in seconds past J2000DU::Real
: canonical distance unitGMs::Vector{Float64}
: gravitational constants of the bodies, in km^3/s^2naif_ids::Vector{String}
: NAIF IDs of the bodiesnaif_frame::String
: inertial frame in which dynamics is integratedabcorr::String
: aberration correction for querying ephemerides of third bodiesfilepath_spherical_harmonics::Union{Nothing,String}
: path to spherical harmonics data filenmax::Int
: maximum degree of spherical harmonics to be includedframe_PCPF::Union{Nothing,String}
: NAIF frame of planet-centered planet-fixed frameget_jacobian_func::Bool
: whether to construct symbolic Jacobian function (only forNbody
dynamics)interpolate_ephem_span::Union{Nothing,Vector{Float64}}
: span of epochs to interpolate ephemeridesinterpolation_time_step::Real
: time step for interpolation
Base.show
— MethodOverload method for showing InterpolatedEphemeris
Main.HighFidelityEphemerisModel.get_pos
— MethodInterpolate ephemeris position at a given epoch
Arguments
ephem::InterpolatedEphemeris
: interpolated ephemeris structet::Float64
: epoch to interpolate
Main.HighFidelityEphemerisModel.get_state
— MethodInterpolate ephemeris state at a given epoch
Arguments
ephem::InterpolatedEphemeris
: interpolated ephemeris structet::Float64
: epoch to interpolate
Main.HighFidelityEphemerisModel.InterpolatedEphemeris
— TypeInterpolatedEphemeris struct
Fields
naif_id::String
: NAIF ID of the bodyet_range::Tuple{Float64, Float64}
: span of epochs to interpolatesplines::Array{Spline1D, 1}
: splines for the interpolated ephemerisrescale_epoch::Bool
: whether to rescale the epoch to the canonical time unittstar::Float64
: canonical time unit
Arguments
naif_id::String
: NAIF ID of the bodyets::Vector{Float64}
: epochs to interpolatervs::Array{Float64, 2}
: position and velocity vectors of the body, in km and km/srescale_epoch::Bool
: whether to rescale the epoch to the canonical time unittstar::Float64
: canonical time unitspline_order::Int
: order of the spline
Base.show
— MethodOverload method for showing InterpolatedTransformation
Main.HighFidelityEphemerisModel.get_euler_angles
— MethodInterpolate Euler angles at a given epoch
Arguments
transformation::InterpolatedTransformation
: interpolated transformation structet::Float64
: epoch to interpolate
Main.HighFidelityEphemerisModel.pxform
— MethodInterpolate transformation matrix at a given epoch
Arguments
transformation::InterpolatedTransformation
: interpolated transformation structet::Float64
: epoch to interpolate
Main.HighFidelityEphemerisModel.InterpolatedTransformation
— TypeInterpolatedTransformation struct
Fields
et_range::Tuple{Float64, Float64}
: span of epochs to interpolateframe_from::String
: frame from which the transformation is computedframe_to::String
: frame to which the transformation is computedaxis_sequence::Tuple{Int, Int, Int}
: sequence of axes for the Euler anglessplines::Array{Spline1D, 1}
: splines for the interpolated transformationrescale_epoch::Bool
: whether to rescale the epoch to the canonical time unitTU::Float64
: canonical time unit
Arguments
ets::Vector{Float64}
: epochs to interpolateframe_from::String
: frame from which the transformation is computedframe_to::String
: frame to which the transformation is computedrescale_epoch::Bool
: whether to rescale the epoch to the canonical time unitTU::Float64
: canonical time unitspline_order::Int
: order of the spline
Equations of motion
Main.HighFidelityEphemerisModel.dfdx_Nbody_Interp
— Methoddfdx_Nbody_Interp(x, u, params, t)
Evaluate Jacobian of N-body problem
Main.HighFidelityEphemerisModel.eom_Nbody_Interp!
— Methodeom_Nbody_Interp!(dx, x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_Nbody_Interp
— Methodeom_Nbody_Interp(x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_stm_Nbody_Interp!
— Methodeom_stm_Nbody_Interp!(dx_stm, x_stm, params, t)
Right-hand side of N-body equations of motion with STM compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_stm_Nbody_Interp_fd!
— Methodeom_stm_Nbody_Interp_fd!(dx_stm, x_stm, params, t)
Right-hand side of N-body equations of motion with STM compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.dfdx_Nbody_SPICE
— Methoddfdx_Nbody_SPICE(x, u, params, t)
Evaluate Jacobian of N-body problem
Main.HighFidelityEphemerisModel.eom_Nbody_SPICE!
— Methodeom_Nbody_SPICE!(dx, x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_Nbody_SPICE
— Methodeom_Nbody_SPICE(x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_stm_Nbody_SPICE!
— Methodeom_stm_Nbody_SPICE!(dx_stm, x_stm, params, t)
Right-hand side of N-body equations of motion with STMcompatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_stm_Nbody_SPICE_fd!
— Methodeom_stm_Nbody_SPICE_fd!(dx_stm, x_stm, params, t)
Right-hand side of N-body equations of motion with STM compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.dfdx_NbodySH_Interp_fd
— Methoddfdx_NbodySH_Interp_fd(x, u, params, t)
Evaluate Jacobian of N-body problem
Main.HighFidelityEphemerisModel.eom_NbodySH_Interp!
— Methodeom_NbodySH_Interp!(dx, x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_NbodySH_Interp
— Methodeom_NbodySH_Interp(x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_stm_NbodySH_Interp_fd!
— Methodeom_stm_NbodySH_Interp_fd!(dx_stm, x_stm, params, t)
Right-hand side of N-body equations of motion with STM compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_NbodySH_SPICE!
— Methodeom_NbodySH_SPICE!(dx, x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_NbodySH_SPICE
— Methodeom_NbodySH_SPICE(x, params, t)
Right-hand side of N-body equations of motion compatible with DifferentialEquations.jl
Main.HighFidelityEphemerisModel.eom_stm_NbodySH_SPICE_fd!
— Methodeom_stm_NbodySH_SPICE_fd!(dx_stm, x_stm, params, t)
Right-hand side of N-body equations of motion with STM compatible with DifferentialEquations.jl
Jacobians & Hessians
Main.HighFidelityEphemerisModel.eom_hessian_fd
— Methodeom_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 signatureeom(x, params, t)
–>dx
x
: State vectoru
: Control inputparams
: Parameterst
: Time
Returns
hess
: Hessian of the equations of motion
Main.HighFidelityEphemerisModel.eom_jacobian_fd
— Methodeom_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 signatureeom(x, params, t)
–>dx
x
: State vectoru
: Control inputparams
: Parameterst
: Time
Returns
jac
: Jacobian of the equations of motion
Perturbations
Main.HighFidelityEphemerisModel.get_VWnm
— MethodGet multipliers into gravity potential c.f. Montenbruck & Gill pg.66 eqn (3.27)
Arguments
phi::Real
: latitudelambda::Real
: longitudeR::Real
: radius of Earthr::Real
: radius of spacecraftn::Int
: degreem::Int
: order
Returns
Vnm::Real
: Vnm componentWnm::Real
: Wnm component
Main.HighFidelityEphemerisModel.legendre
— MethodCompute Legendre function of degree n order m
Arguments
n::Int
: degree (must be >= m)m::Int
: ordert::Real
: argument (cosine of latitude)
Returns
Pnm::Real
: Legendre function value
Main.HighFidelityEphemerisModel.spherical_harmonics_accel_PCPF
— MethodGet acceleration due to potential up to degree nmax in planet-centered planet-fixed frame
Arguments
phi::Real
: latitudelambda::Real
: longitudeCS::Matrix
: coefficient matrixGM::Real
: gravitational parameterR::Real
: reference radiusr::Real
: radius of spacecraftnmax::Int
: maximum degree
Returns
accel_PCPF::Vector{Real}
: acceleration in planet-centered planet-fixed frame
Main.HighFidelityEphemerisModel.spherical_harmonics_nm_accel_PCPF
— MethodGet acceleration due to (n,m) potential in planet-centered planet-fixed frame c.f. Montenbruck & Gill pg.68 eqn (3.33)
Arguments
phi::Real
: latitudelambda::Real
: longitudeCnm_dict::Dict
: dictionary of Cnm coefficients (de-normalized)Snm_dict::Dict
: dictionary of Snm coefficients (de-normalized)GM::Real
: gravitational parameterR::Real
: reference radiusr::Real
: radius of spacecraftn::Int
: degreem::Int
: order
Returns
accel_PCPF::Vector{Real}
: acceleration in planet-centered planet-fixed frame
Main.HighFidelityEphemerisModel.third_body_accel
— Methodthird_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.
Events
Main.HighFidelityEphemerisModel.angle_difference
— Methodangle_difference(ϕ_fwd::Real, ϕ_bck::Real)
Compute angle difference for periodic angles between 0 and 2π
Main.HighFidelityEphemerisModel.cart2trueanomaly
— MethodCompute 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 parameterto2pi::Bool
: if true, return the true anomaly in the range [0, 2π]
Main.HighFidelityEphemerisModel.get_trueanomaly_event
— MethodGet event function to detect target osculating true anomaly
Arguments
θ_target::Real
: target osculating true anomalyt_bounds::Tuple{Real,Real}
: time boundsradius_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}
Main.HighFidelityEphemerisModel.mod_custom
— Methodmod_custom(a, n)
Custom modulo function Ref: https://stackoverflow.com/questions/1878907/how-can-i-find-the-difference-between-two-angles