opengnc.propagators package

Submodules

opengnc.propagators.base module

Abstract base class for orbit propagators.

class opengnc.propagators.base.Propagator[source]

Bases: ABC

Abstract base class for orbit propagators.

abstractmethod propagate(r_i: ndarray, v_i: ndarray, dt: float, **kwargs: Any) tuple[ndarray, ndarray][source]

Propagates the state vector (position and velocity) forward in time.

Parameters:
  • r_i (np.ndarray) – Initial position vector (m).

  • v_i (np.ndarray) – Initial velocity vector (m/s).

  • dt (float) – Time duration for propagation (s).

  • **kwargs (dict) – Additional arguments specific to the propagator implementation.

Returns:

(r_f, v_f). r_f : Final position vector (m). v_f : Final velocity vector (m/s).

Return type:

tuple[np.ndarray, np.ndarray]

opengnc.propagators.cowell module

Numerical Cowell Propagator.

class opengnc.propagators.cowell.CowellPropagator(integrator: Integrator | None = None, mu: float = 398600441800000.0)[source]

Bases: Propagator

Numerical Cowell Propagator.

Integrates the equations of motion numerically, allowing for perturbations.

Parameters:
  • integrator (Integrator, optional) – Numerical integrator instance. Defaults to RK4.

  • mu (float, optional) – Gravitational parameter (m^3/s^2). Default is Earth’s.

propagate(r_i: ndarray, v_i: ndarray, dt: float, perturbation_acc_fn: Callable[[float, ndarray, ndarray], ndarray] | None = None, **kwargs: Any) tuple[ndarray, ndarray][source]

Propagates the state using numerical integration.

Parameters:
  • r_i (np.ndarray) – Initial position vector (m).

  • v_i (np.ndarray) – Initial velocity vector (m/s).

  • dt (float) – Propagation duration (s).

  • perturbation_acc_fn (callable, optional) – Function that returns perturbation accelerations. Signature: acc_pert = f(t, r, v) -> np.ndarray.

  • **kwargs (dict) – Additional arguments, e.g., ‘dt_step’ for integration step size.

Returns:

(r_f, v_f). r_f : Final position vector (m). v_f : Final velocity vector (m/s).

Return type:

tuple[np.ndarray, np.ndarray]

opengnc.propagators.encke module

Encke’s Method Propagator.

class opengnc.propagators.encke.EnckePropagator(integrator: Integrator | None = None, mu: float = 398600441800000.0, rect_tol: float = 1e-06)[source]

Bases: Propagator

Encke’s Method Propagator.

Integrates the deviation from a reference Keplerian orbit. Suitable for orbit propagation with small disturbances.

Parameters:
  • integrator (Integrator, optional) – Numerical integrator for deviations. Defaults to RK4.

  • mu (float, optional) – Gravitational parameter (m^3/s^2). Default is Earth’s.

  • rect_tol (float, optional) – Rectification tolerance (|dr| / r_ref). Default 1e-6.

propagate(r_i: ndarray, v_i: ndarray, dt: float, perturbation_acc_fn: Callable[[float, ndarray, ndarray], ndarray] | None = None, **kwargs: Any) tuple[ndarray, ndarray][source]

Propagate state using Encke’s method.

Parameters:
  • r_i (np.ndarray) – Initial position vector (m).

  • v_i (np.ndarray) – Initial velocity vector (m/s).

  • dt (float) – Propagation duration (s).

  • perturbation_acc_fn (callable, optional) – Function that returns perturbation accelerations. Signature: acc_pert = f(t, r, v) -> np.ndarray.

  • **kwargs (dict) – Additional arguments, e.g., ‘dt_step’ for integration step size.

Returns:

(r_f, v_f). r_f : Final position vector (m). v_f : Final velocity vector (m/s).

Return type:

tuple[np.ndarray, np.ndarray]

opengnc.propagators.gve module

Gauss-Variational Equations (GVE) Propagator.

class opengnc.propagators.gve.GVEPropagator(mu: float = 398600441800000.0)[source]

Bases: Propagator

Propagator based on Gauss-Variational Equations.

