Library

Documentation for SatelliteAnalysis.jl.

SatelliteAnalysis._F_and_∂F_l0pMethod
_F_and_∂F_l0p(l::Integer, p::Integer, i::Number) -> BigFloat, BigFloat

Compute the inclination function F_{l,0,p}(i) and its derivative ∂F_{l,0,p} / ∂ias defined in [1, p. 642].

Note

Internally, we must use BigFloat to allow calculations on high degrees.

References

  • [1] Vallado, D. A (2013). Fundamentals of Astrodynamics and Applications. 4th ed. Microcosm Press, Hawthorne, CA.
source
SatelliteAnalysis.beta_angleMethod
beta_angle(orb::KerplerianElements{Tepoch, T}, Δjd::Number; kwargs...) -> Float64

Compute the beta angle [rad] for the orbit orb after Δjd days from its epoch.

The algorithm was obtained from [1].

Note

It is expected that the input elements are represented in the TOD reference frame. If it is not the case, they can be converted using the function orb_eci_to_eci of SatelliteToolboxTransformations.jl.

Keywords

  • perturbation::Symbol: Select the perturbation terms that must be used when propagating the right ascencion of the ascending node. The possible values are:
    • :J0: Consider a Keplerian orbit.
    • :J2: Consider the perturbation terms up to J₂.
    • :J4: Consider the perturbation terms J₂, J₂², and J₄.
    (Default: :J2)

References

  • [1]: Mortari, D., Wilkins, M. P., and Bruccoleri, C. On Sun-Synchronous Orbits and Associated Constellations

Extended Help

The beta angle is the angle between the orbit plane and the Sun. The positive direction is defined as that of the orbit angular momentum.

The beta angle is useful when computing the mean amount of solar radiation a satellite receives in a particular orbit.

Examples

julia> using SatelliteAnalysis, UnicodePlots

julia> jd₀ = date_to_jd(2021, 1, 1, 0, 0, 0)

julia> orb = KeplerianElements(
            jd₀,
            7130.982e3,
            0.001111,
            98.405 |> deg2rad,
            ltdn_to_raan(10.5, jd₀),
            90     |> deg2rad,
            0
        )

julia> β = beta_angle.(orb, 0:1:364)

julia> lineplot(0:1:364, rad2deg.(β), xlabel = "Day", ylabel = "Beta angle [°]")
                     ┌────────────────────────────────────────┐
                  30 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⣠⠞⠉⠙⢦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⡴⠁⠀⠀⠀⠀⠹⡄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠘⡄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡰⠃⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠹⡄⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⢀⡼⠁⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠹⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⣀⣀⣀⠀⠀⠀⠀⠀⠀⣠⠞⠀⠀⠀⠀⠀⠀│
   Beta angle [°]    │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠱⡀⠀⠀⠀⠀⠀⠀⠀⠀⣠⠖⠋⠀⠀⠀⠉⠓⠲⠤⠴⠚⠁⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠹⣄⠀⠀⠀⠀⠀⣠⠞⠁⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠦⣄⣀⡠⠞⠁⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                  10 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
                     └────────────────────────────────────────┘
                     ⠀0⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀400⠀
                     ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀Day⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀
source
SatelliteAnalysis.design_sun_sync_ground_repeating_orbitMethod
design_sun_sync_ground_repeating_orbit(minimum_repetition::Int, maximum_repetition::Int; kwargs...) -> DataFrame

List all the Sun synchronous, ground repeating orbits in which their repetition period is in the interval [minimum_repetition, maximum_repetition] days.

This function returns a DataFrame with the following columns:

  • semi_major_axis: Orbit semi-major axis.
  • altitude: Orbit altitude above the Equator (a - R0).
  • inclination: Orbit inclination.
  • period: Orbital period.
  • rev_per_days: If the keyword pretify_rev_per_days is false, this column contains Tuples with the integer and rational parts of the number of revolutions per day. Otherwise, it contains a string with a prety representation of the number of revolutions per day.
  • adjacent_gt_distance: Distance between two adjacent ground tracks at Equator.
  • adjacent_gt_angle: Angle between two adjacent ground tracks at Equator measured from the satellite position.
Note

The units of those columns depends on the keywords.

Keywords

  • angle_unit::Symbol: Unit for all the angles in the output DataFrame. It can be :deg for degrees or :rad for radians. (Default: :deg)
  • distance_unit::Symbol: The unit for all the distances in the output DataFrame. It can be :m for meters or :km for kilometers. (Default: :km)
  • eccentricity::Number: Orbit eccentricity. (Default: 0)
  • int_rev_per_day::Tuple: Tuple with the integer parts of the number of revolutions per day to be analyzed. (Default = (13, 14, 15, 16, 17))
  • pretity_rev_per_days::Bool: If true, the column with the revolutions per day will be conveted to a string with a pretty representation of this information. (Default: true)
  • maximum_altitude::Union{Nothing, Number}: Maximum altitude [m] of the orbits in the output DataFrame. If it is nothing, the algorithm will not apply a higher limit to the orbital altitude. (Default = nothing)
  • minimum_altitude::Union{Nothing, Number}: Minimum altitude [m] of the orbits in the output DataFrame. If it is nothing, the algorithm will not apply a lower limit to the orbital altitude. (Default = nothing)
  • time_unit::Symbol: Unit for all the time values in the output DataFrame. It can be :s for seconds, :m for minutes, or :h for hours. (Default = :h)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default = GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default = EGM_2008_J2)
  • R0::Number: Earth's equatorial radius [m]. (Default = EARTH_EQUATORIAL_RADIUS)
  • we::Number: Earth's angular speed [rad / s]. (Default: EARTH_ANGULAR_SPEED)
source
SatelliteAnalysis.eclipse_time_summaryMethod
eclipse_time_summary(orbp::OrbitPropagator; kwargs...) -> DataFrame

Compute the eclipse time summary for the orbit propagator orbp. The summary is computed as the total time the object stays in the sunlight, penumbra, and umbra regions per orbit at each day.

Keywords

  • num_days::Number: Number of days in which the analysis will be performed. (Default = 365)
  • step::Number: The step in which the propagation will occur. Notice that this function has a crossing estimation to accurately estimate the transition between the regions. However, if this step is very large, we may miss some small regions. If it is negative, it will be selected as the time in which the mean anomaly advances 0.5°. (Default = -1)
  • unit::Symbol: Select the unit in which the results will be generated. The possible values are:
    • :s for seconds (Default);
    • :m for minutes; or
    • :h for hours.

Returns

  • DataFrame: The function returns a DataFrame with three columns:
    • sunlight: Total sunlight time per orbit at each day [unit].
    • penumbra: Total penumbra time per orbit at each day [unit].
    • umbra: Total umbra time per orbit at each day [unit].
    The unit of each column is stored in the DataFrame using metadata.

Extended Help

Examples

julia> using SatelliteAnalysis

julia> orb = KeplerianElements(
           date_to_jd(2021, 1, 1, 0, 0, 0),
           7130.982e3,
           0.001111,
           98.405 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           90     |> deg2rad,
           0
       )

julia> orbp = Propagators.init(Val(:J2), orb)

julia> df = eclipse_time_summary(orbp; num_days = 5)
5×4 DataFrame
 Row │ date        sunlight  penumbra  umbra
     │ Date        Float64   Float64   Float64
