Library
Documentation for SatelliteToolboxPropagators.jl
.
SatelliteToolboxPropagators.J2OsculatingPropagator
— Typemutable struct J2OsculatingPropagator{Tepoch<:Number, T<:Number}
J2 osculating orbit propagator structure.
SatelliteToolboxPropagators.J2Propagator
— Typemutable struct J2Propagator{Tepoch<:Number, T<:Number}
J2 orbit propagator structure.
SatelliteToolboxPropagators.J2PropagatorConstants
— Typestruct J2PropagatorConstants{T<:Number}
Constants for the J2 orbit propagator.
Fields
R0::T
: Earth equatorial radius [m].μm::T
: √(GM / R0^3) [er/s]^(3/2).J2::T
: The second gravitational zonal harmonic of the Earth.
SatelliteToolboxPropagators.J4OsculatingPropagator
— Typemutable struct J4OsculatingPropagator{Tepoch<:Number, T<:Number}
J4 osculating orbit propagator structure.
SatelliteToolboxPropagators.J4Propagator
— TypeJ4Propagator{Tepoch, T}
J4 orbit propagator structure.
SatelliteToolboxPropagators.J4PropagatorConstants
— Typestruct J4PropagatorConstants{T<:Number}
Constants for the J4 orbit propagator.
Fields
R0::T
: Earth equatorial radius [m].μm::T
: √(GM / R0^3) [er/s]^(3/2).J2::T
: The second gravitational zonal harmonic of the Earth.J4::T
: The fourth gravitational zonal harmonic of the Earth.
SatelliteToolboxPropagators.OrbitPropagatorJ2
— TypeOrbitPropagatorJ2{Tepoch, T} <: OrbitPropagator{Tepoch, T}
J2 orbit propagator.
Fields
j2d
: Structure that stores the J2 orbit propagator data (seeJ2Propagator
).
SatelliteToolboxPropagators.OrbitPropagatorJ2Osculating
— TypeOrbitPropagatorJ2Osculating{Tepoch, T} <: OrbitPropagator{Tepoch, T}
J2 osculating orbit propagator.
Fields
j2oscd
: Structure that stores the J2 osculating orbit propagator data (seeJ2OsculatingPropagator
).
SatelliteToolboxPropagators.OrbitPropagatorJ4
— TypeOrbitPropagatorJ4{Tepoch, T} <: OrbitPropagator{Tepoch, T}
J4 orbit propagator.
Fields
j4d
: Structure that stores the J4 orbit propagator data (seeJ4Propagator
).
SatelliteToolboxPropagators.OrbitPropagatorJ4Osculating
— TypeOrbitPropagatorJ4Osculating{Tepoch, T} <: OrbitPropagator{Tepoch, T}
J4 osculating orbit propagator.
Fields
j4oscd
: Structure that stores the J4 osculating orbit propagator data (seeJ4OsculatingPropagator
).
SatelliteToolboxPropagators.OrbitPropagatorSgp4
— TypeOrbitPropagatorSgp4{Tepoch, T} <: OrbitPropagator{Tepoch, T}
SGP4 orbit propagator.
Fields
sgp4d
: Structure that stores the SGP4 orbit propagator data.
SatelliteToolboxPropagators.OrbitPropagatorTwoBody
— TypeOrbitPropagatorTwoBody{Tepoch, T} <: OrbitPropagator{Tepoch, T}
Two body orbit propagator.
Fields
tbd
: Structure that stores the two body orbit propagator data (seeTwoBodyPropagator
).
SatelliteToolboxPropagators.TwoBodyPropagator
— Typemutable struct TwoBodyPropagator{Tepoch<:Number, T<:Number}
Two body orbit propagator structure.
SatelliteToolboxPropagators.Propagators.fit_mean_elements!
— MethodPropagators.fit_mean_elements!(orbp::OrbitPropagatorJ2, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector}
Fit a set of mean Keplerian elements for the J2 orbit propagator orbp
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The orbit propagator orbp
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements!
— MethodPropagators.fit_mean_elements!(orbp::OrbitPropagatorJ2Osculating, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector}
Fit a set of mean Keplerian elements for the J2 osculating orbit propagator orbp
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The orbit propagator orbp
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements!
— MethodPropagators.fit_mean_elements!(orbp::OrbitPropagatorJ4, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector}
Fit a set of mean Keplerian elements for the J4 orbit propagator orbp
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The orbit propagator orbp
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements!
— MethodPropagators.fit_mean_elements!(orbp::OrbitPropagatorJ4Osculating, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector}
Fit a set of mean Keplerian elements for the J4 osculating orbit propagator orbp
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The orbit propagator orbp
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements!
— MethodPropagators.fit_mean_elements!(orbp::OrbitPropagatorSgp4, vjd::AbstractVector{Tjd}, vr_teme::AbstractVector{Tv}, vv_teme::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector}
Fit a Two-Line Element set (TLE
) for the SGP4 orbit propagator orbp
using the osculating elements represented by a set of position vectors vr_teme
[m] and a set of velocity vectors vv_teme
[m / s] represented in the True-Equator, Mean-Equinox reference frame (TEME) at instants in the array vjd
[Julian Day].
This algorithm was based on [1].
The SGP4 orbit propagator orbp
will be initialized with the TLE returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)estimate_bstar::Bool
: Iftrue
, the algorithm will try to estimate the B* parameter. Otherwise, it will be set to 0 or to the value in initial guess (see section Initial Guess). (Default = true)initial_guess::Union{Nothing, AbstractVector, TLE}
: Initial guess for the TLE fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_teme
andvv_teme
. For more information, see the section Initial Guess. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted TLE. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)classification::Char
: Satellite classification character for the output TLE. (Default = 'U')element_set_number::Int
: Element set number for the output TLE. (Default = 0)international_designator::String
: International designator string for the output TLE. (Default = "999999")name::String
: Satellite name for the output TLE. (Default = "UNDEFINED")revolution_number::Int
: Revolution number for the output TLE. (Default = 0)satellite_number::Int
: Satellite number for the output TLE. (Default = 9999)
Returns
TLE
: The fitted TLE.SMatrix{7, 7, T}
: Final covariance matrix of the least-square algorithm.
Initial Guess
This algorithm uses a least-square algorithm to fit a TLE based on a set of osculating state vectors. Since the system is chaotic, a good initial guess is paramount for algorithm convergence. We can provide an initial guess using the keyword initial_guess
.
If initial_guess
is a TLE
, we update the TLE epoch using the function update_sgp4_tle_epoch!
to the desired one in mean_elements_epoch
. Afterward, we use this new TLE as the initial guess.
If initial_guess
is an AbstractVector
, we use this vector as the initial mean state vector for the algorithm. It must contain 7 elements as follows:
┌ ┐
│ IDs 1 to 3: Mean position [km] │
│ IDs 4 to 6: Mean velocity [km / s] │
│ ID 7: Bstar [1 / er] │
└ ┘
If initial_guess
is nothing
, the algorithm takes the closest osculating state vector to the mean_elements_epoch
and uses it as the initial mean state vector. In this case, the epoch is set to the same epoch of the osculating data in vjd
. When the fitted TLE is obtained, the algorithm uses the function update_sgp4_tle_epoch!
to change its epoch to mean_elements_epoch
.
If initial_guess
is not nothing
, the B* initial estimate is obtained from the TLE or the state vector. Hence, if estimate_bstar
is false
, it will be kept constant with this initial value.
References
- [1] Vallado, D. A., Crawford, P (2008). SGP4 Orbit Determination. American Institute of Aeronautics ans Astronautics.
SatelliteToolboxPropagators.Propagators.fit_mean_elements
— MethodPropagators.fit_mean_elements(::Val{:J2osc}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J2 osculating orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J2 osculating propagator with the default constants j2c_egm2008
. If another set of constants are required, use the function Propagators.fit_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements
— MethodPropagators.fit_mean_elements(::Val{:J2}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J2 orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J2 propagator with the default constants j2c_egm2008
. If another set of constants are required, use the function Propagators.fit_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements
— MethodPropagators.fit_mean_elements(::Val{:J4osc}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J4 osculating orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J4 osculating propagator with the default constants j4c_egm2008
. If another set of constants are required, use the function Propagators.fit_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements
— MethodPropagators.fit_mean_elements(::Val{:J4}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J4 orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J4 propagator with the default constants j4c_egm2008
. If another set of constants are required, use the function Propagators.fit_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
SatelliteToolboxPropagators.Propagators.fit_mean_elements
— MethodPropagators.fit_mean_elements(::Val{:SGP4}, vjd::AbstractVector{Tjd}, vr_teme::AbstractVector{Tv}, vv_teme::AbstractVector{Tv}; kwargs...) -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a Two-Line Element set (TLE
) for the SGP4 orbit propagator using the osculating elements represented by a set of position vectors vr_teme
[m] and a set of velocity vectors vv_teme
[m / s] represented in the True-Equator, Mean-Equinox reference frame (TEME) at instants in the array vjd
[Julian Day].
This algorithm was based on [1].
This algorithm version will allocate a new SGP4 propagator with the default constants sgp4c_wgs84
. If another set of constants are required, use the function Propagators.fit_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)estimate_bstar::Bool
: Iftrue
, the algorithm will try to estimate the B* parameter. Otherwise, it will be set to 0 or to the value in initial guess (see section Initial Guess). (Default = true)initial_guess::Union{Nothing, AbstractVector, TLE}
: Initial guess for the TLE fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_teme
andvv_teme
. For more information, see the section Initial Guess. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted TLE. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)classification::Char
: Satellite classification character for the output TLE. (Default = 'U')element_set_number::Int
: Element set number for the output TLE. (Default = 0)international_designator::String
: International designator string for the output TLE. (Default = "999999")name::String
: Satellite name for the output TLE. (Default = "UNDEFINED")revolution_number::Int
: Revolution number for the output TLE. (Default = 0)satellite_number::Int
: Satellite number for the output TLE. (Default = 9999)
Returns
TLE
: The fitted TLE.SMatrix{7, 7, T}
: Final covariance matrix of the least-square algorithm.
Initial Guess
This algorithm uses a least-square algorithm to fit a TLE based on a set of osculating state vectors. Since the system is chaotic, a good initial guess is paramount for algorithm convergence. We can provide an initial guess using the keyword initial_guess
.
If initial_guess
is a TLE
, we update the TLE epoch using the function update_sgp4_tle_epoch!
to the desired one in mean_elements_epoch
. Afterward, we use this new TLE as the initial guess.
If initial_guess
is an AbstractVector
, we use this vector as the initial mean state vector for the algorithm. It must contain 7 elements as follows:
┌ ┐
│ IDs 1 to 3: Mean position [km] │
│ IDs 4 to 6: Mean velocity [km / s] │
│ ID 7: Bstar [1 / er] │
└ ┘
If initial_guess
is nothing
, the algorithm takes the closest osculating state vector to the mean_elements_epoch
and uses it as the initial mean state vector. In this case, the epoch is set to the same epoch of the osculating data in vjd
. When the fitted TLE is obtained, the algorithm uses the function update_sgp4_tle_epoch!
to change its epoch to mean_elements_epoch
.
If initial_guess
is not nothing
, the B* initial estimate is obtained from the TLE or the state vector. Hence, if estimate_bstar
is false
, it will be kept constant with this initial value.
References
- [1] Vallado, D. A., Crawford, P (2008). SGP4 Orbit Determination. American Institute of Aeronautics ans Astronautics.
SatelliteToolboxPropagators.Propagators.init!
— MethodPropagators.init!(orbp::OrbitPropagatorJ2, orb₀::KeplerianElements) -> Nothing
Initialize the J2 orbit propagator structure orbp
using the mean Keplerian elements orb₀
[SI units].
The propagation constants j2c::J2PropagatorConstants
in orbp.j2d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.Propagators.init!
— MethodPropagators.init!(orbp::OrbitPropagatorJ2Osculating, orb₀::KeplerianElements) -> Nothing
Initialize the J2 osculating orbit propagator structure orbp
using the mean Keplerian elements orb₀
[SI units].
The propagation constants j2c::J2PropagatorConstants
in orbp.j2d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.Propagators.init!
— MethodPropagators.init!(orbp::OrbitPropagatorJ4, orb₀::KeplerianElements) -> Nothing
Initialize the J4 orbit propagator structure orbp
using the mean Keplerian elements orb₀
.
The propagation constants j4c::J4PropagatorConstants
in j4d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.Propagators.init!
— MethodPropagators.init!(orbp::OrbitPropagatorJ4Osculating, orb₀::KeplerianElements) -> Nothing
Initialize the J4 osculating orbit propagator structure orbp
using the mean Keplerian elements orb₀
.
The propagation constants j4c::J4PropagatorConstants
in orbp.j4d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.Propagators.init!
— MethodPropagators.init!(orbp::OrbitPropagatorSgp4, epoch::Number, n₀::Number, e₀::Number, i₀::Number, Ω₀::Number, ω₀::Number, M₀::Number, bstar::Number; kwargs...) -> Nothing
Propagators.init!(orbp::OrbitPropagatorSgp4, tle::TLE; kwargs...) -> Nothing
Initialize the SGP4 orbit propagator structure orbp
using the initial orbit specified by the arguments.
The propagation constants sgp4c::Sgp4Constants
in orbp.sgp4d
will not be changed. Hence, they must be initialized.
Arguments
epoch::Number
: Epoch of the orbital elements [Julian Day].n₀::Number
: SGP type "mean" mean motion at epoch [rad/s].e₀::Number
: "Mean" eccentricity at epoch.i₀::Number
: "Mean" inclination at epoch [rad].Ω₀::Number
: "Mean" longitude of the ascending node at epoch [rad].ω₀::Number
: "Mean" argument of perigee at epoch [rad].M₀::Number
: "Mean" mean anomaly at epoch [rad].bstar::Number
: Drag parameter (B*).tle::TLE
: Two-line elements used for the initialization.
SatelliteToolboxPropagators.Propagators.init!
— MethodPropagators.init!(orbp::OrbitPropagatorTwoBody, orb₀::KeplerianElements) -> Nothing
Initialize the two-body orbit propagator structure orbp
using the mean Keplerian elements orb₀
.
The propagation constant m0::Number
in tbd
will not be changed. Hence, it must be initialized.
Arguments
orb₀::KeplerianElements
: Initial mean Keplerian elements [SI units].
SatelliteToolboxPropagators.Propagators.init
— MethodPropagators.init(Val(:J2osc), orb₀::KeplerianElements; kwargs...) -> OrbitPropagatorJ2Osculating
Create and initialize the J2 osculating orbit propagator structure using the mean Keplerian elements orb₀
[SI units].
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
Keywords
j2c::J2PropagatorConstants
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
SatelliteToolboxPropagators.Propagators.init
— MethodPropagators.init(Val(:J2), orb₀::KeplerianElements; kwargs...) -> OrbitPropagatorJ2
Create and initialize the J2 orbit propagator structure using the mean Keplerian elements orb₀
[SI units].
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
Keywords
j2c::J2PropagatorConstants
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
SatelliteToolboxPropagators.Propagators.init
— MethodPropagators.init(Val(:J4osc), orb₀::KeplerianElements; kwargs...) -> OrbitPropagatorJ4Osculating
Create and initialize the J4 osculating orbit propagator structure using the mean Keplerian elements orb₀
.
The type used in the propagation will be the same as used to define the constants in the structure j4c
.
Keywords
j4c::J4PropagatorConstants
: J4 orbit propagator constants (seeJ4PropagatorConstants
). (Default =j4c_egm2008
)
SatelliteToolboxPropagators.Propagators.init
— MethodPropagators.init(Val(:J4), orb₀::Orbit; kwargs...) -> OrbitPropagatorJ4
Create and initialize the J4 orbit propagator structure using the mean Keplerian elements orb₀
.
The type used in the propagation will be the same as used to define the constants in the structure j4c
.
Keywords
j4c::J4PropagatorConstants
: J4 orbit propagator constants (seeJ4PropagatorConstants
). (Default =j4c_egm2008
)
SatelliteToolboxPropagators.Propagators.init
— MethodPropagators.init(Val(:SGP4), epoch::Number, n₀::Number, e₀::Number, i₀::Number, Ω₀::Number, ω₀::Number, M₀::Number, bstar::Number; kwargs...) -> OrbitPropagatorSgp4
Propagators.init(Val(:SGP4), tle::TLE; kwargs...) -> OrbitPropagatorSgp4
Create and initialize the SGP4 orbit propagator structure using the initial orbit specified by the arguments.
The type used in the propagation will be the same as used to define the constants in the structure sgp4c
.
Arguments
epoch::Number
: Epoch of the orbital elements [Julian Day].n₀::Number
: SGP type "mean" mean motion at epoch [rad/s].e₀::Number
: "Mean" eccentricity at epoch.i₀::Number
: "Mean" inclination at epoch [rad].Ω₀::Number
: "Mean" longitude of the ascending node at epoch [rad].ω₀::Number
: "Mean" argument of perigee at epoch [rad].M₀::Number
: "Mean" mean anomaly at epoch [rad].bstar::Number
: Drag parameter (B*).tle::TLE
: Two-line elements used for the initialization.
Keywords
sgp4c::Sgp4Constants
: SGP4 orbit propagator constants (seeSgp4Constants
). (Default =sgp4c_wgs84
)
SatelliteToolboxPropagators.Propagators.init
— MethodPropagators.init(Val(:TwoBody), orb₀::KeplerianElements; kwargs...) -> OrbitPropagatorTwoBody
Create and initialize the two-body orbit propagator structure using the mean Keplerian elements orb₀
.
The type used in the propagation will be the same as used to define the gravitational constant m0
.
Keywords
m0::T
: Standard gravitational parameter of the central body [m³ / s²]. (Default =tbc_m0
)
SatelliteToolboxPropagators._j2_jacobian
— Method_j2_jacobian(j2d::J2OsculatingPropagator{Tepoch, T}, Δt::Number, x₁::SVector{6, T}, y₁::SVector{6, T}; kwargs...)) where {T<:Number, Tepoch<:Number} -> SMatrix{6, 6, T}
Compute the J2 orbit propagator Jacobian by finite-differences using the propagator j2d
at instant Δt
considering the input mean elements x₁
that must provide the output vector y₁
. Hence:
∂j2(x, Δt) │
J = ────────── │
∂x │ x = x₁
Keywords
perturbation::T
: Initial state perturbation to compute the finite-difference:Δx = x * perturbation
. (Default = 1e-3)perturbation_tol::T
: Tolerance to accept the perturbation. If the computed perturbation is lower thanperturbation_tol
, we increase it until it absolute value is higher thanperturbation_tol
. (Default = 1e-7)
SatelliteToolboxPropagators._j2osc_jacobian
— Method_j2osc_jacobian(j2oscd::J2OsculatingPropagator{Tepoch, T}, Δt::Number, x₁::SVector{6, T}, y₁::SVector{6, T}; kwargs...)) where {T<:Number, Tepoch<:Number} -> SMatrix{6, 6, T}
Compute the J2 osculating orbit propagator Jacobian by finite-differences using the propagator j2oscd
at instant Δt
considering the input mean elements x₁
that must provide the output vector y₁
. Hence:
∂j2osc(x, Δt) │
J = ───────────── │
∂x │ x = x₁
Keywords
perturbation::T
: Initial state perturbation to compute the finite-difference:Δx = x * perturbation
. (Default = 1e-3)perturbation_tol::T
: Tolerance to accept the perturbation. If the computed perturbation is lower thanperturbation_tol
, we increase it until it absolute value is higher thanperturbation_tol
. (Default = 1e-7)
SatelliteToolboxPropagators._j4_jacobian
— Method_j4_jacobian(j4d::J4OsculatingPropagator{Tepoch, T}, Δt::Number, x₁::SVector{6, T}, y₁::SVector{6, T}; kwargs...)) where {T<:Number, Tepoch<:Number} -> SMatrix{6, 6, T}
Compute the J4 orbit propagator Jacobian by finite-differences using the propagator j4d
at instant Δt
considering the input mean elements x₁
that must provide the output vector y₁
. Hence:
∂j4(x, Δt) │
J = ────────── │
∂x │ x = x₁
Keywords
perturbation::T
: Initial state perturbation to compute the finite-difference:Δx = x * perturbation
. (Default = 1e-3)perturbation_tol::T
: Tolerance to accept the perturbation. If the computed perturbation is lower thanperturbation_tol
, we increase it until it absolute value is higher thanperturbation_tol
. (Default = 1e-7)
SatelliteToolboxPropagators._j4osc_jacobian
— Method_j4osc_jacobian(j4oscd::J4OsculatingPropagator{Tepoch, T}, Δt::Number, x₁::SVector{6, T}, y₁::SVector{6, T}; kwargs...)) where {T<:Number, Tepoch<:Number} -> SMatrix{6, 6, T}
Compute the J4 osculating orbit propagator Jacobian by finite-differences using the propagator j4oscd
at instant Δt
considering the input mean elements x₁
that must provide the output vector y₁
. Hence:
∂j4osc(x, Δt) │
J = ───────────── │
∂x │ x = x₁
Keywords
perturbation::T
: Initial state perturbation to compute the finite-difference:Δx = x * perturbation
. (Default = 1e-3)perturbation_tol::T
: Tolerance to accept the perturbation. If the computed perturbation is lower thanperturbation_tol
, we increase it until it absolute value is higher thanperturbation_tol
. (Default = 1e-7)
SatelliteToolboxPropagators.fit_j2_mean_elements!
— Methodfit_j2_mean_elements!(j2d::J2Propagator{Tepoch, T}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {T<:Number, Tepoch<:Number, Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Tepoch, T}, SMatrix{6, 6, T}
Fit a set of mean Keplerian elements for the J2 orbit propagator j2d
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The J2 orbit propagator j2d
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Tepoch, T}
: Fitted Keplerian elements.SMatrix{6, 6, T}
: Final covariance matrix of the least-square algorithm.
Examples
# Allocate a new J2 orbit propagator using a dummy Keplerian elements.
julia> j2d = j2_init(KeplerianElements{Float64, Float64}(0, 7000e3, 0, 0, 0, 0, 0));
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-1781.214419290065, 1619.7795321872854, 6707.771633846665] .* 1000,
[ 5693.643675547716, -1192.342828671633, 4123.976025977494] .* 1000,
[ 5291.613719530499, -2354.5417593130833, -4175.561367156414] .* 1000,
[-2416.3705905186903, -268.74923235392623, -6715.411357310478] .* 1000,
[-6795.043410709359, 2184.4414321930635, -0.4327055325971031] .* 1000,
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[6.875680282038698, -1.864319399615942, 2.270603214569518] .* 1000,
[3.8964090757666496, -2.1887896252945875, -5.9960180359219075] .* 1000,
[-4.470258022565413, 0.5119576359985208, -5.9608372367141635] .* 1000,
[-6.647358060413909, 2.495415251255861, 2.292118747543002] .* 1000,
[0.3427096905434428, 1.040125572862349, 7.3936887585116855] .* 1000,
];
julia> vjd = [
2.46002818657856e6
2.460028200467449e6
2.460028214356338e6
2.4600282282452267e6
2.4600282421341157e6
2.4600282560230047e6
];
julia> orb, P = fit_j2_mean_elements!(j2d, vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J2 propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 23 4.3413 0.00540076 4341.31 -2.90721e-06 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T18:08:40.388), [0.16604846233666615 0.06643574144803302 … -3.85541368066423e-5 0.000124032254172014; 0.0664357414470214 0.26633448262787296 … -1.7943612056596758e-5 -1.9567956793856743e-5; … ; -3.855413680493565e-5 -1.7943612058983507e-5 … 4.3971984339116176e-7 -8.092704691911699e-8; 0.0001240322541726098 -1.9567956793417353e-5 … -8.092704692135158e-8 1.2451922454639337e-7])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T18:08:40.388)
Semi-major axis : 7131.63 km
Eccentricity : 0.00114299
Inclination : 98.4366 °
RAAN : 162.177 °
Arg. of Perigee : 101.286 °
True Anomaly : 258.689 °
SatelliteToolboxPropagators.fit_j2_mean_elements
— Methodfit_j2_mean_elements(vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J2 orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J2 propagator with the default constants j2c_egm2008
. If another set of constants are required, use the function fit_j2_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
Examples
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-1781.214419290065, 1619.7795321872854, 6707.771633846665] .* 1000,
[ 5693.643675547716, -1192.342828671633, 4123.976025977494] .* 1000,
[ 5291.613719530499, -2354.5417593130833, -4175.561367156414] .* 1000,
[-2416.3705905186903, -268.74923235392623, -6715.411357310478] .* 1000,
[-6795.043410709359, 2184.4414321930635, -0.4327055325971031] .* 1000,
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[6.875680282038698, -1.864319399615942, 2.270603214569518] .* 1000,
[3.8964090757666496, -2.1887896252945875, -5.9960180359219075] .* 1000,
[-4.470258022565413, 0.5119576359985208, -5.9608372367141635] .* 1000,
[-6.647358060413909, 2.495415251255861, 2.292118747543002] .* 1000,
[0.3427096905434428, 1.040125572862349, 7.3936887585116855] .* 1000,
];
julia> vjd = [
2.46002818657856e6
2.460028200467449e6
2.460028214356338e6
2.4600282282452267e6
2.4600282421341157e6
2.4600282560230047e6
];
julia> orb, P = fit_j2_mean_elements(vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J2 propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 23 4.3413 0.00540076 4341.31 -2.90721e-06 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T18:08:40.388), [0.16604846233666615 0.06643574144803302 … -3.85541368066423e-5 0.000124032254172014; 0.0664357414470214 0.26633448262787296 … -1.7943612056596758e-5 -1.9567956793856743e-5; … ; -3.855413680493565e-5 -1.7943612058983507e-5 … 4.3971984339116176e-7 -8.092704691911699e-8; 0.0001240322541726098 -1.9567956793417353e-5 … -8.092704692135158e-8 1.2451922454639337e-7])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T18:08:40.388)
Semi-major axis : 7131.63 km
Eccentricity : 0.00114299
Inclination : 98.4366 °
RAAN : 162.177 °
Arg. of Perigee : 101.286 °
True Anomaly : 258.689 °
SatelliteToolboxPropagators.fit_j2osc_mean_elements!
— Methodfit_j2osc_mean_elements!(j2oscd::J2OsculatingPropagator{Tepoch, T}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {T<:Number, Tepoch<:Number, Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Tepoch, T}, SMatrix{6, 6, T}
Fit a set of mean Keplerian elements for the J2 osculating orbit propagator j2oscd
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The J2 osculating orbit propagator j2oscd
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Tepoch, T}
: Fitted Keplerian elements.SMatrix{6, 6, T}
: Final covariance matrix of the least-square algorithm.
Examples
# Allocate a new J2 osculating orbit propagator using a dummy Keplerian elements.
julia> j2oscd = j2osc_init(KeplerianElements{Float64, Float64}(0, 7000e3, 0, 0, 0, 0, 0));
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-6357.88873265975, 2391.9476768911686, 2181.838771262736] .* 1000
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[2.5285015912807003, 0.27812476784300005, 7.030323100703928] .* 1000
];
julia> vjd = [
2.46002818657856e6,
2.460028190050782e6
];
julia> orb, P = fit_j2osc_mean_elements!(j2oscd, vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J2 osculating propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 4 1.69161e-05 0.00260193 2.60198 -2.06772e-09 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T16:33:40.388), [0.9999427882949705 -0.0049011518111948425 … -0.00012021282848110985 -0.00011550570919063845; -0.00490115181520327 1.0007325785722665 … 0.0032661568999667063 3.8567115793316056e-5; … ; -0.00012021282849269182 0.003266156899966125 … 2.204826778451109e-5 5.382733068487529e-8; -0.00011550570919086411 3.8567115786674836e-5 … 5.382733066340212e-8 2.1657420198753813e-5])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T16:33:40.388)
Semi-major axis : 7135.8 km
Eccentricity : 0.00135383
Inclination : 98.4304 °
RAAN : 162.113 °
Arg. of Perigee : 64.9256 °
True Anomaly : 313.085 °
SatelliteToolboxPropagators.fit_j2osc_mean_elements
— Methodfit_j2osc_mean_elements(vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J2 osculating orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J2 osculating propagator with the default constants j2c_egm2008
. If another set of constants are required, use the function fit_j2osc_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
Examples
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-6357.88873265975, 2391.9476768911686, 2181.838771262736] .* 1000
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[2.5285015912807003, 0.27812476784300005, 7.030323100703928] .* 1000
];
julia> vjd = [
2.46002818657856e6,
2.460028190050782e6
];
julia> orb, P = fit_j2osc_mean_elements(vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J2 osculating propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 4 1.69161e-05 0.00260193 2.60198 -2.06772e-09 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T16:33:40.388), [0.9999427882949705 -0.0049011518111948425 … -0.00012021282848110985 -0.00011550570919063845; -0.00490115181520327 1.0007325785722665 … 0.0032661568999667063 3.8567115793316056e-5; … ; -0.00012021282849269182 0.003266156899966125 … 2.204826778451109e-5 5.382733068487529e-8; -0.00011550570919086411 3.8567115786674836e-5 … 5.382733066340212e-8 2.1657420198753813e-5])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T16:33:40.388)
Semi-major axis : 7135.8 km
Eccentricity : 0.00135383
Inclination : 98.4304 °
RAAN : 162.113 °
Arg. of Perigee : 64.9256 °
True Anomaly : 313.085 °
SatelliteToolboxPropagators.fit_j4_mean_elements!
— Methodfit_j4_mean_elements!(j4d::J4Propagator{Tepoch, T}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {T<:Number, Tepoch<:Number, Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Tepoch, T}, SMatrix{6, 6, T}
Fit a set of mean Keplerian elements for the J4 orbit propagator j4d
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The J4 orbit propagator j4d
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Tepoch, T}
: Fitted Keplerian elements.SMatrix{6, 6, T}
: Final covariance matrix of the least-square algorithm.
Examples
# Allocate a new J4 orbit propagator using a dummy Keplerian elements.
julia> j4d = j4_init(KeplerianElements{Float64, Float64}(0, 7000e3, 0, 0, 0, 0, 0));
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-1781.214419290065, 1619.7795321872854, 6707.771633846665] .* 1000,
[ 5693.643675547716, -1192.342828671633, 4123.976025977494] .* 1000,
[ 5291.613719530499, -2354.5417593130833, -4175.561367156414] .* 1000,
[-2416.3705905186903, -268.74923235392623, -6715.411357310478] .* 1000,
[-6795.043410709359, 2184.4414321930635, -0.4327055325971031] .* 1000,
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[6.875680282038698, -1.864319399615942, 2.270603214569518] .* 1000,
[3.8964090757666496, -2.1887896252945875, -5.9960180359219075] .* 1000,
[-4.470258022565413, 0.5119576359985208, -5.9608372367141635] .* 1000,
[-6.647358060413909, 2.495415251255861, 2.292118747543002] .* 1000,
[0.3427096905434428, 1.040125572862349, 7.3936887585116855] .* 1000,
];
julia> vjd = [
2.46002818657856e6
2.460028200467449e6
2.460028214356338e6
2.4600282282452267e6
2.4600282421341157e6
2.4600282560230047e6
];
julia> orb, P = fit_j4_mean_elements!(j4d, vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J4 propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 23 4.33863 0.00539962 4338.63 -2.90777e-06 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T18:08:40.388), [0.16604866262321472 0.06643593039980161 … -3.8553206984036995e-5 0.00012403204426258087; 0.06643593040175688 0.26633435614867085 … -1.7942563356579497e-5 -1.95681107792859e-5; … ; -3.855320698470177e-5 -1.794256335557519e-5 … 4.3972013198895673e-7 -8.092682623169118e-8; 0.000124032044261828 -1.9568110780915995e-5 … -8.092682623078798e-8 1.2451901466558624e-7])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T18:08:40.388)
Semi-major axis : 7131.64 km
Eccentricity : 0.00114298
Inclination : 98.4366 °
RAAN : 162.177 °
Arg. of Perigee : 101.282 °
True Anomaly : 258.693 °
SatelliteToolboxPropagators.fit_j4_mean_elements
— Methodfit_j4_mean_elements(vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J4 orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J4 propagator with the default constants j4c_egm2008
. If another set of constants are required, use the function fit_j4_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
Examples
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-1781.214419290065, 1619.7795321872854, 6707.771633846665] .* 1000,
[ 5693.643675547716, -1192.342828671633, 4123.976025977494] .* 1000,
[ 5291.613719530499, -2354.5417593130833, -4175.561367156414] .* 1000,
[-2416.3705905186903, -268.74923235392623, -6715.411357310478] .* 1000,
[-6795.043410709359, 2184.4414321930635, -0.4327055325971031] .* 1000,
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[6.875680282038698, -1.864319399615942, 2.270603214569518] .* 1000,
[3.8964090757666496, -2.1887896252945875, -5.9960180359219075] .* 1000,
[-4.470258022565413, 0.5119576359985208, -5.9608372367141635] .* 1000,
[-6.647358060413909, 2.495415251255861, 2.292118747543002] .* 1000,
[0.3427096905434428, 1.040125572862349, 7.3936887585116855] .* 1000,
];
julia> vjd = [
2.46002818657856e6
2.460028200467449e6
2.460028214356338e6
2.4600282282452267e6
2.4600282421341157e6
2.4600282560230047e6
];
julia> orb, P = fit_j4_mean_elements(vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J4 propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 23 4.33863 0.00539962 4338.63 -2.90777e-06 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T18:08:40.388), [0.16604866262321472 0.06643593039980161 … -3.8553206984036995e-5 0.00012403204426258087; 0.06643593040175688 0.26633435614867085 … -1.7942563356579497e-5 -1.95681107792859e-5; … ; -3.855320698470177e-5 -1.794256335557519e-5 … 4.3972013198895673e-7 -8.092682623169118e-8; 0.000124032044261828 -1.9568110780915995e-5 … -8.092682623078798e-8 1.2451901466558624e-7])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T18:08:40.388)
Semi-major axis : 7131.64 km
Eccentricity : 0.00114298
Inclination : 98.4366 °
RAAN : 162.177 °
Arg. of Perigee : 101.282 °
True Anomaly : 258.693 °
SatelliteToolboxPropagators.fit_j4osc_mean_elements!
— Methodfit_j4osc_mean_elements!(j4oscd::J4OsculatingPropagator{Tepoch, T}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {T<:Number, Tepoch<:Number, Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Tepoch, T}, SMatrix{6, 6, T}
Fit a set of mean Keplerian elements for the J4 osculating orbit propagator j4oscd
using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
The J4 osculating orbit propagator j4oscd
will be initialized with the Keplerian elements returned by the function.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Tepoch, T}
: Fitted Keplerian elements.SMatrix{6, 6, T}
: Final covariance matrix of the least-square algorithm.
Examples
# Allocate a new J4 osculating orbit propagator using a dummy Keplerian elements.
julia> j4oscd = j4osc_init(KeplerianElements{Float64, Float64}(0, 7000e3, 0, 0, 0, 0, 0));
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-6357.88873265975, 2391.9476768911686, 2181.838771262736] .* 1000
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[2.5285015912807003, 0.27812476784300005, 7.030323100703928] .* 1000
];
julia> vjd = [
2.46002818657856e6,
2.460028190050782e6
];
julia> orb, P = fit_j4osc_mean_elements!(j4oscd, vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J4 osculating propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 4 1.68852e-05 0.00259716 2.59722 6.52143e-10 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T16:33:40.388), [0.999942775712293 -0.004901178682591695 … -0.00012021347649120565 -0.0001155014732440418; -0.004901178674245216 1.0007325984965127 … 0.0032661543862499806 3.8581997167487534e-5; … ; -0.00012021347646254671 0.003266154386250846 … 2.204822797002821e-5 5.3982040158043245e-8; -0.00011550147323899825 3.8581997169916946e-5 … 5.398204016417114e-8 2.1657607378172235e-5])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T16:33:40.388)
Semi-major axis : 7135.79 km
Eccentricity : 0.00135314
Inclination : 98.4304 °
RAAN : 162.113 °
Arg. of Perigee : 64.9687 °
True Anomaly : 313.042 °
SatelliteToolboxPropagators.fit_j4osc_mean_elements
— Methodfit_j4osc_mean_elements(vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector} -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
Fit a set of mean Keplerian elements for the J4 osculating orbit propagator using the osculating elements represented by a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] represented in an inertial reference frame at instants in the array vjd
[Julian Day].
This algorithm version will allocate a new J4 osculating propagator with the default constants j4c_egm2008
. If another set of constants are required, use the function fit_j4osc_mean_elements!
instead.
Keywords
atol::Number
: Tolerance for the residue absolute value. If the residue is lower thanatol
at any iteration, the computation loop stops. (Default = 2e-4)rtol::Number
: Tolerance for the relative difference between the residues. If the relative difference between the residues in two consecutive iterations is lower thanrtol
, the computation loop stops. (Default = 2e-4)initial_guess::Union{Nothing, KeplerianElements}
: Initial guess for the mean elements fitting process. If it isnothing
, the algorithm will obtain an initial estimate from the osculating elements invr_i
andvv_i
. (Default = nothing)jacobian_perturbation::Number
: Initial state perturbation to compute the finite-difference when calculating the Jacobian matrix. (Default = 1e-3)jacobian_perturbation_tol::Number
: Tolerance to accept the perturbation when calculating the Jacobian matrix. If the computed perturbation is lower thanjacobian_perturbation_tol
, we increase it until it absolute value is higher thanjacobian_perturbation_tol
. (Default = 1e-7)max_iterations::Int
: Maximum number of iterations allowed for the least-square fitting. (Default = 50)mean_elements_epoch::Number
: Epoch for the fitted mean elements. (Default = vjd[end])verbose::Bool
: Iftrue
, the algorithm prints debugging information tostdout
. (Default = true)weight_vector::AbstractVector
: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrixW
as a diagonal matrix with the elements inweight_vector
at its diagonal. (Default =@SVector(ones(Bool, 6))
)
Returns
KeplerianElements{Float64, Float64}
: Fitted Keplerian elements.SMatrix{6, 6, Float64}
: Final covariance matrix of the least-square algorithm.
Examples
julia> vr_i = [
[-6792.402703741442, 2192.6458461287293, 0.18851758695295118] .* 1000,
[-6357.88873265975, 2391.9476768911686, 2181.838771262736] .* 1000
];
julia> vv_i = [
[0.3445760107690598, 1.0395135806993514, 7.393686131436984] .* 1000,
[2.5285015912807003, 0.27812476784300005, 7.030323100703928] .* 1000
];
julia> vjd = [
2.46002818657856e6,
2.460028190050782e6
];
julia> orb, P = fit_j4osc_mean_elements(vjd, vr_i, vv_i)
ACTION: Fitting the mean elements for the J4 osculating propagator.
Iteration Position RMSE Velocity RMSE Total RMSE RMSE Variation
[km] [km / s] [ ]
PROGRESS: 4 1.68852e-05 0.00259716 2.59722 6.52143e-10 %
(KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T16:33:40.388), [0.999942775712293 -0.004901178682591695 … -0.00012021347649120565 -0.0001155014732440418; -0.004901178674245216 1.0007325984965127 … 0.0032661543862499806 3.8581997167487534e-5; … ; -0.00012021347646254671 0.003266154386250846 … 2.204822797002821e-5 5.3982040158043245e-8; -0.00011550147323899825 3.8581997169916946e-5 … 5.398204016417114e-8 2.1657607378172235e-5])
julia> orb
KeplerianElements{Float64, Float64}:
Epoch : 2.46003e6 (2023-03-24T16:33:40.388)
Semi-major axis : 7135.79 km
Eccentricity : 0.00135314
Inclination : 98.4304 °
RAAN : 162.113 °
Arg. of Perigee : 64.9687 °
True Anomaly : 313.042 °
SatelliteToolboxPropagators.j2!
— Methodj2!(j2d::J2Propagator{Tepoch, T}, t::Number) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
Propagate the orbit defined in j2d
(see J2Propagator
) to t
[s] after the epoch of the input mean elements in j2d
.
The internal values in j2d
will be modified.
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j2
— Methodj2(Δt::Number, orb₀::KeplerianElements; kwargs...) -> SVector{3, T}, SVector{3, T}, J2Propagator
Initialize the J2 propagator structure using the input elements orb₀
[SI units] and propagate the orbit until the time Δt [s].
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
Keywords
j2c::J2PropagatorConstants
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.J2Propagator
: Structure with the initialized parameters.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j2_init!
— Methodj2_init!(j2d::J2Propagator, orb₀::KeplerianElements) -> Nothing
Initialize the J2 orbit propagator structure j2d
using the mean Keplerian elements orb₀
[SI units].
The propagation constants j2c::J2PropagatorConstants
in j2d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.j2_init
— Methodj2_init(orb₀::KeplerianElements; kwargs...) -> J2Propagator
Create and initialize the J2 orbit propagator structure using the mean Keplerian elements orb₀
[SI units].
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
Keywords
j2c::J2PropagatorConstants
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
SatelliteToolboxPropagators.j2osc!
— Methodj2osc!(j2oscd::J2OsculatingPropagator{Tepoch, T}, t::Number) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
Propagate the orbit defined in j2oscd
(see J2OsculatingPropagator
) to t
[s] after the epoch of the input mean elements in j2d
.
The internal values in j2oscd
will be modified.
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j2osc
— Methodj2osc(Δt::Number, orb₀::KeplerianElements; kwargs...) -> SVector{3, T}, SVector{3, T}, J2OsculatingPropagator
Initialize the J2 osculating propagator structure using the input elements orb₀
[SI units] and propagate the orbit until the time Δt [s].
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
Keywords
j2c::J2PropagatorConstants{T}
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.J2OsculatingPropagator
: Structure with the initialized parameters.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j2osc_init!
— Methodj2osc_init!(j2oscd::J2OsculatingPropagator, orb₀::KeplerianElements) -> Nothing
Initialize the J2 osculating orbit propagator structure j2oscd
using the mean Keplerian elements orb₀
[SI units].
The propagation constants j2c::J2PropagatorConstants
in j2oscd.j2d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.j2osc_init
— Methodj2osc_init(orb₀::KeplerianElements; kwargs...) -> J2OsculatingPropagator
Create and initialize the J2 osculating orbit propagator structure using the mean Keplerian elements orb₀
[SI units].
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
Keywords
j2c::J2PropagatorConstants
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
SatelliteToolboxPropagators.j4!
— Methodj4!(j4d::J4Propagator{Tepoch, T}, t::Number) where {Tepoch<:Number, T<:Number} -> SVector{3, T}, SVector{3, T}
Propagate the orbit defined in j4d
(see J4Propagator
) to t
[s] after the epoch of the input mean elements in j4d
.
The internal values in j4d
will be modified.
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j4
— Methodj4(Δt::Number, orb₀::KeplerianElements; kwargs...) -> SVector{3, T}, SVector{3, T}, J4Propagator
Initialize the J4 propagator structure using the input elements orb₀
and propagate the orbit until the time Δt [s].
The type used in the propagation will be the same as used to define the constants in the structure j4c
.
Keywords
j4c::J4PropagatorConstants
: J4 orbit propagator constants (seeJ4PropagatorConstants
). (Default =j4c_egm2008
)
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.J4Propagator
: Structure with the initialized parameters.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j4_init!
— Methodj4_init!(j4d::J4Propagator, orb₀::KeplerianElements) -> Nothing
Initialize the J4 orbit propagator structure j4d
using the mean Keplerian elements orb₀
.
The propagation constants j4c::J4PropagatorConstants
in j4d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.j4_init
— Methodj4_init(orb₀::KeplerianElements; kwargs...) -> J4Propagator
Create and initialize the J4 orbit propagator structure using the mean Keplerian elements orb₀
.
The type used in the propagation will be the same as used to define the constants in the structure j4c
.
Keywords
j4c::J4PropagatorConstants
: J4 orbit propagator constants (seeJ4PropagatorConstants
). (Default =j4c_egm2008
)
SatelliteToolboxPropagators.j4osc!
— Methodj4osc!(j4oscd::J4OsculatingPropagator{Tepoch, T}, t::Number) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
Propagate the orbit defined in j4oscd
(see J4OsculatingPropagator
) to t
[s] after the epoch of the input mean elements in j4d
.
The internal values in j4oscd
will be modified.
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j4osc
— Methodj4osc(Δt::Number, orb₀::KeplerianElements; kwargs...) -> SVector{3, T}, SVector{3, T}, J4OsculatingPropagator
Initialize the J4 osculating propagator structure using the input elements orb₀
and propagate the orbit until the time Δt [s].
The type used in the propagation will be the same as used to define the constants in the structure j4c
.
Keywords
j4c::J4PropagatorConstants{T}
: J4 orbit propagator constants (seeJ4PropagatorConstants
). (Default =j4c_egm2008
)
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.J4OsculatingPropagator
: Structure with the initialized parameters.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters. Notice that the perturbation theory requires an inertial frame with true equator.
SatelliteToolboxPropagators.j4osc_init!
— Methodj4osc_init!(j4oscd::J4OsculatingPropagator, orb₀::KeplerianElements) -> Nothing
Initialize the J4 osculating orbit propagator structure j4oscd
using the mean Keplerian elements orb₀
.
The propagation constants j4c::J4PropagatorConstants
in j4oscd.j4d
will not be changed. Hence, they must be initialized.
SatelliteToolboxPropagators.j4osc_init
— Methodj4osc_init(orb₀::KeplerianElements; kwargs...) -> J4OsculatingPropagator
Create and initialize the J4 osculating orbit propagator structure using the mean Keplerian elements orb₀
.
The type used in the propagation will be the same as used to define the constants in the structure j4c
.
Keywords
j4c::J4PropagatorConstants
: J4 orbit propagator constants (seeJ4PropagatorConstants
). (Default =j4c_egm2008
)
SatelliteToolboxPropagators.twobody!
— Methodtwobody!(tbd::TwoBodyPropagator{Tepoch, T}, t::Number) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
Propagate the orbit defined in tbd
(see TwoBodyPropagator
) to t
[s] after the epoch of the input mean elements in tbd
.
The internal values in tbd
will be modified.
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters.
SatelliteToolboxPropagators.twobody
— Methodtwobody(Δt::Number, orb₀::KeplerianElements; kwargs...) -> SVector{3, T}, SVector{3, T}, TwoBodyPropagator
Initialize the two-body propagator structure using the input elements orb₀
and propagate the orbit until the time Δt [s].
The type used in the propagation will be the same as used to define the gravitational constant m0
.
Keywords
m0::T
: Standard gravitational parameter of the central body [m³ / s²]. (Default =tbc_m0
)
Returns
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.TwoBodyPropagator
: Structure with the initialized parameters.
Remarks
The inertial frame in which the output is represented depends on which frame it was used to generate the orbit parameters.
SatelliteToolboxPropagators.twobody_init!
— Methodtwobody_init!(tbd::TwoBodyPropagator, orb₀::KeplerianElements; kwargs...) -> Nothing
Initialize the two-body propagator structure tbd
using the mean Keplerian elements orb₀
.
The propagation constant μ::Number
in tbd
will not be changed. Hence, it must be initialized.
SatelliteToolboxPropagators.twobody_init
— Methodtwobody_init(orb₀::KeplerianElements; kwargs...) -> TwoBodyPropagator
Create and initialize the two-body propagator structure using the mean Keplerian elements orb₀
.
The type used in the propagation will be the same as used to define the gravitational constant m0
.
Keywords
m0::T
: Standard gravitational parameter of the central body [m³ / s²]. (Default =tbc_m0
)
SatelliteToolboxPropagators.update_j2_mean_elements_epoch!
— Methodupdate_j2_mean_elements_epoch!(j2d::J2Propagator, orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using the propagator j2d
to new_epoch
, which can be represented by a Julian Day or a DateTime
.
The J2 orbit propagator j2d
will be initialized with the Keplerian elements returned by the function.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
# Allocate a new J2 orbit propagator using the created Keplerian elements. Notice that any
# set of Keplerian elements can be used here.
julia> j2d = j2_init(orb);
julia> update_j2_mean_elements_epoch!(j2d, orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9565 °
Arg. of Perigee : 197.078 °
True Anomaly : 127.291 °
SatelliteToolboxPropagators.update_j2_mean_elements_epoch
— Methodupdate_j2_mean_elements_epoch(orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using a J2 orbit propagator to new_epoch
, which can be represented by a Julian Day or a DateTime
.
This algorithm version will allocate a new J2 propagator with the default constants j2c_egm2008
. If another set of constants are required, use the function update_j2osc_mean_elements_epoch!
instead.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
julia> update_j2_mean_elements_epoch(orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9565 °
Arg. of Perigee : 197.078 °
True Anomaly : 127.291 °
SatelliteToolboxPropagators.update_j2osc_mean_elements_epoch!
— Methodupdate_j2osc_mean_elements_epoch!(j2oscd::J2OsculatingPropagator, orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using the propagator j2oscd
to new_epoch
, which can be represented by a Julian Day or a DateTime
.
The J2 osculating orbit propagator j2oscd
will be initialized with the Keplerian elements returned by the function.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
# Allocate a new J2 osculating orbit propagator using the created Keplerian elements. Notice
# that any set of Keplerian elements can be used here.
julia> j2oscd = j2osc_init(orb);
julia> update_j2osc_mean_elements_epoch!(j2oscd, orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9565 °
Arg. of Perigee : 197.078 °
True Anomaly : 127.291 °
SatelliteToolboxPropagators.update_j2osc_mean_elements_epoch
— Methodupdate_j2osc_mean_elements_epoch(orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using a J2 osculating orbit propagator to new_epoch
, which can be represented by a Julian Day or a DateTime
.
This algorithm version will allocate a new J2 osculating propagator with the default constants j2c_egm2008
. If another set of constants are required, use the function update_j2osc_mean_elements_epoch!
instead.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
julia> update_j2osc_mean_elements_epoch(orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9565 °
Arg. of Perigee : 197.078 °
True Anomaly : 127.291 °
SatelliteToolboxPropagators.update_j4_mean_elements_epoch!
— Methodupdate_j4_mean_elements_epoch!(j4d::J4Propagator, orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using the propagator j4d
to new_epoch
, which can be represented by a Julian Day or a DateTime
.
The J4 orbit propagator j4d
will be initialized with the Keplerian elements returned by the function.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
# Allocate a new J4 orbit propagator using the created Keplerian elements. Notice that any
# set of Keplerian elements can be used here.
julia> j4d = j4_init(orb);
julia> update_j4_mean_elements_epoch!(j4d, orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9555 °
Arg. of Perigee : 197.079 °
True Anomaly : 127.293 °
SatelliteToolboxPropagators.update_j4_mean_elements_epoch
— Methodupdate_j4_mean_elements_epoch(orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using a J4 orbit propagator to new_epoch
, which can be represented by a Julian Day or a DateTime
.
This algorithm version will allocate a new J4 propagator with the default constants j4c_egm2008
. If another set of constants are required, use the function update_j4osc_mean_elements_epoch!
instead.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
julia> update_j4_mean_elements_epoch(orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9555 °
Arg. of Perigee : 197.079 °
True Anomaly : 127.293 °
SatelliteToolboxPropagators.update_j4osc_mean_elements_epoch!
— Methodupdate_j4osc_mean_elements_epoch!(j4oscd::J4OsculatingPropagator, orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using the propagator j4oscd
to new_epoch
, which can be represented by a Julian Day or a DateTime
.
The J4 osculating orbit propagator j4oscd
will be initialized with the Keplerian elements returned by the function.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
# Allocate a new J4 osculating orbit propagator using the created Keplerian elements. Notice
# that any set of Keplerian elements can be used here.
julia> j4oscd = j4osc_init(orb);
julia> update_j4osc_mean_elements_epoch!(j4oscd, orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9555 °
Arg. of Perigee : 197.079 °
True Anomaly : 127.293 °
SatelliteToolboxPropagators.update_j4osc_mean_elements_epoch
— Methodupdate_j4osc_mean_elements_epoch(orb::KeplerianElements, new_epoch::Union{Number, DateTime}) -> KepleriranElements
Update the epoch of the mean elements orb
using a J4 osculating orbit propagator to new_epoch
, which can be represented by a Julian Day or a DateTime
.
This algorithm version will allocate a new J4 osculating propagator with the default constants j4c_egm2008
. If another set of constants are required, use the function update_j4osc_mean_elements_epoch!
instead.
Examples
julia> orb = KeplerianElements(
DateTime("2023-01-01") |> datetime2julian,
7190.982e3,
0.001111,
98.405 |> deg2rad,
90 |> deg2rad,
200 |> deg2rad,
45 |> deg2rad
)
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-01T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.0 °
Arg. of Perigee : 200.0 °
True Anomaly : 45.0 °
julia> update_j4osc_mean_elements_epoch(orb, DateTime("2023-01-02"))
KeplerianElements{Float64, Float64}:
Epoch : 2.45995e6 (2023-01-02T00:00:00)
Semi-major axis : 7190.98 km
Eccentricity : 0.001111
Inclination : 98.405 °
RAAN : 90.9555 °
Arg. of Perigee : 197.079 °
True Anomaly : 127.293 °
SatelliteToolboxPropagators.Propagators.OrbitPropagator
— Typeabstract type OrbitPropagator{Tepoch<:Number, T<:Number}
Abstract type for the orbit propagators.
SatelliteToolboxPropagators.Propagators.epoch
— Functionepoch(orbp::OrbitPropagator{Tepoch, T}) where {Tepoch<:Number, T<:Number} -> Tepoch
Return the initial elements' epoch of the propagator orbp
[JD].
SatelliteToolboxPropagators.Propagators.fit_mean_elements
— Functionfit_mean_elements(::Val{:propagator}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector} -> <Mean elements>
fit_mean_elements(::Val{:propagator}, vsv::OrbitStateVector{Tepoch, T}; kwargs...) where {Tepoch<:Number, T<:Number} -> <Mean elements>
Fit a set of mean elements for the propagator
using the osculating state vector represented in an intertial reference frame. The state vector can be represented using a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] obtained at the instants in the array vjd
[Julian Day], or an array of OrbitStateVector
vsv
[SI], containing the same information. The keywords kwargs
depends on the propagator type.
This function returns the set of mean elements used to initialize the propagator
.
SatelliteToolboxPropagators.Propagators.fit_mean_elements!
— Functionfit_mean_elements!(orbp::OrbitPropagator, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) where {Tjd<:Number, Tv<:AbstractVector} -> <Mean elements>
fit_mean_elements!(orbp::OrbitPropagator, vsv::Vector{OrbitStateVector{Tepoch, T}}; kwargs...) where {Tepoch<:Number, T<:Number} -> <Mean elements>
Fit a set of mean elements for the propagator orbp
using the osculating state vector represented in an intertial reference frame. The state vector can be represented using a set of position vectors vr_i
[m] and a set of velocity vectors vv_i
[m / s] obtained at the instants in the array vjd
[Julian Day], or an array of OrbitStateVector
vsv
[SI], containing the same information. The keywords kwargs
depends on the propagator type.
This function returns the set of mean elements used to initialize the propagator
and also initializes orbp
with the fitted mean elements.
SatelliteToolboxPropagators.Propagators.init
— Functioninit(::Val{:propagator}, args...; kwargs...) -> OrbitPropagator
Create and initialize the orbit propagator
. The arguments args
and keywords kwargs
depends of the propagator type.
SatelliteToolboxPropagators.Propagators.init!
— Functioninit!(orbp::OrbitPropagator, args...; kwargs...) -> Nothing
Initialize the orbit propagator orbp
. The arguments args
and keywords kwargs
depends of the propagator type.
SatelliteToolboxPropagators.Propagators.last_instant
— Functionlast_instant(orbp::OrbitPropagator{Tepoch, T}) where {Tepoch<:Number, T<:Number} -> T
Return the last propagation instant [s] measured from the epoch.
SatelliteToolboxPropagators.Propagators.mean_elements
— Methodmean_elements(orbp::OrbitPropagator) -> Union{Nothing, KeplerianElements}
Return the mean elements using the structure KeplerianElements
of the latest propagation performed by orbp
. This is an optinal function in the API. It will return nothing
if the propagator does not support it.
SatelliteToolboxPropagators.Propagators.name
— Methodname(orbp::OrbitPropagator) -> String
Return the name of the orbit propagator orbp
. If this function is not defined, the structure name is used: typeof(orbp) |> string
.
SatelliteToolboxPropagators.Propagators.propagate!
— Functionpropagate!(orbp::OrbitPropagator{Tepoch, T}, Δt::Number[, sink = Tuple]) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
propagate!(orbp::OrbitPropagator{Tepoch, T}, p::Union{Dates.Period, Dates.CompoundPeriod}[, sink = Tuple]) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
propagate!(orbp::OrbitPropagator{Tepoch, T}, Δt::Number, sink = OrbitStateVector) where {Tepoch, T} -> OrbitStateVector{Tepoch, T}
propagate!(orbp::OrbitPropagator{Tepoch, T}, p::Union{Dates.Period, Dates.CompoundPeriod}, sink = OrbitStateVector) where {Tepoch, T} -> OrbitStateVector{Tepoch, T}
Propagate the orbit using orbp
by Δt
[s] or by the period defined by p
from the initial orbit epoch. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
and the output is a tuple with the position and velocity vectors.
Returns
If sink
is Tuple
:
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
If sink
is OrbitStateVector
:
OrbitStateVector{Tepoch, T}
: Structure with the orbit state vector [SI] at the propagation instant.
SatelliteToolboxPropagators.Propagators.propagate!
— Methodpropagate!(orbp::OrbitPropagator{Tepoch, T}, vt::AbstractVector[, sink = Tuple]; kwargs...) where {Tepoch <: Number, T <: Number} -> Vector{SVector{3, T}}, Vector{SVector{3, T}}
propagate!(orbp::OrbitPropagator{Tepoch, T}, vp::AbstractVector{Union{Dates.Period, Dates.CompundPeriod}}[, sink = Tuple]; kwargs...) where {Tepoch <: Number, T <: Number} -> Vector{SVector{3, T}}, Vector{SVector{3, T}}
propagate!(orbp::OrbitPropagator{Tepoch, T}, vt::AbstractVector, sink = OrbitStateVector; kwargs...) where {Tepoch <: Number, T <: Number} -> Vector{OrbitStateVector{Tepoch, T}}
propagate!(orbp::OrbitPropagator{Tepoch, T}, vp::AbstractVector{Union{Dates.Period, Dates.CompundPeriod}}, sink = OrbitStateVector; kwargs...) where {Tepoch <: Number, T <: Number} -> Vector{OrbitStateVector{Tepoch, T}}
Propagate the orbit using orbp
for every instant defined in vt
[s] or for every period defined in vp
from the initial orbit epoch. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
, and the output is a tuple with the arrays containing the position and velocity vectors.
Keywords
ntasks::Integer
: Number of parallel tasks to propagate the orbit. If it is set to a number equal or lower than 1, the function will propagate the orbit sequentially. (Default =Threads.nthreads()
)
Returns
If sink
is Tuple
:
Vector{SVector{3, T}}
: Array with the position vectors [m] in the inertial frame at each propagation instant defined invt
orvp
.Vector{SVector{3, T}}
: Array with the velocity vectors [m / s] in the inertial frame at each propagation instant defined invt
orvp
.
If sink
is OrbitStateVector
:
Vector{OrbitStateVector{Tepoch, T}}
: Array with the orbit state vectors [SI] at each propagation instant defined invt
orvp
.
SatelliteToolboxPropagators.Propagators.propagate
— Methodpropagate([sink = Tuple, ]::Val{:propagator}, Δt::Number, args...; kwargs...) -> SVector{3, T}, SVector{3, T}, OrbitPropagator{Tepoch, T}
propagate([sink = Tuple, ]::Val{:propagator}, p::Union{Dates.Period, Dates.CompundPeriod}, args...; kwargs...) -> SVector{3, T}, SVector{3, T}, OrbitPropagator{Tepoch, T}
propagate(sink = OrbitStateVector, ::Val{:propagator}, Δt::Number, args...; kwargs...) -> OrbitStateVector{Tepoch, T}, OrbitPropagator{Tepoch, T}
propagate(sink = OrbitStateVector, ::Val{:propagator}, p::Union{Dates.Period, Dates.CompundPeriod}, args...; kwargs...) -> OrbitStateVector{Tepoch, T}, OrbitPropagator{Tepoch, T}
Initialize the orbit propagator
and propagate the orbit by Δt
[s] or by the period defined by p
from the initial orbit epoch. The initialization arguments args...
and kwargs...
are the same as in the initialization function Propagators.init
. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
and the output is a tuple with the position and velocity vectors.
T
is the propagator number type. For more information, see Propagators.init
.
Returns
If sink
is Tuple
:
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.OrbitPropagator{Tepoch, T}
: Structure with the initialized propagator.
If sink
is OrbitStateVector
:
OrbitStateVector{Tepoch, T}
: Structure with the orbit state vector [SI] at the propagation instant.OrbitPropagator{Tepoch, T}
: Structure with the initialized propagator.
SatelliteToolboxPropagators.Propagators.propagate
— Methodpropagate([sink = Tuple, ]::Val{:propagator}, vt::AbstractVector, args...; kwargs...) -> Vector{SVector{3, T}}, Vector{SVector{3, T}}, OrbitPropagator{Tepoch, T}
propagate([sink = Tuple, ]::Val{:propagator}, vp::AbstractVector{Union{Dates.Period, Dates.CompundPeriod}}, args...; kwargs...) -> Vector{SVector{3, T}}, Vector{SVector{3, T}}, OrbitPropagator{Tepoch, T}
propagate(sink = OrbitStateVector, ::Val{:propagator}, vt::AbstractVector, args...; kwargs...) -> Vector{OrbitStateVector{Tepoch, T}}, OrbitPropagator{Tepoch, T}
propagate(sink = OrbitStateVector, ::Val{:propagator}, vp::AbstractVector{Union{Dates.Period, Dates.CompundPeriod}}, args...; kwargs...) -> Vector{OrbitStateVector{Tepoch, T}}, OrbitPropagator{Tepoch, T}
Initialize the orbit propagator
and propagate the orbit for every instant defined in vt
[s] or for every period defined in vp
from the initial orbit epoch. The initialization arguments args...
and kwargs...
(except for ntasks
) are the same as in the initialization function Propagators.init
. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
, and the output is a tuple with the arrays containing the position and velocity vectors.
T
is the propagator number type. For more information, see Propagators.init
.
Keywords
ntasks::Integer
: Number of parallel tasks to propagate the orbit. If it is set to a number equal or lower than 1, the function will propagate the orbit sequentially. (Default =Threads.nthreads()
)
Returns
If sink
is Tuple
:
Vector{SVector{3, T}}
: Array with the position vectors [m] in the inertial frame at each propagation instant defined invt
orvp
.Vector{SVector{3, T}}
: Array with the velocity vectors [m / s] in the inertial frame at each propagation instant defined invt
orvp
.OrbitPropagator{Tepoch, T}
: Structure with the initialized propagator.
If sink
is OrbitStateVector
:
Vector{OrbitStateVector{Tepoch, T}}
: Array with the orbit state vectors [SI] at each propagation instant defined invt
orvp
.OrbitPropagator{Tepoch, T}
: Structure with the initialized propagator.
SatelliteToolboxPropagators.Propagators.propagate_to_epoch!
— Methodpropagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, jd::Number[, sink = Tuple]) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
propagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, dt::DateTime[, sink = Tuple]) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
propagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, jd::Number, sink = OrbitStateVector) where {Tepoch, T} -> OrbitStateVector{Tepoch, T}
propagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, dt::DateTime, sink = OrbitStateVector) where {Tepoch, T} -> OrbitStateVector{Tepoch, T}
Propagate the orbit using orbp
until the epoch defined either by the Julian Day jd
[UTC] or by the DateTime
object dt
[UTC]. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
and the output is a tuple with the position and velocity vectors.
Returns
If sink
is Tuple
:
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
If sink
is OrbitStateVector
:
OrbitStateVector{Tepoch, T}
: Structure with the orbit state vector [SI] at the propagation instant.
SatelliteToolboxPropagators.Propagators.propagate_to_epoch!
— Methodpropagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, vjd::AbstractVector[, sink = Tuple]; kwargs...) where {Tepoch, T} -> Vector{SVector{3, T}}, Vector{SVector{3, T}}
propagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, vdt::AbstractVector{DateTime}[, sink = Tuple]; kwargs...) where {Tepoch, T} -> Vector{SVector{3, T}}, Vector{SVector{3, T}}
propagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, vjd::AbstractVector, sink = OrbitStateVector; kwargs...) where {Tepoch, T} -> Vector{OrbitStateVector{Tepoch, T}}
propagate_to_epoch!(orbp::OrbitPropagator{Tepoch, T}, vdt::AbstractVector{DateTime}, sink = OrbitStateVector; kwargs...) where {Tepoch, T} -> Vector{OrbitStateVector{Tepoch, T}}
Propagate the orbit using orbp
for every epoch defined in the vector of Julian Days vjd
[UTC] or in the vector of DateTime
objects vdt
[UTC]. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
, and the output is a tuple with the arrays containing the position and velocity vectors.
Keywords
ntasks::Integer
: Number of parallel tasks to propagate the orbit. If it is set to a number equal or lower than 1, the function will propagate the orbit sequentially. (Default =Threads.nthreads()
)
Returns
If sink
is Tuple
:
Vector{SVector{3, T}}
: Array with the position vectors [m] in the inertial frame at each propagation instant defined invjd
orvdt
.Vector{SVector{3, T}}
: Array with the velocity vectors [m / s] in the inertial frame at each propagation instant defined invjd
orvdt
.
If sink
is OrbitStateVector
:
Vector{OrbitStateVector{Tepoch, T}}
: Array with the orbit state vectors [SI] at each propagation instant defined invjd
orvdt
.
SatelliteToolboxPropagators.Propagators.propagate_to_epoch
— Methodpropagate_to_epoch([sink = Tuple, ]::Val{:propagator}, jd::Number, args...; kwargs...) -> SVector{3, T}, SVector{3, T}, OrbitPropagator{Tepoch, T}
propagate_to_epoch([sink = Tuple, ]::Val{:propagator}, dt::DateTime, args...; kwargs...) -> SVector{3, T}, SVector{3, T}, OrbitPropagator{Tepoch, T}
propagate_to_epoch(sink = OrbitStateVector, ::Val{:propagator}, jd::Number, args...; kwargs...) -> OrbitStateVector{Tepoch, T}, OrbitPropagator{Tepoch, T}
propagate_to_epoch(sink = OrbitStateVector, ::Val{:propagator}, dt::DateTime, args...; kwargs...) -> OrbitStateVector{Tepoch, T}, OrbitPropagator{Tepoch, T}
Initialize the orbit propagator
and propagate the orbit until the epoch defined by either the Julian Day jd
[UTC] or by a DateTime
object dt
[UTC] from the initial orbit epoch. The initialization arguments args...
and kwargs...
are the same as in the initialization function Propagators.init
. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
and the output is a tuple with the position and velocity vectors.
T
is the propagator number type. For more information, see Propagators.init
.
Returns
If sink
is Tuple
:
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.OrbitPropagator{Tepoch, T}
: Structure with the initialized propagator.
If sink
is OrbitStateVector
:
OrbitStateVector{Tepoch, T}
: Structure with the orbit state vector [SI] at the propagation instant.OrbitPropagator{Tepoch, T}
: Structure with the initialized propagator.
SatelliteToolboxPropagators.Propagators.propagate_to_epoch
— Methodpropagate_to_epoch([sink = Tuple, ]::Val{:propagator}, vjd::AbstractVector, args...; kwargs...) -> Vector{SVector{3, T}}, Vector{SVector{3, T}}, OrbitPropagator{Tepoch, T}
propagate_to_epoch([sink = Tuple, ]::Val{:propagator}, vdt::DateTime, args...; kwargs...) -> Vector{SVector{3, T}}, Vector{SVector{3, T}}, OrbitPropagator{Tepoch, T}
propagate_to_epoch(sink = OrbitStateVector, ::Val{:propagator}, vjd::AbstractVector, args...; kwargs...) -> OrbitStateVector{Tepoch, T}, OrbitPropagator{Tepoch, T}
propagate_to_epoch(sink = OrbitStateVector, ::Val{:propagator}, vdt::DateTime, args...; kwargs...) -> OrbitStateVector{Tepoch, T}, OrbitPropagator{Tepoch, T}
Initialize the orbit propagator
and propagate the orbit for every epoch defined in the vector of Julian Days vjd
[UTC] or in the vector of DateTime
objects vdt
[UTC]. The initialization arguments args...
and kwargs...
(except for ntasks
) are the same as in the initialization function Propagators.init
. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
, and the output is a tuple with the arrays containing the position and velocity vectors.
T
is the propagator number type. For more information, see Propagators.init
.
Keywords
ntasks::Integer
: Number of parallel tasks to propagate the orbit. If it is set to a number equal or lower than 1, the function will propagate the orbit sequentially. (Default =Threads.nthreads()
)
Returns
If sink
is Tuple
:
Vector{SVector{3, T}}
: Array with the position vectors [m] in the inertial frame at each propagation instant defined invjd
orvdt
.Vector{SVector{3, T}}
: Array with the velocity vectors [m / s] in the inertial frame at each propagation instant defined invjd
orvdt
.OrbitPropagator{Tepoch, T}
: Structure with the initialized parameters.
If sink
is OrbitStateVector
:
Vector{OrbitStateVector{Tepoch, T}}
: Array with the orbit state vectors [SI] at each propagation instant defined invjd
orvdt
.OrbitPropagator{Tepoch, T}
: Structure with the initialized parameters.
SatelliteToolboxPropagators.Propagators.step!
— Methodstep!(orbp::OrbitPropagator{Tepoch, T}, Δt::Number[, sink = Tuple]) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
step!(orbp::OrbitPropagator{Tepoch, T}, p::Union{Dates.Period, Dates.CompoundPeriod}[, sink = Tuple]) where {Tepoch, T} -> SVector{3, T}, SVector{3, T}
step!(orbp::OrbitPropagator{Tepoch, T}, Δt::Number, sink = OrbitStateVector) where {Tepoch, T} -> OrbitStateVector{Tepoch, T}
step!(orbp::OrbitPropagator{Tepoch, T}, p::Union{Dates.Period, Dates.CompoundPeriod}, sink = OrbitStateVector) where {Tepoch, T} -> OrbitStateVector{Tepoch, T}
Propagate the orbit using orbp
by Δt
[s] or by the period defined by p
from the current orbit epoch. The output type depends on the parameter sink
. If it is omitted, it defaults to Tuple
and the output is a tuple with the position and velocity vectors.
Returns
If sink
is Tuple
:
SVector{3, T}
: Position vector [m] represented in the inertial frame at propagation instant.SVector{3, T}
: Velocity vector [m / s] represented in the inertial frame at propagation instant.
If sink
is OrbitStateVector
:
OrbitStateVector{Tepoch, T}
: Structure with the orbit state vector [SI] at the propagation instant.