J4 Analytical Orbit Propagator

The J4 analytical orbit propagator considers the perturbation caused by the terms $J_2$, $J_2^2$, and $J_4$. It is a slightly more precise version of the J2 analytical orbit propagator useful for mission design and analysis.

Algorithm

The algorithm implemented here is based on [1, 2]. This version takes into account only the secular effects in the orbit elements.

Considering only the terms $J_2$, $J_2^2$, and $J_4$, we obtain that only the right accession of the ascending node (RAAN), the argument of perigee, and mean motion are perturbed by the geopotential field with secular effects. Hence, the algorithm propagates the mean elements using the following equations:

\[\begin{aligned} & p_0 = a_0 \cdot (1 - e_0^2) \\ & \beta = \sqrt{1 - e_0^2} \\ & \bar{n} = n_0 \cdot \left[ \right. \\ & \begin{aligned} \quad 1 + \\ \quad \frac{ 3}{ 4} & \cdot J_2 \cdot \left(\frac{R_0}{p_0}\right)^2 \cdot \beta \cdot (2 - 3 \sin^2 i_0) + \\ \quad \frac{ 3}{128} & \cdot J_2^2 \cdot \left(\frac{R_0}{p_0}\right)^4 \cdot \beta \cdot \left(120 + 64\beta - 40\beta^2 + (-240 - 192\beta + 40\beta^2) \sin{i_0}^2 + (105 + 144\beta + 25\beta^2) \sin{i_0}^4\right) -\\ \quad \frac{45}{128} & \cdot J_4 \cdot \left(\frac{R_0}{p_0}\right)^4 \cdot \beta \cdot e_0^2 \cdot \left(-8 + 40\sin{i_0}^2 - 35\sin{i_0}^4\right) \end{aligned} \\ & \left. \right] \\ & \Omega(t) = \Omega_0 + \left[ \right. \\ & \begin{aligned} \quad -\frac{ 3}{ 2} & \cdot \bar{n} \cdot J_2 \left(\frac{R_0}{p_0}\right)^2 \cdot \cos i_0 + \\ \quad \frac{ 3}{32} & \cdot \bar{n} \cdot J_2^2 \left(\frac{R_0}{p_0}\right)^4 \cdot \cos i_0 \cdot \left(-36 - 4e_0^2 + 48\beta + (40 - 5e_0^2 - 72\beta) \sin{i_0}^2\right) + \\ \quad \frac{15}{32} & \cdot n_0 \cdot J_4 \left(\frac{R_0}{p_0}\right)^4 \cdot \cos i_0 \cdot \left(8 + 12e_0^2 - (14 + 21e_0^2) \sin{i_0}^2 \right) \end{aligned} \\ & \left. \right] \cdot \left(t - t_0\right) \\ & \omega(t) = \omega_0 + \left[ \right. \\ & \begin{aligned} \quad \frac{ 3}{ 4} & \cdot \bar{n} \cdot J_2 \left(\frac{R_0}{p_0}\right)^2 \cdot (4 - 5\sin^2 i_0) + \\ \quad \frac{ 3}{128} & \cdot \bar{n} \cdot J_2^2 \left(\frac{R_0}{p_0}\right)^4 \cdot \left(384 + 96e_0^2 - 384\beta + (-824 - 116e_0^2 + 1056\beta) \sin{i_0}^2 + (430 - 5e_0^2 - 720\beta) \sin{i_0}^4\right) - \\ \quad \frac{15}{ 16} & \cdot n_0 \cdot J_2^2 \left(\frac{R_0}{p_0}\right)^4 \cdot e_0^2 \cdot \cos{i_0}^4 -\\ \quad \frac{15}{128} & \cdot n_0 \cdot J_4 \left(\frac{R_0}{p_0}\right)^4 \cdot (64 + 72e_0^2 - (248 + 252e_0^2) \sin{i_0}^2 + (196 + 189e_0^2) \sin{i_0}^4) \\ \end{aligned} \\ & \left. \right] \cdot \left(t - t_0\right) \\ & M(t) = M_0 + \bar{n} \cdot \left(t - t_0\right) \\ \end{aligned}\]

where the subscript $_0$ indicates the initial element, $t_0$ is the initial mean elements' epoch, $a$ is the mean semi-major axis, $n_0$ is mean motion, $\bar{n}$ is the perturbed mean motion, $e$ is the eccentricity, $p$ is the semi-latus rectum, $R_0$ is the Earth's equatorial radius, $\Omega$ is the RAAN, and $M$ is the mean anomaly.

Warning

This propagator provides the same results as the ANSYS STK® for all Keplerian elements except for the right ascension of the ascending node (RAAN). In this case, we would need to flip the sign of the $J_4$ perturbation term to obtain similar values. However, this modification does not seem right given the references [1] and [2], and the version implemented here leads to lower errors when fitting orbits. Hence, further investigation is required. For more information, see this issue.

Initialization

We can initialize the J4 analytical orbit propagator with the following function:

Propagators.init(Val(:J4), orb₀::KeplerianElements; kwargs...) -> OrbitPropagatorJ4