─────┼─────────────────────────────────────────
   1 │ 2021-01-01   3972.63   20.4117  2006.96
   2 │ 2021-01-02   3973.85   20.4376  2005.71
   3 │ 2021-01-03   3974.77   20.4575  2004.77
   4 │ 2021-01-04   3975.74   20.4758  2003.79
   5 │ 2021-01-05   3976.94   20.5022  2002.55

julia> df = eclipse_time_summary(orbp; num_days = 5, unit = :m)
5×4 DataFrame
 Row │ date        sunlight  penumbra  umbra
     │ Date        Float64   Float64   Float64
─────┼─────────────────────────────────────────
   1 │ 2021-01-01   66.2105  0.340195  33.4493
   2 │ 2021-01-02   66.2308  0.340627  33.4285
   3 │ 2021-01-03   66.2461  0.340958  33.4129
   4 │ 2021-01-04   66.2623  0.341263  33.3964
   5 │ 2021-01-05   66.2824  0.341704  33.3759

julia> colmetadata(df)
Dict{Symbol, Dict{String, Symbol}} with 3 entries:
  :penumbra => Dict("Unit"=>:m)
  :sunlight => Dict("Unit"=>:m)
  :umbra    => Dict("Unit"=>:m)
source
SatelliteAnalysis.fetch_country_polygonsFunction
fetch_country_polygons(url = "https://pkgstore.datahub.io/core/geo-countries/countries/archive/23f420f929e0e09c39d916b8aaa166fb/countries.geojson"; kwargs...) -> String

Fetch the GeoJSON file with the country polygons in url. The algorithm stores the file in a scratch space. The function returns a String with the file path.

Note

If the file has already been downloaded, this function only returns its path. However, if the keyword force_download is true, the file is downloaded again from url.

Keywords

  • force_download::Bool: Download the file from url even if it already exists. (Default = false)
source
SatelliteAnalysis.find_crossingMethod
find_crossing(f::Function, t₀::Number, t₁::Number, s₀, s₁, vargs...; Δ = 1e-3, max = 100, kwargs...) -> T

Return the crossing time tc in which the function f(t) goes from the state s₀ to the state s₁. It is assumed that f(t₀) = s₀ and f(t₁) = s₁.

The parameters in vargs... are passed to the function f after t, and the keywords kwargs... are also passed to f. Hence, it will always be called as f(t, vargs...; kwargs...).

If the computed interval is smaller than Δ, or if the number of iterations is higher than max, the algorithm stops.

Note

The output type T is obtained by the type of (t₁ + t₂) / 2.

Examples

julia> SatelliteAnalysis.find_crossing(
    t -> (sin(t) > 0),
    -0.3,
    0.3,
    false,
    true;
    Δ = 1e-10
)
6.984919309616089e-11
source
SatelliteAnalysis.frozen_orbitMethod
frozen_orbit(a::Number, i::Number; kwargs...) -> Float64, Float64

Compute the eccentricity [ ] and argument of perigee [rad] to obtain a frozen orbit when the orbit has semi-major axis a [m] and inclination i [rad]. This function uses the theory in [1].

Note

This function uses BigFloat internally to perform all computations, allowing very high degrees. However, the user must ensure that the default precision is enough for the required degree. Refer to the function setprecision for more information.

Keywords

  • gravity_model::Union{Nothing, AbstractGravityModel}: Gravity model used to compute the frozen eccentricity. Refer to the object AbstractGravityModel of the package SatelliteToolboxGravityModels.jl for more information. If it is nothing, the system will automatically fetch and load the EGM96 gravity model. However, loading a gravity model can significantly decrease the performance. Thus, it is advisable to pass a gravity model here. (Default = nothing)
  • max_degree: Maximum gravity model degree used to compute the frozen eccentricity. If it is equal to or lower than 0, the maximum degree in grav_model will be used. Otherwise, if it is lower than 3 or higher than the grav_model maximum degree, it will be clamped accordingly. (Default = 53)

References

  • [1] Rosborough, G. W.; Ocampo, C. A (1991). Influence of higher degree zonals on the frozen orbit geometry. Proceedings of the AAS/AIAA Astrodynamics Conference, Durango, CO.

Extended Help

Due to the Earth's gravitational perturbation, the orbit of a salite will experience secular changes in the argument of perigee. Hence, the satellite mean altitude per latitude will differ during the mission. This effect can be problematic, especially if we must compare images by a camera onboard the satellite in different periods. The altitude variation will change the resolution, leading to some problems when comparing the data.

We can avoid this problem if we compute an eccentricity and the argument of perigee that yields theoretically:

de      dω
── = 0, ── = 0
dt      dt

This orbit is called frozen. Refer to [1] for more information.

Examples

julia> using SatelliteAnalysis

julia> frozen_orbit(7130.982e3, 98.410 |> deg2rad)
(0.0011641853028456078, 1.5707963267948966)

