J2 Osculating Analytical Orbit Propagator
This algorithm uses the J2 propagator to obtain the secular effects of the Keplerian elements caused only by the J2 term of the geopotential field. Afterward, it adds short-term perturbations. This model is useful when fitting an orbit to a set of mean elements for the J2 orbit propagator. Hence, we can use it to verify, for example, how close a satellite is from a Sun-Synchronous orbit.
Algorithm
The algorithm implemented here is based on [1].
Initialization
We can initialize the J2 osculating analytical orbit propagator with the following function:
Propagators.init(Val(:J2osc), orb₀::KeplerianElements; kwargs...) -> OrbitPropagatorJ2Osculating
which creates a J2 osculating propagator structure OrbitPropagatorJ2Osculating
with the mean Keplerian elements orb₀
.
The following keyword selects the gravitational constants for the propagation algorithm:
j2c::J2PropagatorConstants
: J2 orbit propagator constants (seeJ2PropagatorConstants
). (Default =j2c_egm2008
)
This package contains some pre-built propagation constants for this propagator:
J2 Propagator Constant | Description | Type |
---|---|---|
j2c_egm2008 | EGM-2008 gravitational constants | Float64 |
j2c_egm2008_f32 | EGM-2008 gravitational constants | Float32 |
j2c_egm1996 | EGM-1996 gravitational constants | Float64 |
j2c_egm1996_f32 | EGM-1996 gravitational constants | Float32 |
j2c_jgm02 | JGM-02 gravitational constants | Float64 |
j2c_jgm02_f32 | JGM-02 gravitational constants | Float32 |
j2c_jgm03 | JGM-03 gravitational constants | Float64 |
j2c_jgm03_f32 | JGM-03 gravitational constants | Float32 |
The type used in the propagation will be the same as used to define the constants in the structure j2c
.
julia> orb = KeplerianElements( date_to_jd(2023, 1, 1, 0, 0, 0), 7190.982e3, 0.001111, 98.405 |> deg2rad, 100 |> deg2rad, 90 |> deg2rad, 19 |> 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 : 100.0 ° Arg. of Perigee : 90.0 ° True Anomaly : 19.0 °
julia> orbp = Propagators.init(Val(:J2osc), orb)
OrbitPropagatorJ2Osculating{Float64, Float64}: Propagator name : J2 Osculating Orbit Propagator Propagator epoch : 2023-01-01T00:00:00 Last propagation : 2023-01-01T00:00:00
Fitting Mean Elements
We can use the function:
Propagators.fit_mean_elements(::Val{:J2osc}, vjd::AbstractVector{Tjd}, vr_i::AbstractVector{Tv}, vv_i::AbstractVector{Tv}; kwargs...) -> KeplerianElements{Float64, Float64}, SMatrix{6, 6, Float64}
to 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].
It returns the fitted Keplerian elements and the final covariance matrix of the least-square algorithm.
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.
The following keywords are available to configure the fitting process:
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))
)
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 = Propagators.fit_mean_elements(Val(:J2), 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: 1 2.13455 0.00396508 2134.56 --- PROGRESS: 2 7.5425e-05 0.00975016 9.75045 -99.5432 % PROGRESS: 3 6.15798e-05 0.00975023 9.75043 -0.0002039 % (KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T16:33:40.388), [0.9999771881586103 3.3330538199638824e-7 … -9.720588376900117e-5 -7.124383306089359e-5; 3.3330159687496843e-7 0.9999779181614062 … 0.0032646452140992063 1.701808112036195e-5; … ; -9.720588377946891e-5 0.0032646452140986325 … 2.2082152622482682e-5 1.7944993096216286e-8; -7.124383307731214e-5 1.701808111721979e-5 … 1.7944993088607555e-8 2.172473646254285e-5])
julia> orb
KeplerianElements{Float64, Float64}: Epoch : 2.46003e6 (2023-03-24T16:33:40.388) Semi-major axis : 7155.99 km Eccentricity : 0.00307382 Inclination : 98.4357 ° RAAN : 162.113 ° Arg. of Perigee : 33.0518 ° True Anomaly : 344.956 °
References
- [1] Vallado, D. A (2013). Fundamentals of Astrodynamics and Applications. 4th ed. Microcosm Press, Hawthorn, CA, USA.