Suitable for analyzing the effects of small perturbations on Keplerian elements.

Parameters:

mu (float) – Gravitational parameter (m^3/s^2).

gve_derivatives(t: float, state: ndarray, perturbation_func: Callable[[float, ndarray], ndarray]) ndarray[source]

Calculate GVE state derivatives.

Parameters:
  • t (float) – Time.

  • state (np.ndarray) – $[a, e, i, Omega, omega, theta]$.

  • perturbation_func (Callable) – RIC acceleration function.

Returns:

$[da, de, di, dOmega, domega, dtheta]$.

Return type:

np.ndarray

propagate(t0: float, state0: ndarray, dt: float, perturbation_func: Callable[[float, ndarray], ndarray]) ndarray[source]

Propagate orbital elements one step using RK4 on GVE.

Parameters:
  • t0 (float) – Initial time (s).

  • state0 (np.ndarray) – Initial Keplerian elements $[a, e, i, Omega, omega, theta]$.

  • dt (float) – Time step (s).

  • perturbation_func (Callable) – Function returning acceleration vector $[a_r, a_s, a_w]$ in RIC frame.

Returns:

New orbital elements.

Return type:

np.ndarray

opengnc.propagators.kepler module

Analytical Two-Body Propagator using Kepler’s Equation.

class opengnc.propagators.kepler.KeplerPropagator(mu: float = 398600441800000.0)[source]

Bases: Propagator

Analytical Two-Body Propagator.

Solves Kepler’s Problem using Universal Variables: $M = E - e sin E$

Parameters:

mu (float, optional) – Gravitational parameter ($m^3/s^2$).

propagate(r_i: ndarray, v_i: ndarray, dt: float, n_iter: int = 100, **kwargs: Any) tuple[ndarray, ndarray][source]

Solves the two-body problem using the Kepler equation (Vallado implementation).

Parameters:
  • r_i (np.ndarray) – Initial position vector (m).

  • v_i (np.ndarray) – Initial velocity vector (m/s).

  • dt (float) – Propagation duration (s).

  • n_iter (int, optional) – Maximum number of iterations for Newton-Raphson. Default is 100.

  • **kwargs (dict) – Additional arguments.

Returns:

(r_f, v_f). r_f : Final position vector (m). v_f : Final velocity vector (m/s).

Return type:

tuple[np.ndarray, np.ndarray]

opengnc.propagators.sgp4_propagator module

SGP4/SDP4 Analytical Propagator.

class opengnc.propagators.sgp4_propagator.Sgp4Propagator(line1: str, line2: str)[source]

Bases: Propagator

SGP4/SDP4 Analytical Propagator.

Wraps the sgp4 library to propagate orbits using Two-Line Elements (TLEs). State output is typically in the TEME frame.

Parameters:
  • line1 (str) – First line of the TLE.

  • line2 (str) – Second line of the TLE.

propagate(r_i: ndarray, v_i: ndarray, dt: float, **kwargs: Any) tuple[ndarray, ndarray][source]

Propagates the satellite state from TLE epoch forward by dt.

Note: initial position/velocity input (r_i, v_i) is IGNORED, as SGP4 is fully determined by the initialized TLE.

Parameters:
  • r_i (np.ndarray) – Ignored initial position.

  • v_i (np.ndarray) – Ignored initial velocity.

  • dt (float) – Propagation time relative to TLE epoch (s).

  • **kwargs (dict) – Additional arguments.

Returns:

(r_f, v_f). r_f : Position in TEME frame (m). v_f : Velocity in TEME frame (m/s).

Return type:

tuple[np.ndarray, np.ndarray]

propagate_to_jd(jd_f: float, jd_f_frac: float = 0.0) tuple[ndarray, ndarray][source]

Propagate to a specific Julian Date.

Parameters:
  • jd_f (float) – Final Julian Date (integer part).

  • jd_f_frac (float, optional) – Final Julian Date (fractional part). Default is 0.0.

Returns:

(r_f, v_f) in TEME frame (m, m/s).

Return type:

tuple[np.ndarray, np.ndarray]

Raises:

RuntimeError – If SGP4 propagation fails.

Module contents