julia> jgm3 = GravityModels.load(IcgemFile, fetch_icgem_file(:JGM3))
[ Info: Downloading the ICGEM file 'JGM3.gfc' from 'http://icgem.gfz-potsdam.de/getmodel/gfc/a3375e01a717ac162962138a5e94f10
466b71aa4a130d7f7d5b18ab3d5f90c3d/JGM3.gfc'...
IcgemFile{Float64}:
      Product type : gravity_field
       Model name  : JGM3
  Gravity constant : 3.986004415e14
            Radius : 6.3781363e6
    Maximum degree : 70
            Errors : formal
       Tide system : unknown
              Norm : fully_normalized
         Data type : Float64

julia> frozen_orbit(7130.982e3, 98.410 |> deg2rad; gravity_model = jgm3)
(0.001163484769069545, 1.5707963267948966)
source
SatelliteAnalysis.ground_facility_accessesMethod
ground_facility_accesses(orbp, [(WGS84)]; kwargs...) -> DataFrame

Compute the accesses of a satellite with orbit propagator orbp (see Propagators.init) to the ground facilities defined in the vector [(WGS84)]. The analysis interval begins in the propagator epoch plus initial_time and lasts for duration [s], where both are keywords.

The ground facilities are specified using a vector of tuples with three numbers:

Tuple{T1, T2, T3} where {T1 <: Number, T2 <: Number, T3 <: Number}

containing the WGS84 position of each ground facility [(WGS84)]:

(latitude [rad], longitude [rad], altitude [m])

Those geodetic information are transformed to an ECEF vector using the function geodetic_to_ecef.

Warning

This function computes the accesses using multiple threads. Hence, the function f_eci_to_ecef must be thread safe.

Keywords

  • duration::Number: Duration of the analysis [s]. (Default = 86400)

  • f_eci_to_ecef::Function: Function to convert the orbit propagator position represented in the Earth-centered inertial (ECI) reference frame to the Earth-centered, Earth-fixed (ECEF) reference frame. The signature must be

    f_eci_to_ecef(r_i::AbstractVector, jd::Number) -> AbstractVector

    and it must return the position vector r_i represented in the ECEF at the instant jd [Julian Day]. By default, we use TEME as the ECI and PEF as the ECEF. (Default: _ground_facility_default_eci_to_ecef)

  • initial_time::Number: Initial time of the analysis after the propagator epoch [s]. (Default = 0)

  • minimum_elevation::Number: Minimum elevation angle for communication between the satellite and the ground facilities [rad]. (Default = 10°)

  • num_chunks::Number: Number of chunks the algorithm will divide the time vector to compute the accesses. (Default = Threads.nthreads())

  • reduction::Function: A function that receives a boolean vector with the visibility between the satellite and each ground facility. It must return a boolean value indicating if the access must be computed or not. This is useful to merge access time between two or more facilities. (Default = v -> |(v...) i.e. compute the access if at least one ground facilities is visible)

  • step::Number: The step [s] used to propagate the orbit. Notice that we perform a cross tuning to accurately obtain the access time. However, if an access is lower than the step, it can be neglected. (Default = 60)

  • unit::Symbol: Select the unit in which the duration will be computed. The possible values are:

    • :s for seconds (Default);
    • :m for minutes; or
    • :h for hours.

Returns

  • DataFrame: The function returns a DataFrame with three columns:
    • access_beginning: Time of the access beginning [UTC] encoded using DateTime.
    • access_end: Time of the access end [UTC] encoded using DateTime.
    • duration: Duration of the access [s].
    The unit of the column duration is stored in the DataFrame using metadata.

Extended Help

Examples

julia> using SatelliteAnalysis

julia> jd₀ = date_to_jd(2024, 1, 1);

julia> orb = KeplerianElements(
           jd₀,
           7130.982e3,
           0.001111,
           98.405 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           π / 2,
           0
       );

julia> orbp = Propagators.init(Val(:J2), orb);

julia> ground_facility_accesses(orbp, (0, 0, 0))
2×3 DataFrame
 Row │ access_beginning         access_end               duration 
     │ DateTime                 DateTime                 Float64  
─────┼────────────────────────────────────────────────────────────
   1 │ 2024-01-01T10:20:03.136  2024-01-01T10:30:02.971   599.835
   2 │ 2024-01-01T22:49:55.910  2024-01-01T22:59:23.470   567.56

julia> ground_facility_accesses(orbp, (0, 0, 0); unit = :m)
2×3 DataFrame
 Row │ access_beginning         access_end               duration 
     │ DateTime                 DateTime                 Float64  
─────┼────────────────────────────────────────────────────────────
   1 │ 2024-01-01T10:20:03.136  2024-01-01T10:30:02.971   9.99725
   2 │ 2024-01-01T22:49:55.910  2024-01-01T22:59:23.470   9.45933
source
SatelliteAnalysis.ground_facility_gapsMethod
ground_facility_gaps(orbp, args...; duration::Number = 86400, initial_time::Number = 0, kwargs...) -> DataFrame

Compute the gaps between the accesses of ground facilities. The arguments and keywords are the same as the ones used in the function ground_facility_accesses.

Notice that the gap analysis starts in the orbit propagator epoch plus initial_time and lasts for duration [s].

Returns

  • DataFrame: The function returns a DataFrame with three columns:
    • gap_beginning: Time of the access beginning [UTC] encoded using DateTime.
    • gap_end: Time of the access end [UTC] encoded using DateTime.
    • duration: Duration of the access [s].
    The unit of the column duration is stored in the DataFrame using metadata.

Extended Help

Examples

julia> using SatelliteAnalysis

julia> jd₀ = date_to_jd(2024, 1, 1);

julia> orb = KeplerianElements(
           jd₀,
           7130.982e3,
           0.001111,
           98.405 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           π / 2,
           0
       );

julia> orbp = Propagators.init(Val(:J2), orb);

julia> ground_facility_gaps(orbp, (0, 0, 0))
3×3 DataFrame
 Row │ gap_beginning            gap_end                  duration 
     │ DateTime                 DateTime                 Float64  
─────┼────────────────────────────────────────────────────────────
   1 │ 2024-01-01T00:00:00      2024-01-01T10:20:03.136  37203.1
   2 │ 2024-01-01T10:30:02.971  2024-01-01T22:49:55.910  44392.9
   3 │ 2024-01-01T22:59:23.470  2024-01-02T00:00:00       3636.53

julia> ground_facility_gaps(orbp, (0, 0, 0); unit = :m)
3×3 DataFrame
 Row │ gap_beginning            gap_end                  duration 
     │ DateTime                 DateTime                 Float64  
─────┼────────────────────────────────────────────────────────────
   1 │ 2024-01-01T00:00:00      2024-01-01T10:20:03.136  620.052
   2 │ 2024-01-01T10:30:02.971  2024-01-01T22:49:55.910  739.882
   3 │ 2024-01-01T22:59:23.470  2024-01-02T00:00:00       60.6088
source
SatelliteAnalysis.ground_facility_visibility_circleMethod
ground_facility_visibility_circle(gf_wgs84::Tuple, satellite_position_norm::Number; kwargs...) -> Vector{NTuple{2, Float64}}

Compute the ground facility visibility circle from the position gf_wgs84 (WGS84) to a satellite in which its distance from the Earth's center is satellite_position_norm [m]. It returns a vector of NTuple{2, Float64} where the first element is the latitude [rad] and the second is the longitude [rad] of each point in the visibility circle.

The ground facility is specified using a tuple with its WGS84 position:

(latitude [rad], longitude [rad], altitude [m])

Keywords

  • azimuth_step::Number: The step in the azimuth used to compute the visibility circle. (Default: 0.1 |> deg2rad)
  • minimum_elevation::Number: Minimum elevation angle for communication between the satellite and the ground facility [rad]. (Default: 10 |> deg2rad)

Extended Help

Examples

julia> using SatelliteAnalysis, UnicodePlots

julia> gfv = ground_facility_visibility_circle((0, 0, 0), EARTH_EQUATORIAL_RADIUS + 700e3);

julia> lineplot(last.(gfv) .|> rad2deg, first.(gfv) .|> rad2deg; xlim = (-180, 180), ylim = (-90, 90))
       ┌────────────────────────────────────────┐ 
    90 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡼⠉⡏⢧⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⡧⠤⡧⢼⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⢳⣀⣇⡞⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
   -90 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       └────────────────────────────────────────┘ 
       ⠀-180⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀180⠀ 

julia> gfv = ground_facility_visibility_circle((-40 |> deg2rad, -60 |> deg2rad, 700), EARTH_EQUATORIAL_RADIUS + 700e3);

julia> lineplot(last.(gfv) .|> rad2deg, first.(gfv) .|> rad2deg; xlim = (-180, 180), ylim = (-90, 90))
       ┌────────────────────────────────────────┐ 
    90 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⡧⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤⠤│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⣠⠖⠒⢦⡀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⢰⠃⠀⠀⠀⢹⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠸⣄⠀⠀⠀⣸⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠙⠒⠋⠁⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
   -90 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│ 
       └────────────────────────────────────────┘ 
       ⠀-180⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀180⠀ 
source
SatelliteAnalysis.ground_repeating_orbit_adjacent_track_angleMethod
ground_repeating_orbit_adjacent_track_angle(a::T1, e::T2, i::T3, orbit_cycle::Integer; kwargs...) where {T1 <: Number, T2 <: Number, T3 <: Number}

Compute the adjacent track angle [rad] at Equator in a ground repeating orbit measured from the satellite position. The orbit is described by its semi-major axis a [m], eccentricity [ ], inclination i [rad], and orbit cycle orbit_cyle [day].

Warning

The code does not check if the orbit is ground-repeating with orbit_cycle [day].

Note

Internally, this function uses the precision obtained by promoting T1, T2, and T3 to a float-pointing number T.

Keywords

  • perturbation::Symbol: Symbol to select the perturbation terms that will be used. It can be :J0, :J2, or :J4. (Default: :J2)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default: GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default: EGM_2008_J2)
  • J4::Number: J₄ perturbation term. (Default: EGM_2008_J4)
  • R0::Number: Earth's equatorial radius [m]. (Default: EARTH_EQUATORIAL_RADIUS)
  • we::Number: Earth's angular speed [rad / s]. (Default: EARTH_ANGULAR_SPEED)

