Constellations
The Constellation class provides factory methods for designing satellite
constellations using several well-known algorithms.
Walker Delta
A Walker Delta constellation distributes satellites evenly across orbital planes with a RAAN spread of 360 degrees.
import lox_space as lox
epoch = lox.Time("2025-01-01T00:00:00.000 TAI")
constellation = lox.Constellation.walker_delta(
"iridium",
epoch,
lox.Origin("Earth"),
nsats=66,
nplanes=6,
semi_major_axis=7159 * lox.km,
inclination=53 * lox.deg,
phasing=1,
)
print(len(constellation)) # 66
print(constellation.name) # "iridium"
Walker Star
Same as Walker Delta but with a RAAN spread of 180 degrees.
constellation = lox.Constellation.walker_star(
"polar",
epoch,
lox.Origin("Earth"),
nsats=8,
nplanes=4,
semi_major_axis=7000 * lox.km,
inclination=90 * lox.deg,
)
Street-of-Coverage
Optimizes satellite placement for continuous coverage using the method of Huang, Colombo, and Bernelli-Zazzera (2021).
constellation = lox.Constellation.street_of_coverage(
"coverage",
epoch,
lox.Origin("Earth"),
nsats=24,
nplanes=4,
semi_major_axis=7159 * lox.km,
inclination=53 * lox.deg,
coverage_fold=1,
)
Flower
Flower constellations produce repeating ground tracks. The orbital shape can be computed from a perigee altitude (the mean radius, gravitational parameter, and rotation rate are derived from the origin automatically) or provided directly as semi-major axis and eccentricity.
constellation = lox.Constellation.flower(
"flower14",
epoch,
lox.Origin("Earth"),
n_petals=14,
n_days=1,
nsats=28,
phasing_numerator=1,
phasing_denominator=28,
inclination=53 * lox.deg,
perigee_altitude=780 * lox.km,
)
Using with Scenarios
Constellations can be added to a Scenario for visibility analysis:
scenario = lox.Scenario(t0, t1).with_constellation(constellation)
Each satellite is converted to a Spacecraft using the propagator specified
at constellation creation time (default: "vallado", also available: "j2").
Constellation
A named collection of satellites produced by a constellation design algorithm.
Use the classmethods to create constellations of different types.
Methods:
-
flower–Create a Flower constellation (repeating ground tracks).
-
street_of_coverage–Create a Street-of-Coverage constellation.
-
walker_delta–Create a Walker Delta constellation (RAAN spread = 360 deg).
-
walker_star–Create a Walker Star constellation (RAAN spread = 180 deg).
Attributes:
-
name(str) –Return the constellation name.
-
satellites(list[ConstellationSatellite]) –Return the list of satellites.
flower
classmethod
flower(
name: str,
time: Time,
origin: Origin,
*,
n_petals: int,
n_days: int,
nsats: int,
phasing_numerator: int,
phasing_denominator: int,
inclination: Angle,
perigee_altitude: Distance | None = None,
semi_major_axis: Distance | None = None,
eccentricity: float | None = None,
argument_of_periapsis: Angle | None = None,
propagator: str = "vallado"
) -> Constellation
Create a Flower constellation (repeating ground tracks).
street_of_coverage
classmethod
street_of_coverage(
name: str,
time: Time,
origin: Origin,
*,
nsats: int,
nplanes: int,
semi_major_axis: Distance,
inclination: Angle,
eccentricity: float = 0.0,
coverage_fold: int = 1,
argument_of_periapsis: Angle | None = None,
propagator: str = "vallado"
) -> Constellation
Create a Street-of-Coverage constellation.
walker_delta
classmethod
walker_delta(
name: str,
time: Time,
origin: Origin,
*,
nsats: int,
nplanes: int,
semi_major_axis: Distance,
inclination: Angle,
eccentricity: float = 0.0,
phasing: int = 0,
argument_of_periapsis: Angle | None = None,
propagator: str = "vallado"
) -> Constellation
Create a Walker Delta constellation (RAAN spread = 360 deg).
walker_star
classmethod
walker_star(
name: str,
time: Time,
origin: Origin,
*,
nsats: int,
nplanes: int,
semi_major_axis: Distance,
inclination: Angle,
eccentricity: float = 0.0,
phasing: int = 0,
argument_of_periapsis: Angle | None = None,
propagator: str = "vallado"
) -> Constellation
Create a Walker Star constellation (RAAN spread = 180 deg).
ConstellationSatellite
A single satellite in a constellation.
Attributes:
-
index_in_plane(int) –Return the index within the plane (0-based).
-
plane(int) –Return the orbital plane index (0-based).