which creates a J4 propagator structure OrbitPropagatorJ4 with the mean Keplerian elements orb₀.

The following keyword selects the gravitational constants for the propagation algorithm:

  • j4c::J4PropagatorConstants: J4 orbit propagator constants (see J4PropagatorConstants). (Default = j4c_egm2008)

This package contains some pre-built propagation constants for this propagator:

J4 Propagator ConstantDescriptionType
j4c_egm2008EGM-2008 gravitational constantsFloat64
j4c_egm2008_f32EGM-2008 gravitational constantsFloat32
j4c_egm1996EGM-1996 gravitational constantsFloat64
j4c_egm1996_f32EGM-1996 gravitational constantsFloat32
j4c_jgm02JGM-02 gravitational constantsFloat64
j4c_jgm02_f32JGM-02 gravitational constantsFloat32
j4c_jgm03JGM-03 gravitational constantsFloat64
j4c_jgm03_f32JGM-03 gravitational constantsFloat32
Note

The type used in the propagation will be the same as used to define the constants in the structure j4c.

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(:J4), orb; j4c = j4c_jgm03)OrbitPropagatorJ4{Float64, Float64}: Propagator name : J4 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{:J4}, 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 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].

It returns the fitted Keplerian elements and the final covariance matrix of the least-square algorithm.

Note

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.

The following keywords are available to configure the fitting process:

  • atol::Number: Tolerance for the residue absolute value. If the residue is lower than atol 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 than rtol, the computation loop stops. (Default = 2e-4)
  • initial_guess::Union{Nothing, KeplerianElements}: Initial guess for the mean elements fitting process. If it is nothing, the algorithm will obtain an initial estimate from the osculating elements in vr_i and vv_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 than jacobian_perturbation_tol, we increase it until it absolute value is higher than jacobian_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: If true, the algorithm prints debugging information to stdout. (Default = true)
  • weight_vector::AbstractVector: Vector with the measurements weights for the least-square algorithm. We assemble the weight matrix W as a diagonal matrix with the elements in weight_vector at its diagonal. (Default = @SVector(ones(Bool, 6)))
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 = Propagators.fit_mean_elements(Val(:J4), 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: 1 53.2399 0.0547615 53239.9 --- PROGRESS: 2 6.08696 0.00629814 6086.96 -88.5669 % PROGRESS: 3 6.27894 0.00643005 6278.95 3.15403 % PROGRESS: 4 6.21357 0.00639258 6213.57 -1.04121 % PROGRESS: 5 6.14231 0.00635162 6142.31 -1.14685 % PROGRESS: 6 6.06484 0.00630732 6064.85 -1.26114 % PROGRESS: 7 5.98081 0.00625956 5980.82 -1.3855 % PROGRESS: 8 5.88987 0.00620821 5889.88 -1.52054 % PROGRESS: 9 5.79173 0.00615322 5791.73 -1.66635 % PROGRESS: 10 5.68617 0.00609456 5686.17 -1.82253 % PROGRESS: 11 5.57314 0.00603234 5573.14 -1.9879 % PROGRESS: 12 5.45275 0.00596676 5452.76 -2.16006 % PROGRESS: 13 5.32544 0.0058982 5325.45 -2.33476 % PROGRESS: 14 5.19205 0.00582729 5192.05 -2.50488 % PROGRESS: 15 5.05399 0.00575492 5053.99 -2.659 % PROGRESS: 16 4.91352 0.00568241 4913.52 -2.77942 % PROGRESS: 17 4.774 0.00561156 4774.01 -2.83941 % PROGRESS: 18 4.64033 0.00554484 4640.34 -2.79997 % PROGRESS: 19 4.51938 0.00548552 4519.38 -2.60666 % PROGRESS: 20 4.42046 0.00543786 4420.47 -2.18861 % PROGRESS: 21 4.35577 0.0054073 4355.77 -1.46358 % PROGRESS: 22 4.33863 0.00539962 4338.63 -0.393498 % PROGRESS: 23 4.33863 0.00539961 4338.63 -2.90694e-06 % (KeplerianElements{Float64, Float64}: Epoch = 2.46003e6 (2023-03-24T18:08:40.388), [0.16604866252575995 0.066435930408688 … -3.8553206810206474e-5 0.0001240320441360566; 0.06643593040643207 0.26633435614589746 … -1.7942563352684816e-5 -1.9568110768822313e-5; … ; -3.855320681002824e-5 -1.7942563355231136e-5 … 4.3972013142494105e-7 -8.092682604708755e-8; 0.0001240320441369116 -1.9568110767091483e-5 … -8.092682604801845e-8 1.2451901450868635e-7])
julia> orbKeplerianElements{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 °

References

  • [1] Vallado, D. A (2013). Fundamentals of Astrodynamics and Applications. 4th ed. Microcosm Press, Hawthorn, CA, USA.
  • [2] Kozai, Y (1959). The Motion of a Close Earth Satellite. The Astronomical Journal, v. 64, no. 1274, pp. 367 – 377.