Extended help

A ground repeating orbit is any orbit that the number of revolutions per day is a rational number. Hence, this type of orbit repeats its ground trace after a finite number of days.

The information orbit_cyle is redundant given that we have a, e, and i. However, it is necessary to improve the algorithm precision. Otherwise, the orbit_cycle must be obtained by computing the orbit period using a, e, and i and then converting it to a rational number, leading to numerical problems.

source
SatelliteAnalysis.ground_repeating_orbit_adjacent_track_distanceMethod
ground_repeating_orbit_adjacent_track_distance(orbit_period::T1, i::T2, orbit_cycle::Integer; kwargs...) where {T1 <: Number, T2 <: Number} -> T

Compute the adjacent track distance [m] at Equator in a ground repeating orbit. The orbit is described by its orbital period orbit_period [s], inclination i [rad], and orbit cycle orbit_cycle [day].

Note

Internally, this function uses the precision obtained by promoting T1 and T2 to a float-pointing number T.

Keywords

  • perturbation::Symbol: Symbol to select the perturbation terms that will be used. It can be :J0, :J2, or :J4. (Default: :J2)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default: GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default: EGM_2008_J2)
  • J4::Number: J₄ perturbation term. (Default: EGM_2008_J4)
  • R0::Number: Earth's equatorial radius [m]. (Default: EARTH_EQUATORIAL_RADIUS)
  • we::Number: Earth's angular speed [rad / s]. (Default: EARTH_ANGULAR_SPEED)

Extended help

A ground repeating orbit is any orbit that the number of revolutions per day is a rational number. Hence, this type of orbit repeats its ground trace after a finite number of days.

The information orbit_period and orbit_cyle is redundant. However, they are necessary to improve the algorithm precision. Otherwise, the orbit_cycle must be obtained by converting the floating-point number orbit_period to a rational number, leading to numerical problems.

source
SatelliteAnalysis.ground_trackMethod
ground_track(orbp::OrbitPropagator; kwargs...) -> Vector{NTuple{2, Float64}}

Compute the satellite ground track using the orbit propagator orbp. It returns a vector of NTuple{2, Float64} where the first element is the latitude [rad] and the second is the longitude [rad] of each point in the ground track.

Keywords

  • add_nans::Bool: If true, we add NaN if there is a discontinuity in the ground track to improve plotting. (Default: true)

  • duration::Number: Duration of the analysis. (Default: 86400)

  • initial_time::Number: Initial time regarding the orbit propagator orbp epoch [s]. (Default: 0)

  • f_eci_to_ecef::Function: Function to convert the orbit propagator position represented in the Earth-centered inertial (ECI) reference frame to the Earth-centered, Earth-fixed (ECEF) reference frame. The signature must be

    f_eci_to_ecef(r_i::AbstractVector, jd::Number) -> AbstractVector

    and it must return the position vector r_i represented in the ECEF at the instant jd [Julian Day]. By default, we use TEME as the ECI and PEF as the ECEF. (Default: _ground_track_default_eci_to_ecef)

  • step::Union{Nothing, Number}: Step for the computation. If nothing, we will roughly compute the step to approximate 1° in the mean anomaly. (Default: nothing)

  • track_types::Symbol: A symbol describing what kind of track types we must add to the output vector. It can be :ascending for only ascending passages, :descending for only descending passages, or :all for both. (Default: :all)

Extended Help

Examples

julia> using SatelliteAnalysis, UnicodePlots

julia> jd₀ = date_to_jd(2024, 1, 1);

julia> orb = KeplerianElements(
           jd₀,
           7130.982e3,
           0.001111,
           98.405 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           π / 2,
           0
       );

julia> orbp = Propagators.init(Val(:J2), orb);

julia> gt = ground_track(orbp);

julia> lineplot(last.(gt), first.(gt))
      ┌────────────────────────────────────────┐
    2 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⣸⡿⠿⣻⠿⣿⠿⣟⡯⢟⡿⢿⣿⡿⣿⡟⣿⡿⢿⣿⢿⣿⢻⣿⠿⣿⡟⣿⡿⣿⡿⢿⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠞⠹⡼⠙⣶⠋⢶⠋⢣⠏⢳⡞⠀⣿⡁⡽⡇⣹⡇⢨⢏⢈⢯⠀⣿⡀⡽⡁⣹⡇⢸⣯⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢦⢰⢳⢠⢿⡀⡟⡄⡞⡇⡸⣇⢸⠁⢷⠃⣷⡇⢸⡏⠸⡼⠈⣾⠀⣷⠃⢣⠇⢹⣏⠏⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢸⡞⠘⣾⠀⣿⠁⢧⠃⢹⡇⢸⡏⠀⣼⠀⣿⡆⢰⡇⢀⣇⠀⣷⠀⣾⠀⣸⡀⢸⣿⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢀⡇⠀⣇⠀⣿⠀⢸⠀⢸⡁⢨⡇⠀⡇⡇⡏⡇⣸⢳⢸⢸⢰⢻⣀⡏⡆⡇⡇⡜⣯⠀⠀⠀⠀⠀│
      │⠤⠤⠤⠤⢼⢧⢤⢿⠤⣿⠤⡿⡦⡼⡧⢼⣧⢴⠥⢧⡧⢼⡧⢼⡾⠼⣼⠤⣿⠤⣧⠧⢷⣧⢿⠤⠤⠤⠤⠤│
      │⠀⠀⠀⠀⡼⢸⢸⠸⣼⠁⣇⠇⡇⡇⢳⡞⢸⢸⠀⢸⡇⢸⡇⠈⡇⠀⡏⠀⣿⠀⢸⠀⢸⢹⠘⡆⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠁⠈⡏⠀⡿⠀⢿⠀⢹⠃⢸⡇⠘⡇⠀⡜⡇⣸⡇⢸⢧⢠⢿⠀⣿⠀⡾⡄⡸⡟⠀⣷⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢆⢰⢳⢀⢿⠀⡿⡄⡼⡆⣸⣇⢰⢧⢠⠇⣧⡇⢹⡞⠸⡼⠘⣾⠁⣷⠃⢧⢧⢷⢀⡟⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢘⣎⢈⣞⠀⣿⡁⣳⡃⣹⡇⢸⣏⢘⣞⢀⡿⢆⡼⢧⣰⢳⣠⠿⣀⠾⣄⡞⣆⢈⣿⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢹⣾⣿⣾⣿⣶⣿⣧⣿⣷⣿⣷⣾⣿⣼⣿⣷⣾⣷⣾⣵⣲⣽⣶⣿⣶⣯⣶⣼⣿⣾⣿⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
   -2 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      └────────────────────────────────────────┘
      ⠀-4⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀4⠀

julia> gt = ground_track(orbp, track_types = :ascending);

julia> lineplot(last.(gt), first.(gt))
      ┌────────────────────────────────────────┐
    2 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢸⡛⠽⡛⠿⣟⠫⣟⠫⢍⠛⠻⢿⡻⢽⡛⡿⡛⠿⣟⠯⣟⠫⢟⠻⢿⡛⢽⡛⠽⡛⠯⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠹⡄⠘⡆⠈⢆⠈⢣⠀⢳⡀⠀⢳⡀⠹⡇⠙⡆⠈⢆⠈⢧⠀⢳⡀⠱⡀⠙⡄⠘⣆⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢦⠀⢱⠀⢸⡀⠘⡄⠈⡇⠀⣇⠀⠀⢧⠀⣷⠀⢸⠀⠸⡄⠈⡆⠀⣇⠀⢣⠀⢹⠀⠈⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢸⡀⠘⡆⠀⡇⠀⢇⠀⢹⠀⢸⠀⠀⢸⠀⡟⡆⠀⡇⠀⣇⠀⢳⠀⢸⠀⠸⡀⠈⡇⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⡇⠀⣇⠀⢳⠀⢸⠀⠸⡀⠈⡇⠀⠀⡇⡇⡇⠀⢳⠀⢸⠀⢸⡀⠈⡆⠀⡇⠀⢧⠀⠀⠀⠀⠀│
      │⠤⠤⠤⠤⠤⢧⠤⢼⠤⢼⠤⠼⡦⠤⡧⠤⣧⠤⠤⢧⡧⢼⠤⢼⠤⠼⡤⠤⡧⠤⡧⠤⢷⠤⢼⠤⠤⠤⠤⠤│
      │⠀⠀⠀⠀⠀⢸⠀⠸⡄⠀⡇⠀⡇⠀⢳⠀⢸⠀⠀⢸⡇⠸⡄⠈⡇⠀⡇⠀⢧⠀⢸⠀⢸⠀⠘⡆⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠈⡇⠀⡇⠀⢧⠀⢹⠀⢸⡀⠘⡆⠀⠘⡇⠀⡇⠀⢧⠀⢹⠀⢸⠀⠸⡄⠈⡇⠀⣇⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢆⠀⢳⠀⢸⠀⠸⡄⠘⡆⠀⣇⠀⢧⠀⠀⣧⠀⢹⠀⠸⡀⠘⡆⠀⣇⠀⢧⠀⢳⠀⠘⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠘⣆⠈⢆⠀⢧⡀⢳⡀⠹⡄⠘⣄⠘⣆⠀⡏⢆⠀⢧⠀⢳⡀⠹⡀⠸⣄⠘⣆⠈⢧⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠐⣬⣗⣮⣷⣦⣵⣢⣽⣲⣽⣶⣬⣗⣬⣗⣧⣬⣓⣦⣵⣢⣽⣲⣽⣶⣬⣖⣬⣗⣮⡵⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
   -2 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      └────────────────────────────────────────┘
      ⠀-4⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀4⠀

julia> gt = ground_track(orbp, track_types = :descending);

julia> lineplot(last.(gt), first.(gt))
      ┌────────────────────────────────────────┐
    2 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⣸⡽⠟⣻⠽⣻⠿⢛⡯⢛⡯⢟⡻⠝⣻⠝⣿⠽⢛⠯⢛⡯⢛⡿⠟⡻⠝⣻⠽⣻⠿⢓⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠞⠀⡼⠁⣰⠃⢰⠋⢠⠏⢀⡞⠀⡜⠁⡼⡇⣰⠃⢠⠏⢀⠎⠀⡞⠀⡼⠁⣰⠃⢰⣫⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⢰⠃⢠⠇⠀⡏⠀⡞⠀⡸⠀⢸⠁⢰⠃⣇⡇⠀⡏⠀⡼⠀⣸⠀⢰⠃⢠⠇⠀⣏⠇⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⡞⠀⣸⠀⢸⠁⢠⠃⠀⡇⠀⡏⠀⡼⠀⣿⠀⢰⠁⢀⡇⠀⡇⠀⡞⠀⣸⠀⢸⢸⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢀⡇⠀⡇⠀⡜⠀⢸⠀⢸⠁⢠⠇⠀⡇⠀⡏⠀⣸⠀⢸⠀⢰⠃⢀⡇⠀⡇⠀⡜⡏⠀⠀⠀⠀⠀│
      │⠤⠤⠤⠤⢼⠤⢤⠧⠤⡧⠤⡯⠤⡼⠤⢼⠤⢴⠥⢤⡧⠤⡧⠤⡾⠤⢼⠤⢼⠤⢤⠧⠤⣧⠧⠤⠤⠤⠤⠤│
      │⠀⠀⠀⠀⡼⠀⢸⠀⢸⠁⢀⠇⠀⡇⠀⡞⠀⢸⠀⢸⡇⢠⠃⠀⡇⠀⡏⠀⡼⠀⢸⠀⢸⢹⠀⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠁⠀⡎⠀⡸⠀⢸⠀⢰⠃⢀⡇⠀⡇⠀⡜⡇⣸⠀⢸⠁⢠⠇⠀⡇⠀⡎⠀⡸⡞⠀⣰⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⢰⠃⢀⠇⠀⡏⠀⡼⠀⣸⠀⢰⠁⢠⠇⣇⡇⠀⡞⠀⡼⠀⢸⠁⢰⠃⢀⢧⠇⢀⡇⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢀⠎⢀⡞⠀⡼⠁⡰⠃⣠⠃⢠⠏⢀⡞⢀⡿⠀⡼⠁⣰⠃⢠⠇⢀⠎⢀⡞⠀⢀⡜⠀⠀⠀⠀⠀│
      │⠀⠀⠀⠀⢩⣖⣯⣔⣮⣴⣾⣥⣺⣥⣲⣥⣖⣯⣔⣯⣷⣾⣵⣺⣥⣲⣥⣶⣯⣖⣋⣤⣔⣯⣔⣎⠀⠀⠀⠀│
      │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
   -2 │⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⡇⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀│
      └────────────────────────────────────────┘
      ⠀-4⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀4⠀
source
SatelliteAnalysis.ground_track_inclinationMethod
ground_track_inclination(a::Number, e::Number, i::Number; kwargs...) -> T
ground_track_inclination(orb::Orbit{Tepoch, T}); kwargs...) where {Tepoch <: Number, T <: Number} -> T

Compute the ground track inclination at the Equator [rad] in an orbit with semi-major axis a [m], eccentricity e [ ], and inclination i [rad]. The orbit can also be specified by orb (see Orbit).

Note

The output type T in the first signature is obtained by promoting the inputs to a float type.

Warning

The algorithm here assumes a small orbit eccentricity.

Keywords

  • perturbation::Symbol: Symbol to select the perturbation terms that will be used. It can be :J0, :J2, or :J4. (Default: :J2)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default: GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default: EGM_2008_J2)
  • J4::Number: J₄ perturbation term. (Default: EGM_2008_J4)
  • R0::Number: Earth's equatorial radius [m]. (Default: EARTH_EQUATORIAL_RADIUS)
  • we::Number: Earth's angular speed [rad / s]. (Default: EARTH_ANGULAR_SPEED)

Extended Help

We define the ground track inclination as the angle that the ground track has with respect to the Equator. This information is important to compute, for example, the required swath for a remote sensing satellite to cover the entire Earth.

The ground track inclination i_gt is given by:

            ┌                       ┐
            │      ω_s ⋅ sin i      │
i_gt = atan │ ───────────────────── │
            │ ω_s ⋅ cos i - ω_e + Ω̇ │
            └                       ┘

where ω_s = n + ω̇ is the satellite angular velocity, n is the perturbed mean motion, i is the orbit inclination, ω is the orbit argument of perigee, Ω is the orbit right ascension of the ascending node, and ω_e is the Earth's angular rate.

Formally, we should use the satellite instantaneous angular speed at the Equator instead of the mean angular speed ω_s. However, given the perturbations caused by the Earth's gravitational potential, the former is not simple to compute. This calculation would required to implement an orbit propagator here. Thus, we simplify it by assuming that the orbit eccentricity is small. This assumption is reasonable given the missions that would benefit from the computation of the ground track inclination. In this case, the orbit angular speed is almost constant and equal to ω_s.

Examples

julia> using SatelliteAnalysis

julia> ground_track_inclination(7130.982e3, 0.00111, 98.410 |> deg2rad) |> rad2deg
102.30052101661899

julia> jd₀ = date_to_jd(2021, 1, 1)
2.4592155e6

julia> orb = KeplerianElements(
           jd₀,
           7130.982e3,
           0.001111,
           98.410 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           π / 2,
           0
       )
KeplerianElements{Float64, Float64}:
           Epoch :    2.45922e6 (2021-01-01T00:00:00)
 Semi-major axis : 7130.98     km
    Eccentricity :    0.001111
     Inclination :   98.41     °
            RAAN :   78.4021   °
 Arg. of Perigee :   90.0      °
    True Anomaly :    0.0      °

julia> ground_track_inclination(orb) |> rad2deg
102.30052101658998
source
SatelliteAnalysis.is_ground_facility_visibleMethod
is_ground_facility_visible(sat_r_e::AbstractVector, gf_lat::Number, gf_lon::Number, gf_h::Number, θ::Number) -> Bool
is_ground_facility_visible(sat_r_e::AbstractVector, gf_r_e::AbstractVector, θ::Number) -> Bool
is_ground_facility_visible(sat_r_ned::AbstractVector, θ::Number) -> Bool

Check if the satellite with position vector sat_r_e (ECEF) is inside the visibility circle of a ground facility with latitude gf_lat [rad], longitude gf_lon [rad], altitude gf_h (WGS-84) or ECEF position gf_r_e [m]. The algorithm considers that the ground station has visibility to the satellite if its elevation angle is larger than θ [rad].

The user can also pass the satellite position represented in the NED (North-East-Down) reference frame sat_r_ned [m] at the ground station location, which increases the performance since the algorithm performs no reference frame conversion.

Returns

  • Bool: true if the satellite is inside the visibility circle, or false otherwise.
source
SatelliteAnalysis.lighting_conditionMethod
lighting_condition(r_i::AbstractVector, s_i::AbstractVector)

Compute the lighting condition at the position r_i [m] considering the Sun position vector s_i [m]. The possible return values are:

  • :sunlight: The point is under direct sunlight.
  • :penumbra: The point is in penumbra region.
  • :umbra: The point is in umbra region.

The algorithm used in this function was based on [1].

Note

The vectors r_i and s_i must be represented in the same reference frame.

References

  • [1] Longo, C. R. O., Rickman, S. L (1995). Method for the Calculation of Spacecraft Umbra and Penumbra Shadow Terminator Points. NASA Technical Paper 3547.
source
SatelliteAnalysis.plot_ground_facility_visibility_circles!Method
plot_ground_facility_visibility_circles!(ax::Axis, vgf_vc::Vector{Vector{NTuple{2, Number}}}; kwargs...) -> Nothing

Plot in the Makie.jl axis ax the ground facility visibility circles in the vector vgf_vc, where each element is computed using the function ground_facility_visibility_circle.

Warning

This function only works after loading the package GeoMakie.jl. Furthermore, the user must also load one Makie.jl back end (CairoMakie.jl or GLMakie.jl, for example) to see the result.

Keywords

  • ground_facility_names::Union{Nothing, Vector{String}}: The user can provide a vector of Strings with the length of vgf_vc to be plotted with the visibility circles. If this parameter is nothing, no ground facility name is added to the figure. (Default = nothing)

Extended Help

Examples

julia> using SatelliteAnalysis, GeoMakie, GLMakie

julia> gfv1 = ground_facility_visibility_circle((0, 0, 0), EARTH_EQUATORIAL_RADIUS + 700e3);

julia> gfv2 = ground_facility_visibility_circle((-40 |> deg2rad, -60 |> deg2rad, 0), EARTH_EQUATORIAL_RADIUS + 700e3);

julia> fig = Figure(size = (1000, 1000))

julia> ax = Axis(fig[1, 1])

julia> plot_ground_facility_visibility_circles!(
           ax,
           [gfv1, gfv2];
           ground_facility_names = ["GF 1", "GF 2"]
       )
source
SatelliteAnalysis.plot_ground_facility_visibility_circlesMethod
plot_ground_facility_visibility_circles(vgf_vc::Vector{Vector{NTuple{2, Number}}}; kwargs...) -> Figure, Axis

Plot the ground facility visibility circles in the vector vgf_vc, where each element is computed using the function ground_facility_visibility_circle. It returns the objects Figure and Axis used to plot the data. For more information, please, refer to Makie.jl documentation.

Note

This function plots the countries' borders in the created figure using the file with the country polygons fetched with the function fetch_country_polygons. Hence, if this files does not exist, the algorithm tries to download it.

Warning

This function only works after loading the package GeoMakie.jl. Furthermore, the user must also load one Makie.jl back end (CairoMakie.jl or GLMakie.jl, for example) to see the result.

Keywords

  • ground_facility_names::Union{Nothing, Vector{String}}: The user can provide a vector of Strings with the length of vgf_vc to be plotted with the visibility circles. If this parameter is nothing, no ground facility name is added to the figure. (Default = nothing)

All other kwargs... are passed to the function plot_world_map.

Extended Help

Examples

julia> using SatelliteAnalysis, GeoMakie, GLMakie

julia> gfv1 = ground_facility_visibility_circle((0, 0, 0), EARTH_EQUATORIAL_RADIUS + 700e3);

julia> gfv2 = ground_facility_visibility_circle((-40 |> deg2rad, -60 |> deg2rad, 0), EARTH_EQUATORIAL_RADIUS + 700e3);

julia> fig, ax = plot_ground_facility_visibility_circles(
           [gfv1, gfv2];
           ground_facility_names = ["GF 1", "GF 2"]
       )
(Scene (1600px, 800px):
  0 Plots
  1 Child Scene:
    └ Scene (1600px, 800px), Axis (7 plots))

julia> fig
source
SatelliteAnalysis.plot_ground_track!Method
plot_ground_track!(ax:Axis, gt::Vector{NTuple{2, Number}}) -> Nothing

Plot in the Makie.jl axis ax the ground track gt computed using the function ground_track.

Warning

This function only works after loading the package GeoMakie.jl. Furthermore, the user must also load one Makie.jl back end (CairoMakie.jl or GLMakie.jl, for example) to see the result.

Extended Help

Examples

julia> using SatelliteAnalysis, GeoMakie, GLMakie

julia> jd₀ = date_to_jd(2021, 1, 1)
2.4592155e6

julia> orb = KeplerianElements(
           jd₀,
           7130.982e3,
           0.001111,
           98.405 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           π / 2,
           0
       )
KeplerianElements{Float64, Float64}:
           Epoch :    2.45922e6 (2021-01-01T00:00:00)
 Semi-major axis : 7130.98     km
    Eccentricity :    0.001111
     Inclination :   98.405    °
            RAAN :   78.4021   °
 Arg. of Perigee :   90.0      °
    True Anomaly :    0.0      °

julia> orbp = Propagators.init(Val(:J2), orb)
OrbitPropagatorJ2{Float64, Float64}:
   Propagator name : J2 Orbit Propagator
  Propagator epoch : 2021-01-01T00:00:00
  Last propagation : 2021-01-01T00:00:00

julia> gt = ground_track(orbp; track_types = :descending, duration = 5 * 86400);


julia> fig = Figure(size = (1000, 1000))

julia> ax = Axis(fig[1, 1])

julia> plot_ground_track!(ax, gt)
source
SatelliteAnalysis.plot_ground_trackMethod
plot_ground_track(gt::Vector{NTuple{2, Number}}; kwargs...) -> Figure, Axis

Plot the ground track gt computed using the function ground_track. It returns the objects Figure and Axis used to plot the data. For more information, please, refer to Makie.jl documentation.

Note

This function plots the countries' borders in the created figure using the file with the country polygons fetched with the function fetch_country_polygons. Hence, if this files does not exist, the algorithm tries to download it.

Warning

This function only works after loading the package GeoMakie.jl. Furthermore, the user must also load one Makie.jl back end (CairoMakie.jl or GLMakie.jl, for example) to see the result.

All kwargs... are passed to the function plot_world_map.

Extended Help

Examples

julia> using SatelliteAnalysis, GeoMakie, GLMakie

julia> jd₀ = date_to_jd(2021, 1, 1)
2.4592155e6

julia> orb = KeplerianElements(
           jd₀,
           7130.982e3,
           0.001111,
           98.405 |> deg2rad,
           ltdn_to_raan(10.5, jd₀),
           π / 2,
           0
       )
KeplerianElements{Float64, Float64}:
           Epoch :    2.45922e6 (2021-01-01T00:00:00)
 Semi-major axis : 7130.98     km
    Eccentricity :    0.001111
     Inclination :   98.405    °
            RAAN :   78.4021   °
 Arg. of Perigee :   90.0      °
    True Anomaly :    0.0      °

julia> orbp = Propagators.init(Val(:J2), orb)
OrbitPropagatorJ2{Float64, Float64}:
   Propagator name : J2 Orbit Propagator
  Propagator epoch : 2021-01-01T00:00:00
  Last propagation : 2021-01-01T00:00:00

julia> gt = ground_track(orbp; track_types = :descending, duration = 5 * 86400);

julia> fig, ax = plot_ground_track(gt; size = (2000, 1000))
(Scene (2000px, 1000px):
  0 Plots
  1 Child Scene:
    └ Scene (2000px, 1000px), Axis (2 plots))

julia> fig
source
SatelliteAnalysis.plot_world_mapMethod
plot_world_map(; kwargs...) -> Figure, Axis

Create a Makie.jl Figure and Axis with the world map. All the kwargs... are passed to the function Figure. For more information, please, refer to Makie.jl documentation.

Note

This function plots the countries' borders in the created figure using the file with the country polygons fetched with the function fetch_country_polygons. Hence, if this files does not exist, the algorithm tries to download it.

source
SatelliteAnalysis.sun_sync_orbit_from_angular_velocityMethod
sun_sync_orbit_from_angular_velocity(angvel::T1, e::T2 = 0; kwargs...) where {T1 <: Number, T2 <: Number} -> T, T, Bool

Compute the Sun-synchronous orbit semi-major axis [m] and inclination [rad] given the angular velocity angvel [rad / s] and the orbit eccentricity e [ ]. If the latter is omitted, the orbit is considered circular, i.e., e = 0.

The algorithm here considers only the perturbation terms up to J₂.

Note

Internally, this function uses the precision obtained by promoting T1 and T2 to a float-pointing number T.

Keywords

  • max_iterations::Number: Maximum number of iterations in the Newton-Raphson method. (Default = 3)
  • no_warnings::Bool: If true, no warnings will be printed. (Default = false)
  • tolerance::Union{Nothing, NTuple{2, Number}}: Residue tolerances to verify if the numerical method has converged. If it is nothing, (√eps(T), √eps(T)) will be used, where T is the internal type for the computations. Notice that the residue function f₁ unit is [deg / day], whereas the f₂ unit is [deg / min]. (Default = 1e-18)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default = GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default = EGM_2008_J2)
  • R0::Number: Earth's equatorial radius [m]. (Default = EARTH_EQUATORIAL_RADIUS)

Returns

  • T: Semi-major axis [m].
  • T: Inclination [rad].
  • Bool: true if the Newton-Raphson algorithm converged, or false otherwise.

Extended help

A Sun-synchronous orbit is defined as an orbit in which the precession of the right ascension of the ascending node (RAAN) equals the Earth's orbit mean motion. In this case, the orbit plane will have the same orientation to the Sun at the ascending node.

The RAAN time-derivative considering only the secular terms up to J₂ is [1, p. 372] is:

∂Ω      3                       n̄
── = - ─── R₀² . J₂ . cos(i) . ─── .
∂t      2                       p²

where:

         ┌                                                ┐
         │      3    R₀²                                  │
n̄ = n₀ . │ 1 + ─── . ─── . J₂ . √(1 - e²) . (2 - 3sin²(i))│.
         │      4     p²                                  │
         └                                                ┘

We can express the orbit angular velocity in terms of its nodal period, i.e., the period it takes for the satellite to cross the ascending node two consecutive times:

          ∂M     ∂ω
angvel = ──── + ────,
          ∂t     ∂t

                  3    R₀²
angvel = n̄ + n̄ . ─── . ─── . J₂ . (4 - 5sin²(i)),
                  4     p²

where n is the perturbed mean motion due to the same consideration as presented for the RAAN time-derivative.

Finally, this function finds the pair (a, i) that simultaneously solves the equations:

∂Ω
── (a, i) = EARTH_ORBIT_MEAN_MOTION,
∂t

 ∂M            ∂ω
──── (a, i) + ──── (a, i) = angvel,
 ∂t            ∂t

using the Newton-Raphson method with the presented equations.

Examples

julia> using SatelliteAnalysis

julia> sun_sync_orbit_from_angular_velocity(0.06 |> deg2rad)
(7.130983932846816e6, 1.7175898375139984, true)

julia> sun_sync_orbit_from_angular_velocity(0.06 |> deg2rad, 0)
(7.130983932846816e6, 1.7175898375139984, true)

julia> sun_sync_orbit_from_angular_velocity(0.06 |> deg2rad, 0.1)
(7.13086251587883e6, 1.7146410689929386, true)

The user can verify some internal information of the solver by turning on the debugging logs:

julia> with_logger(ConsoleLogger(stderr, Logging.Debug)) do
           sun_sync_orbit_from_angular_velocity(0.06 |> deg2rad)
       end
┌ Debug: Iteration #1
│   Estimation :
│     a  = 7136.635455699327 km
│     i  = 81.57099271530629 °
│   Residues :
│     f₁ = 1.9706966881670205 ° / day
│     f₂ = 0.004266929859281898 ° / min
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:394
┌ Debug: Iteration #2
│   Estimation :
│     a  = 7128.856266265137 km
│     i  = 98.46515928332974 °
│   Residues :
│     f₁ = -0.0073785260175135425 ° / day
│     f₂ = -0.0016144146737784304 ° / min
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:394
┌ Debug: Iteration #3
│   Estimation :
│     a  = 7130.983594940013 km
│     i  = 98.41070350863473 °
│   Residues :
│     f₁ = -6.549620124363109e-6 ° / day
│     f₂ = -2.6066638092459016e-7 ° / min
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:394
┌ Debug: Iteration #4
│   Estimation :
│     a  = 7130.983932846698 km
│     i  = 98.41064862339992 °
│   Residues :
│     f₁ = 8.290634845309341e-11 ° / day
│     f₂ = -2.2648549702353193e-14 ° / min
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:394
(7.130983932846816e6, 1.7175898375139984, true)

References

  • [1] Kozai, Y (1959). The Motion of a Close Earth Satellite. The Astronomical Journal, v. 64, no. 1274, pp. 367 – 377.
source
SatelliteAnalysis.sun_sync_orbit_inclinationMethod
sun_sync_orbit_inclination(a::T1, e::T2 = 0; kwargs...) where {T1 <: Number, T2 <: Number} -> T, Bool

Compute the inclination [rad] of the Sun-synchronous orbit with semi-major axis a [m] and the eccentricity e [ ]. If the latter is omitted, the orbit is considered circular, i.e., e = 0.

The algorithm here considers only the perturbation terms up to J₂.

Note

Internally, this function uses the precision obtained by promoting T1 and T2 to a float-pointing number T.

Keywords

  • max_iterations::Number: Maximum number of iterations in the Newton-Raphson method. (Default = 30)
  • tolerance::Union{Nothing, Number}: Residue tolerance to verify if the numerical method has converged. If it is nothing, √eps(T) will be used, where T is the internal type for the computations. Notice that the residue unit is [deg / day]. (Default = nothing)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default = GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default = EGM_2008_J2)
  • R0::Number: Earth's equatorial radius [m]. (Default = EARTH_EQUATORIAL_RADIUS)

Returns

  • T: Inclination [rad] of the Sun-synchronous orbit with semi-major axis a and eccentricity e.
  • Bool: true if the Newton-Raphson algorithm converged, or false otherwise.

Extended help

A Sun-synchronous orbit is defined as an orbit in which the precession of the right ascension of the ascending node (RAAN) equals the Earth's orbit mean motion. In this case, the orbit plane will have the same orientation to the Sun at the ascending node.

The RAAN time-derivative considering only the secular terms up to J₂ is [1, p. 372] is:

∂Ω      3                       n̄
── = - ─── R₀² . J₂ . cos(i) . ─── .
∂t      2                       p²

where:

         ┌                                                ┐
         │      3    R₀²                                  │
n̄ = n₀ . │ 1 + ─── . ─── . J₂ . √(1 - e²) . (2 - 3sin²(i))│.
         │      4     p²                                  │
         └                                                ┘

Finally, this function solves the equation:

∂Ω
── (i) = EARTH_ORBIT_MEAN_MOTION
∂t

for i using the Newton-Raphson method with the presented equations.

Examples

julia> using SatelliteAnalysis

julia> sun_sync_orbit_inclination(7130.982e3)
(1.7175896973066611, true)

julia> sun_sync_orbit_inclination(7130.982e3, 0)
(1.7175896973066611, true)

julia> sun_sync_orbit_inclination(7130.982e3, 0.001111)
(1.7175893324980402, true)

The user can verify some internal information of the solver by turning on the debugging logs:

julia> with_logger(ConsoleLogger(stderr, Logging.Debug)) do
           sun_sync_orbit_inclination(7130.982e3)
       end
┌ Debug: Iteration #1
│   Estimation : 98.41064059121584 °
│   Residue    : 0.0005992085524891833 ° / day
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:686
┌ Debug: Iteration #2
│   Estimation : 98.41064059082426 °
│   Residue    : -4.556321986370904e-11 ° / day
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:686
(1.7175896973066611, true)

References

  • [1] Kozai, Y (1959). The Motion of a Close Earth Satellite. The Astronomical Journal, v. 64, no. 1274, pp. 367 – 377.
source
SatelliteAnalysis.sun_sync_orbit_semi_major_axisMethod
sun_sync_orbit_semi_major_axis(i::T1, e::T2 = 0; kwargs...) where {T1 <: Number, T2 <: Number} -> T, Bool

Compute the semi-major axis [m] of the Sun-synchronous orbit with inclination i [rad] and the eccentricity e [ ]. If the latter is omitted, the orbit is considered circular, i.e., e = 0.

The algorithm here considers only the perturbation terms up to J₂.

Note

Internally, this function uses the precision obtained by promoting T1 and T2 to a float-pointing number T.

Keywords

  • max_iterations::Number: Maximum number of iterations in the Newton-Raphson method. (Default = 30)
  • tolerance::Union{Nothing, Number}: Residue tolerance to verify if the numerical method has converged. If it is nothing, √eps(T) will be used, where T is the internal type for the computations. Notice that the residue unit is [deg / day]. (Default = nothing)
  • m0::Number: Standard gravitational parameter for Earth [m³ / s²]. (Default = GM_EARTH)
  • J2::Number: J₂ perturbation term. (Default = EGM_2008_J2)
  • R0::Number: Earth's equatorial radius [m]. (Default = EARTH_EQUATORIAL_RADIUS)

Returns

  • T: Semi-major axis [m] of the Sun-synchronous orbit with inclination i and eccentricity e.
  • Bool: true if the Newton-Raphson algorithm converged, or false otherwise.

Extended help

A Sun-synchronous orbit is defined as an orbit in which the precession of the right ascension of the ascending node (RAAN) equals the Earth's orbit mean motion. In this case, the orbit plane will have the same orientation to the Sun at the ascending node.

The RAAN time-derivative considering only the secular terms up to J₂ is [1, p. 372] is:

∂Ω      3                       n̄
── = - ─── R₀² . J₂ . cos(i) . ─── .
∂t      2                       p²

where:

         ┌                                                ┐
         │      3    R₀²                                  │
n̄ = n₀ . │ 1 + ─── . ─── . J₂ . √(1 - e²) . (2 - 3sin²(i))│.
         │      4     p²                                  │
         └                                                ┘

Finally, this function solves the equation:

∂Ω
── (a) = EARTH_ORBIT_MEAN_MOTION
∂t

for a using the Newton-Raphson method with the presented equations.

Examples

julia> using SatelliteAnalysis

julia> sun_sync_orbit_semi_major_axis(98.410 |> deg2rad)
(7.130827866508738e6, true)

julia> sun_sync_orbit_semi_major_axis(98.410 |> deg2rad, 0)
(7.130827866508738e6, true)

julia> sun_sync_orbit_semi_major_axis(98.410 |> deg2rad, 0.001111)
(7.1308328955274355e6, true)

The user can verify some internal information of the solver by turning on the debugging logs:

julia> with_logger(ConsoleLogger(stderr, Logging.Debug)) do
           sun_sync_orbit_semi_major_axis(98.41064163374567 |> deg2rad)
       end
┌ Debug: Iteration #1
│   Estimation : 7130.981820550704 km
│   Residue    : 0.0005989504045072862 ° / day
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:559
┌ Debug: Iteration #2
│   Estimation : 7130.982250931794 km
│   Residue    : -2.081337784770338e-7 ° / day
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:559
┌ Debug: Iteration #3
│   Estimation : 7130.982250931845 km
│   Residue    : -2.4312691616901194e-14 ° / day
└ @ SatelliteAnalysis ~/.julia/dev/SatelliteAnalysis/src/sun_synchronous_orbits.jl:559
(7.130982250931845e6, true)

References

  • [1] Kozai, Y (1959). The Motion of a Close Earth Satellite. The Astronomical Journal, v. 64, no. 1274, pp. 367 – 377.
source