Modules documentation

Filters:

Triple Moving Average

Here we take the average of 3 terms x0, A, B where, x0 = The point to be estimated A = weighted average of n terms previous to x0 B = weighted avreage of n terms ahead of x0 n = window size

orbitdeterminator.filters.triple_moving_average.generate_filtered_data(in_data, window)[source]

Apply the filter and generate the filtered data

Parameters:
  • in_data (string) – numpy array containing the positional data
  • window (int) – window size applied into the filter
Returns:

the final filtered array

Return type:

numpy array

orbitdeterminator.filters.triple_moving_average.triple_moving_average(signal_array, window_size)[source]

Apply triple moving average to a signal

Parameters:
  • signal_array (numpy array) – the array of values on which the filter is to be applied
  • window_size (int) – the no. of points before and after x0 which should be considered for calculating A and B
Returns:

a filtered array of size same as that of signal_array

Return type:

numpy array

orbitdeterminator.filters.triple_moving_average.weighted_average(params)[source]

Calculates the weighted average of terms in the input

Parameters:params (list) – a list of numbers
Returns:weighted average of the terms in the list
Return type:list

Savintzky - Golay

Takes a positional data set (time, x, y, z) and applies the Savintzky Golay filter on it based on the polynomial and window parameters we input

orbitdeterminator.filters.sav_golay.golay(data, window, degree)[source]

Apply the Savintzky-Golay filter to a positional data set.

Parameters:
  • data (numpy array) – containing all of the positional data in the format of (time, x, y, z)
  • window (int) – window size of the Savintzky-Golay filter
  • degree (int) – degree of the polynomial in Savintzky-Golay filter
Returns:

filtered data in the same format

Return type:

numpy array

Interpolation:

Lamberts-Kalman Method

Takes a positional data set and produces sets of six keplerian elements using Lambert’s solution for preliminary orbit determination and Kalman filters

orbitdeterminator.kep_determination.lamberts_kalman.check_keplerian(kep)[source]

Checks all the sets of keplerian elements to see if they have wrong values like eccentricity greater that 1 or a negative number for semi major axis

Parameters:kep (numpy array) – all the sets of keplerian elements in [semi major axis (a), eccentricity (e), inclination (i), argument of perigee (ω), right ascension of the ascending node (Ω), true anomaly (v)] format
Returns:the final corrected set of keplerian elements that will be inputed in the kalman filter
Return type:numpy array
orbitdeterminator.kep_determination.lamberts_kalman.create_kep(my_data)[source]

Computes all the keplerian elements for every point of the orbit you provide using Lambert’s solution It implements a tool for deleting all the points that give extremely jittery state vectors

Parameters:data (numpy array) – contains the positional data set in (Time, x, y, z) Format
Returns:array containing all the keplerian elements computed for the orbit given in [semi major axis (a), eccentricity (e), inclination (i), argument of perigee (ω), right ascension of the ascending node (Ω), true anomaly (v)] format
Return type:numpy array
orbitdeterminator.kep_determination.lamberts_kalman.kalman(kep, R)[source]

Takes as an input lots of sets of keplerian elements and produces the fitted value of them by applying kalman filters

Parameters:
  • kep (numpy array) – containing keplerian elements in this format (a, e, i, ω, Ω, v)
  • R – estimate of measurement variance
Returns:

final set of keplerian elements describing the orbit based on kalman filtering

Return type:

numpy array

orbitdeterminator.kep_determination.lamberts_kalman.orbit_trajectory(x1_new, x2_new, time)[source]

Tool for checking if the motion of the sallite is retrogade or counter - clock wise

Parameters:
  • x1 (numpy array) – time and position for point 1 [time1,x1,y1,z1]
  • x2 (numpy array) – time and position for point 2 [time2,x2,y2,z2]
  • time (float) – time difference between the 2 points
Returns:

true if we want to keep retrograde, False if we want counter-clock wise

Return type:

bool

Gibb’s Method

Spline Interpolation

Interpolation using splines for calculating velocity at a point and hence the orbital elements

orbitdeterminator.kep_determination.interpolation.compute_velocity(spline, point)[source]

Calculate the derivative of spline at the point(on the points the given spline corresponds to). This gives the velocity at that point.

Parameters:
  • spline (list) – component wise cubic splines of orbit data points of the format [spline_x, spline_y, spline_z].
  • point (numpy array) – point at which velocity is to be calculated.
Returns:

velocity vector at the given point

Return type:

numpy array

orbitdeterminator.kep_determination.interpolation.cubic_spline(orbit_data)[source]

Compute component wise cubic spline of points of input data

Parameters:
  • orbit_data (numpy array) – array of orbit data points of the
  • [time, x, y, z] (format) –
Returns:

component wise cubic splines of orbit data points of the format [spline_x, spline_y, spline_z]

Return type:

list

orbitdeterminator.kep_determination.interpolation.main(data_points)[source]

Apply the whole process of interpolation for keplerian element computation

Parameters:data_points (numpy array) – positional data set in format of (time, x, y, z)
Returns:computed keplerian elements for every point of the orbit
Return type:numpy array

Ellipse Fit

Finds out the ellipse that best fits to a set of data points and calculates its keplerian elements.

orbitdeterminator.kep_determination.ellipse_fit.determine_kep(data)[source]

Determines keplerian elements that fit a set of points.

Parameters:data (nx3 numpy array) – A numpy array of points in the format [x y z].
Returns:(kep,res) - The keplerian elements and the residuals as a tuple. kep: 1x6 numpy array res: nx3 numpy array

For the keplerian elements: kep[0] - semi-major axis (in whatever units the data was provided in) kep[1] - eccentricity kep[2] - inclination (in degrees) kep[3] - argument of periapsis (in degrees) kep[4] - right ascension of ascending node (in degrees) kep[5] - true anomaly of the first row in the data (in degrees)

For the residuals: (in whatever units the data was provided in) res[0] - residuals in x axis res[1] - residuals in y axis res[2] - residuals in z axis

orbitdeterminator.kep_determination.ellipse_fit.plot_kep(kep, data)[source]

Plots the original data and the orbit defined by the keplerian elements.

Parameters:
  • kep (1x6 numpy array) – keplerian elements
  • data (nx3 numpy array) – original data
Returns:

nothing

Gauss method

Implements Gauss’ method for three topocentric right ascension and declination measurements of celestial bodies. Supports both Earth-centered and Sun-centered orbits.

orbitdeterminator.kep_determination.gauss_method.alpha(x, y, z, u, v, w, mu)[source]

Compute the inverse of the semimajor axis.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
Returns:

alpha = 1/a

Return type:

float

orbitdeterminator.kep_determination.gauss_method.angle_diff_rad(a1, a2)[source]

Computes signed difference between two angles. Input angles are assumed to be in radians. Result is returned in radians. Code adapted from https://rosettacode.org/wiki/Angle_difference_between_two_bearings#Python.

Args:
a1 (float): angle 1 in radians a2 (float): angle 2 in radians
Returns:
r (float): shortest signed difference in radians
orbitdeterminator.kep_determination.gauss_method.argperi(x, y, z, u, v, w, mu)[source]

Compute the argument of pericenter.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
Returns:

argument of pericenter

Return type:

float

orbitdeterminator.kep_determination.gauss_method.earth_ephemeris(t_tdb)[source]

Compute heliocentric position of Earth at Julian date t_tdb (TDB, days), according to SPK kernel defined by astropy.coordinates.solar_system_ephemeris.

Args:
t_tdb (float): TDB instant of requested position
Returns:
(1x3 array): cartesian position in km
orbitdeterminator.kep_determination.gauss_method.eccentricity(x, y, z, u, v, w, mu)[source]

Compute value of eccentricity, e.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
Returns:

eccentricity, e

Return type:

float

orbitdeterminator.kep_determination.gauss_method.gauss_estimate_mpc(mpc_object_data, mpc_observatories_data, inds, r2_root_ind=0)[source]

Gauss method implementation for MPC Near-Earth asteroids ra/dec tracking data.

Parameters:
  • mpc_object_data (string) – path to MPC-formatted observation data file
  • mpc_observatories_data (string) – path to MPC observation sites data file
  • inds (1x3 int array) – indices of requested data
  • r2_root_ind (int) – index of selected Gauss polynomial root
Returns:

updated position at first observation r2 (1x3 array): updated position at second observation r3 (1x3 array): updated position at third observation v2 (1x3 array): updated velocity at second observation D (3x3 array): auxiliary matrix R (1x3 array): three observer position vectors rho1 (1x3 array): LOS vector at first observation rho2 (1x3 array): LOS vector at second observation rho3 (1x3 array): LOS vector at third observation tau1 (float): time interval from second to first observation tau3 (float): time interval from second to third observation f1 (float): Lagrange’s f function value at first observation g1 (float): Lagrange’s g function value at first observation f3 (float): Lagrange’s f function value at third observation g3 (float): Lagrange’s g function value at third observation Ea_hc_pos (1x3 array): cartesian position vectors (Earth wrt Sun) rho_1_sr (float): slant range at first observation rho_2_sr (float): slant range at second observation rho_3_sr (float): slant range at third observation obs_t (1x3 array): three times of observations

Return type:

r1 (1x3 array)

orbitdeterminator.kep_determination.gauss_method.gauss_estimate_sat(iod_object_data, sat_observatories_data, inds, r2_root_ind=0)[source]

Gauss method implementation for Earth-orbiting satellites ra/dec tracking data. Assumes observation data uses IOD format, with angle subformat 2.

Args:
iod_object_data (string): file path to sat tracking observation data of object sat_observatories_data (string): path to file containing COSPAR satellite tracking stations data. inds (1x3 int array): line numbers in data file to be processed r2_root_ind (int): index of selected Gauss polynomial root
Returns:
r1 (1x3 array): updated position at first observation r2 (1x3 array): updated position at second observation r3 (1x3 array): updated position at third observation v2 (1x3 array): updated velocity at second observation D (3x3 array): auxiliary matrix R (1x3 array): three observer position vectors rho1 (1x3 array): LOS vector at first observation rho2 (1x3 array): LOS vector at second observation rho3 (1x3 array): LOS vector at third observation tau1 (float): time interval from second to first observation tau3 (float): time interval from second to third observation f1 (float): Lagrange’s f function value at first observation g1 (float): Lagrange’s g function value at first observation f3 (float): Lagrange’s f function value at third observation g3 (float): Lagrange’s g function value at third observation rho_1_sr (float): slant range at first observation rho_2_sr (float): slant range at second observation rho_3_sr (float): slant range at third observation obs_t_jd (1x3 array): three Julian dates of observations
orbitdeterminator.kep_determination.gauss_method.gauss_iterator_mpc(mpc_object_data, mpc_observatories_data, inds, refiters=0, r2_root_ind=0)[source]

Gauss method iterator for minor planets ra/dec tracking data. Computes a first estimate of the orbit using gauss_estimate_sat function, and then refines this estimate using gauss_refinement. Assumes observation data file follows MPC format.

Args:
mpc_object_data (string): path to MPC-formatted observation data file mpc_observatories_data (string): path to MPC observation sites data file inds (1x3 int array): line numbers in data file to be processed refiters (int): number of refinement iterations to be performed r2_root_ind (int): index of selected Gauss polynomial root
Returns:
r1 (1x3 array): updated position at first observation r2 (1x3 array): updated position at second observation r3 (1x3 array): updated position at third observation v2 (1x3 array): updated velocity at second observation R (1x3 array): three observer position vectors rho1 (1x3 array): LOS vector at first observation rho2 (1x3 array): LOS vector at second observation rho3 (1x3 array): LOS vector at third observation rho_1_sr (float): slant range at first observation rho_2_sr (float): slant range at second observation rho_3_sr (float): slant range at third observation Ea_hc_pos (1x3 array): cartesian position vectors (Earth wrt Sun) obs_t (1x3 array): times of observations
orbitdeterminator.kep_determination.gauss_method.gauss_iterator_sat(iod_object_data, sat_observatories_data, inds, refiters=0, r2_root_ind=0)[source]

Gauss method iterator for Earth-orbiting satellites ra/dec tracking data. Computes a first estimate of the orbit using gauss_estimate_sat function, and then refines this estimate using gauss_refinement. Assumes observation data file is IOD-formatted, with angle subformat 2.

Args:
iod_object_data (string): file path to sat tracking observation data of object sat_observatories_data (string): path to file containing COSPAR satellite tracking stations data. inds (1x3 int array): line numbers in data file to be processed refiters (int): number of refinement iterations to be performed r2_root_ind (int): index of selected Gauss polynomial root
Returns:
r1 (1x3 array): updated position at first observation r2 (1x3 array): updated position at second observation r3 (1x3 array): updated position at third observation v2 (1x3 array): updated velocity at second observation R (1x3 array): three observer position vectors rho1 (1x3 array): LOS vector at first observation rho2 (1x3 array): LOS vector at second observation rho3 (1x3 array): LOS vector at third observation rho_1_sr (float): slant range at first observation rho_2_sr (float): slant range at second observation rho_3_sr (float): slant range at third observation obs_t (1x3 array): times of observations
orbitdeterminator.kep_determination.gauss_method.gauss_method_core(obs_radec, obs_t, R, mu, r2_root_ind=0)[source]

Perform core Gauss method.

Parameters:
  • obs_radec (1x3 SkyCoord array) – three rad/dec observations
  • obs_t (1x3 array) – three times of observations
  • R (1x3 array) – three observer position vectors
  • mu (float) – gravitational parameter of center of attraction
  • r2_root_ind (int) – index of Gauss polynomial root
Returns:

estimated position at first observation r2 (1x3 array): estimated position at second observation r3 (1x3 array): estimated position at third observation v2 (1x3 array): estimated velocity at second observation D (3x3 array): auxiliary matrix rho1 (1x3 array): LOS vector at first observation rho2 (1x3 array): LOS vector at second observation rho3 (1x3 array): LOS vector at third observation tau1 (float): time interval from second to first observation tau3 (float): time interval from second to third observation f1 (float): estimated Lagrange’s f function value at first observation g1 (float): estimated Lagrange’s g function value at first observation f3 (float): estimated Lagrange’s f function value at third observation g3 (float): estimated Lagrange’s g function value at third observation rho_1_sr (float): estimated slant range at first observation rho_2_sr (float): estimated slant range at second observation rho_3_sr (float): estimated slant range at third observation

Return type:

r1 (1x3 array)

orbitdeterminator.kep_determination.gauss_method.gauss_method_mpc(filename, bodyname, obs_arr=None, r2_root_ind_vec=None, refiters=0, plot=True)[source]

Gauss method high-level function for minor planets (asteroids, comets, etc.) orbit determination from MPC-formatted ra/dec tracking data. Roots of 8-th order Gauss polynomial are computed using np.roots function. Note that if r2_root_ind_vec is not specified by the user, then the first positive root returned by np.roots is used by default.

Args:
filename (string): path to MPC-formatted observation data file bodyname (string): user-defined name of minor planet obs_arr (int vector): line numbers in data file to be processed refiters (int): number of refinement iterations to be performed r2_root_ind_vec (1xlen(obs_arr) int array): indices of Gauss polynomial roots. plot (bool): if True, plots data.
Returns:
x (tuple): set of Keplerian orbital elements {(a, e, taup, omega, I, omega, T),t_vec[-1]}
orbitdeterminator.kep_determination.gauss_method.gauss_method_sat(filename, obs_arr=None, bodyname=None, r2_root_ind_vec=None, refiters=0, plot=True)[source]

Gauss method high-level function for orbit determination of Earth satellites from IOD-formatted ra/dec tracking data. IOD angle subformat 2 is assumed. Roots of 8-th order Gauss polynomial are computed using np.roots function. Note that if r2_root_ind_vec is not specified by the user, then the first positive root returned by np.roots is used by default.

Args:
filename (string): path to IOD-formatted observation data file obs_arr (int vector): line numbers in data file to be processed bodyname (string): user-defined name of satellite refiters (int): number of refinement iterations to be performed r2_root_ind_vec (1xlen(obs_arr) int array): indices of Gauss polynomial roots. plot (bool): if True, plots data.
Returns:
x (tuple): set of Keplerian orbital elements (a, e, taup, omega, I, omega, T)
orbitdeterminator.kep_determination.gauss_method.gauss_method_sat_passes(filename, obs_arr=None, bodyname=None, r2_root_ind_vec=None, refiters=10, plot=False)[source]

Gauss method high-level function for orbit determination of Earth satellites from IOD-formatted ra/dec tracking data. Roots of 8-th order Gauss polynomial are computed using np.roots function. Note that if r2_root_ind_vec is not specified by the user, then the first positive root returned by np.roots is used by default.

Args:
filename (string): path to IOD-formatted observation data file obs_arr (int vector): line numbers in data file to be processed bodyname (string): user-defined name of satellite refiters (int): number of refinement iterations to be performed r2_root_ind_vec (1xlen(obs_arr) int array): indices of Gauss polynomial roots. plot (bool): if True, plots data.
Returns:
x (tuple): set of Keplerian orbital elements {(a, e, taup, omega, I, omega, T),t_vec[-1]}
orbitdeterminator.kep_determination.gauss_method.gauss_refinement(mu, tau1, tau3, r2, v2, atol, D, R, rho1, rho2, rho3, f_1, g_1, f_3, g_3)[source]

Perform refinement of Gauss method.

Parameters:
  • mu (float) – gravitational parameter of center of attraction
  • tau1 (float) – time interval from second to first observation
  • tau3 (float) – time interval from second to third observation
  • r2 (1x3 array) – estimated position at second observation
  • v2 (1x3 array) – estimated velocity at second observation
  • atol (float) – absolute tolerance of universal Kepler anomaly computation
  • D (3x3 array) – auxiliary matrix
  • R (1x3 array) – three observer position vectors
  • rho1 (1x3 array) – LOS vector at first observation
  • rho2 (1x3 array) – LOS vector at second observation
  • rho3 (1x3 array) – LOS vector at third observation
  • f_1 (float) – estimated Lagrange’s f function value at first observation
  • g_1 (float) – estimated Lagrange’s g function value at first observation
  • f_3 (float) – estimated Lagrange’s f function value at third observation
  • g_3 (float) – estimated Lagrange’s g function value at third observation
Returns:

updated position at first observation r2 (1x3 array): updated position at second observation r3 (1x3 array): updated position at third observation v2 (1x3 array): updated velocity at second observation rho_1_sr (float): updated slant range at first observation rho_2_sr (float): updated slant range at second observation rho_3_sr (float): updated slant range at third observation f_1_new (float): updated Lagrange’s f function value at first observation g_1_new (float): updated Lagrange’s g function value at first observation f_3_new (float): updated Lagrange’s f function value at third observation g_3_new (float): updated Lagrange’s g function value at third observation

Return type:

r1 (1x3 array)

orbitdeterminator.kep_determination.gauss_method.get_observations_data(mpc_object_data, inds)[source]

Extract three ra/dec observations from MPC observation data file.

Parameters:
  • mpc_object_data (string) – file path to MPC observation data of object
  • inds (int array) – indices of requested data
Returns:

ra/dec observation data obs_t (1x3 Time array): time observation data site_codes (1x3 int array): corresponding codes of observation sites

Return type:

obs_radec (1x3 SkyCoord array)

orbitdeterminator.kep_determination.gauss_method.get_observations_data_sat(iod_object_data, inds)[source]

Extract three ra/dec observations from IOD observation data file.

Parameters:
  • iod_object_data (string) – file path to sat tracking observation data of object
  • inds (int array) – indices of requested data
Returns:

ra/dec observation data obs_t (1x3 Time array): time observation data site_codes (1x3 int array): corresponding codes of observation sites

Return type:

obs_radec (1x3 SkyCoord array)

orbitdeterminator.kep_determination.gauss_method.get_observatory_data(observatory_code, mpc_observatories_data)[source]

Load individual data of MPC observatory corresponding to given observatory code.

Parameters:
  • observatory_code (int) – MPC observatory code.
  • mpc_observatories_data (string) – path to file containing MPC observatories data.
Returns:

observatory data corresponding to code.

Return type:

ndarray

orbitdeterminator.kep_determination.gauss_method.get_observer_pos_wrt_earth(sat_observatories_data, obs_radec, site_codes)[source]

Compute position of observer at Earth’s surface, with respect to the Earth, in equatorial frame, during 3 distinct instants.

Args:
sat_observatories_data (string): path to file containing COSPAR satellite tracking stations data. obs_radec (1x3 SkyCoord array): three rad/dec observations site_codes (1x3 int array): COSPAR codes of observation sites
Returns:
R (1x3 array): cartesian position vectors (observer wrt Earth)
orbitdeterminator.kep_determination.gauss_method.get_observer_pos_wrt_sun(mpc_observatories_data, obs_radec, site_codes)[source]

Compute position of observer at Earth’s surface, with respect to the Sun, in equatorial frame, during 3 distinct instants.

Args:
mpc_observatories_data (string): path to file containing MPC observatories data. obs_radec (1x3 SkyCoord array): three rad/dec observations site_codes (1x3 int array): MPC codes of observation sites
Returns:
R (1x3 array): cartesian position vectors (observer wrt Sun) Ea_hc_pos (1x3 array): cartesian position vectors (Earth wrt Sun)
orbitdeterminator.kep_determination.gauss_method.get_time_of_observation(year, month, day, hour, minute, second, msecond)[source]

creates time variable

Parameters:
  • year
  • month
  • day
  • hour
  • minute
  • second
  • msecond
Returns:

orbitdeterminator.kep_determination.gauss_method.inclination(x, y, z, u, v, w)[source]

Compute value of inclination, I.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
Returns:

inclination, I

Return type:

float

orbitdeterminator.kep_determination.gauss_method.kep_h_norm(x, y, z, u, v, w)[source]

Compute norm of specific angular momentum vector h.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
Returns:

norm of specific angular momentum vector, h.

Return type:

float

orbitdeterminator.kep_determination.gauss_method.kep_h_vec(x, y, z, u, v, w)[source]

Compute specific angular momentum vector h.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
Returns:

specific angular momentum vector, h.

Return type:

float

orbitdeterminator.kep_determination.gauss_method.lagrangef(mu, r2, tau)[source]

Compute 1st order approximation to Lagrange’s f function.

Parameters:
  • mu (float) – gravitational parameter attracting body
  • r2 (float) – radial distance
  • tau (float) – time interval
Returns:

Lagrange’s f function value

Return type:

float

orbitdeterminator.kep_determination.gauss_method.lagrangef_(xi, z, r)[source]

Compute current value of Lagrange’s f function.

Parameters:
  • xi (float) – universal Kepler anomaly
  • z (float) – xi**2/alpha
  • r (float) – radial distance
Returns:

Lagrange’s f function value

Return type:

float

orbitdeterminator.kep_determination.gauss_method.lagrangeg(mu, r2, tau)[source]

Compute 1st order approximation to Lagrange’s g function.

Parameters:
  • mu (float) – gravitational parameter attracting body
  • r2 (float) – radial distance
  • tau (float) – time interval
Returns:

Lagrange’s g function value

Return type:

float

orbitdeterminator.kep_determination.gauss_method.lagrangeg_(tau, xi, z, mu)[source]

Compute current value of Lagrange’s g function.

Parameters:
  • tau (float) – time interval
  • xi (float) – universal Kepler anomaly
  • z (float) – xi**2/alpha
  • r (float) – radial distance
Returns:

Lagrange’s g function value

Return type:

float

orbitdeterminator.kep_determination.gauss_method.load_mpc_data(fname)[source]

Loads minor planet position observation data from MPC-formatted files. MPC format for minor planet observations is described at https://www.minorplanetcenter.net/iau/info/OpticalObs.html TODO: Add support for comets and natural satellites. Add support for radar observations: https://www.minorplanetcenter.net/iau/info/RadarObs.html See also NOTE 2 in: https://www.minorplanetcenter.net/iau/info/OpticalObs.html

Parameters:fname (string) – name of the MPC-formatted text file to be parsed
Returns:array of minor planet position observations following the MPC format.
Return type:x (ndarray)
orbitdeterminator.kep_determination.gauss_method.load_mpc_observatories_data(mpc_observatories_fname)[source]

Load Minor Planet Center observatories data using numpy’s genfromtxt function.

Parameters:mpc_observatories_fname (str) – file name with MPC observatories data.
Returns:data read from the text file (output from numpy.genfromtxt)
Return type:ndarray
orbitdeterminator.kep_determination.gauss_method.longascnode(x, y, z, u, v, w)[source]

Compute value of longitude of ascending node, computed as the angle between x-axis and the vector n = (-hy,hx,0), where hx, hy, are respectively, the x and y components of specific angular momentum vector, h.

Args:
x (float): x-component of position y (float): y-component of position z (float): z-component of position u (float): x-component of velocity v (float): y-component of velocity w (float): z-component of velocity
Returns:
float: longitude of ascending node
orbitdeterminator.kep_determination.gauss_method.losvector(ra_rad, dec_rad)[source]

Compute line-of-sight (LOS) vector for given values of right ascension and declination. Both angles must be provided in radians.

Args:
ra_rad (float): right ascension (rad) dec_rad (float): declination (rad)
Returns:
1x3 numpy array: cartesian components of LOS vector.
orbitdeterminator.kep_determination.gauss_method.object_wrt_sun(t_utc, a, e, taup, omega, I, Omega)[source]

Compute position of celestial object with respect to the Sun, in equatorial frame.

Parameters:
  • t_utc (Time) – UTC time of observation
  • a (float) – semimajor axis
  • e (float) – eccentricity
  • taup (float) – time of pericenter passage
  • omega (float) – argument of pericenter
  • I (float) – inclination
  • Omega (float) – longitude of ascending node
Returns:

cartesian vector

Return type:

(1x3 array)

orbitdeterminator.kep_determination.gauss_method.observer_wrt_sun(long, parallax_s, parallax_c, t_utc)[source]

Compute position of observer at Earth’s surface, with respect to the Sun, in equatorial frame.

Args:
long (float): longitude of observing site parallax_s (float): parallax constant S of observing site parallax_c (float): parallax constant C of observing site t_utc (Time): UTC time of observation
Returns:
(1x3 array): cartesian vector
orbitdeterminator.kep_determination.gauss_method.observerpos_mpc(long, parallax_s, parallax_c, t_utc)[source]

Compute geocentric observer position at UTC instant t_utc, for Sun-centered orbits, at a given observation site defined by its longitude, and parallax constants S and C. Formula taken from top of page 266, chapter 5, Orbital Mechanics book (Curtis). The parallax constants S and C are defined by: S=rho cos phi’ C=rho sin phi’, where rho: slant range phi’: geocentric latitude

Args:
long (float): longitude of observing site parallax_s (float): parallax constant S of observing site parallax_c (float): parallax constant C of observing site t_utc (astropy.time.Time): UTC time of observation
Returns:
1x3 numpy array: cartesian components of observer’s geocentric position
orbitdeterminator.kep_determination.gauss_method.observerpos_sat(lat, long, elev, t_utc)[source]

Compute geocentric observer position at UTC instant t_utc, for Earth-centered orbits, at a given observation site defined by its longitude, geodetic latitude and elevation above reference ellipsoid. Formula taken from bottom of page 265 (Eq. 5.56), chapter 5, Orbital Mechanics book (Curtis).

Args:
lat (float): geodetic latitude (deg) long (float): longitude (deg) elev (float): elevation above reference ellipsoid (m) t_utc (astropy.time.Time): UTC time of observation
Returns:
1x3 numpy array: cartesian components of observer’s geocentric position
orbitdeterminator.kep_determination.gauss_method.radec_obs_vec_mpc(inds, mpc_object_data)[source]

Compute vector of observed ra,dec values for MPC tracking data.

Parameters:
  • inds (int array) – line numbers of data in file
  • mpc_object_data (ndarray) – MPC observation data for object
Returns:

vector of ra/dec observed values

Return type:

rov (1xlen(inds) array)

orbitdeterminator.kep_determination.gauss_method.radec_obs_vec_sat(inds, iod_object_data)[source]

Compute vector of observed ra,dec values for satellite tracking data (IOD-formatted).

Parameters:
  • inds (int array) – line numbers of data in file
  • iod_object_data (ndarray) – observation data
Returns:

vector of ra/dec observed values

Return type:

rov (1xlen(inds) array)

orbitdeterminator.kep_determination.gauss_method.radec_res_vec_rov_mpc(x, inds, mpc_object_data, mpc_observatories_data, rov)[source]

Compute vector of observed minus computed (O-C) residuals for ra/dec MPC-formatted observations of minor planets (asteroids, comets, etc.), with pre-computed observed radec values vector. Assumes ra/dec observed values vector is contained in rov, and they are stored as rov = [ra1, dec1, ra2, dec2, …].

Args:
x (1x6 float array): set of orbital elements (a, e, taup, omega, I, Omega) inds (int array): line numbers of data in file mpc_object_data (ndarray): observation data mpc_observatories_data (ndarray): MPC observatories data rov (1xlen(inds) float-like array): vector of observed ra/dec values
Returns:
rv (1xlen(inds) array): vector of ra/dec (O-C) residuals.
orbitdeterminator.kep_determination.gauss_method.radec_res_vec_rov_sat(x, inds, iod_object_data, sat_observatories_data, rov)[source]

Compute vector of observed minus computed (O-C) residuals for ra/dec Earth-orbiting satellite observations with pre-computed observed radec values vector. Assumes ra/dec observed values vector is contained in rov, and they are stored as rov = [ra1, dec1, ra2, dec2, …].

Args:
x (1x6 float array): set of orbital elements (a, e, taup, omega, I, Omega) inds (int array): line numbers of data in file iod_object_data (ndarray): observation data sat_observatories_data (ndarray): satellite tracking stations data rov (1xlen(inds) float-like array): vector of observed ra/dec values
Returns:
rv (1xlen(inds) array): vector of ra/dec (O-C) residuals.
orbitdeterminator.kep_determination.gauss_method.radec_residual_mpc(x, t_ra_dec_datapoint, long, parallax_s, parallax_c)[source]

Compute observed minus computed (O-C) residual for a given ra/dec datapoint, represented as a SkyCoord object, for MPC observation data.

Args:
x (1x6 array): set of Keplerian elements t_ra_dec_datapoint (SkyCoord): ra/dec datapoint long (float): longitude of observing site parallax_s (float): parallax constant S of observing site parallax_c (float): parallax constant C of observing site
Returns:
(1x2 array): right ascension difference, declination difference
orbitdeterminator.kep_determination.gauss_method.radec_residual_rov_mpc(x, t, ra_obs_rad, dec_obs_rad, long, parallax_s, parallax_c)[source]

Compute right ascension and declination observed minus computed (O-C) residual, using precomputed vector of observed ra/dec values, for MPC observation data.

Args:
x (1x6 array): set of Keplerian elements t (Time): time of observation ra_obs_rad (float): observed right ascension (rad) dec_obs_rad (float): observed declination (rad) long (float): longitude of observing site parallax_s (float): parallax constant S of observing site parallax_c (float): parallax constant C of observing site
Returns:
(1x2 array): right ascension difference, declination difference
orbitdeterminator.kep_determination.gauss_method.rho_vec(long, parallax_s, parallax_c, t_utc, a, e, taup, omega, I, Omega)[source]

Compute slant range vector.

Parameters:
  • long (float) – longitude of observing site
  • parallax_s (float) – parallax constant S of observing site
  • parallax_c (float) – parallax constant C of observing site
  • t_utc (Time) – UTC time of observation
  • a (float) – semimajor axis
  • e (float) – eccentricity
  • taup (float) – time of pericenter passage
  • omega (float) – argument of pericenter
  • I (float) – inclination
  • Omega (float) – longitude of ascending node
Returns:

cartesian vector

Return type:

(1x3 array)

orbitdeterminator.kep_determination.gauss_method.rhovec2radec(long, parallax_s, parallax_c, t_utc, a, e, taup, omega, I, Omega)[source]

Transform slant range vector to ra/dec values.

Parameters:
  • long (float) – longitude of observing site
  • parallax_s (float) – parallax constant S of observing site
  • parallax_c (float) – parallax constant C of observing site
  • t_utc (Time) – UTC time of observation
  • a (float) – semimajor axis
  • e (float) – eccentricity
  • taup (float) – time of pericenter passage
  • omega (float) – argument of pericenter
  • I (float) – inclination
  • Omega (float) – longitude of ascending node
Returns:

right ascension (rad) dec_rad (float): declination (rad)

Return type:

ra_rad (float)

orbitdeterminator.kep_determination.gauss_method.rungelenz(x, y, z, u, v, w, mu)[source]

Compute the cartesian components of Laplace-Runge-Lenz vector.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
Returns:

Laplace-Runge-Lenz vector

Return type:

float

orbitdeterminator.kep_determination.gauss_method.semimajoraxis(x, y, z, u, v, w, mu)[source]

Compute value of semimajor axis, a.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
Returns:

semimajor axis, a

Return type:

float

orbitdeterminator.kep_determination.gauss_method.t_radec_res_vec_mpc(x, inds, mpc_object_data, mpc_observatories_data)[source]

Compute vector of observed minus computed (O-C) residuals for ra/dec MPC-formatted observations of minor planets (asteroids, comets, etc.), with pre-computed observed radec values vector. Assumes ra/dec observed values vector is contained in rov, and they are stored as rov = [ra1, dec1, ra2, dec2, …].

Args:
x (1x6 float array): set of orbital elements (a, e, taup, omega, I, Omega) inds (int array): line numbers of data in file mpc_object_data (ndarray): observation data mpc_observatories_data (ndarray): MPC observatories data rov (1xlen(inds) float-like array): vector of observed ra/dec values
Returns:
rv (1xlen(inds) array): vector of ra/dec (O-C) residuals. tv (1xlen(inds) array): vector of observation times.
orbitdeterminator.kep_determination.gauss_method.t_radec_res_vec_sat(x, inds, iod_object_data, sat_observatories_data, rov)[source]

Compute vector of observed minus computed (O-C) residuals for ra/dec Earth-orbiting satellite observations with pre-computed observed radec values vector. Assumes ra/dec observed values vector is contained in rov, and they are stored as rov = [ra1, dec1, ra2, dec2, …].

Args:
x (1x6 float array): set of orbital elements (a, e, taup, omega, I, Omega) inds (int array): line numbers of data in file iod_object_data (ndarray): observation data sat_observatories_data (ndarray): satellite tracking stations data rov (1xlen(inds) float-like array): vector of observed ra/dec values
Returns:
rv (1xlen(inds) array): vector of ra/dec (O-C) residuals. tv (1xlen(inds) array): vector of observation times.
orbitdeterminator.kep_determination.gauss_method.taupericenter(t, e, f, n)[source]

Compute the time of pericenter passage.

Parameters:
  • t (float) – current time
  • e (float) – eccentricity
  • f (float) – true anomaly
  • n (float) – Keplerian mean motion
Returns:

time of pericenter passage

Return type:

float

orbitdeterminator.kep_determination.gauss_method.trueanomaly5(x, y, z, u, v, w, mu)[source]

Compute the true anomaly from cartesian state.

Parameters:
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
Returns:

true anomaly

Return type:

float

orbitdeterminator.kep_determination.gauss_method.univkepler(dt, x, y, z, u, v, w, mu, iters=5, atol=1e-15)[source]

Compute the current value of the universal Kepler anomaly, xi.

Parameters:
  • dt (float) – time interval
  • x (float) – x-component of position
  • y (float) – y-component of position
  • z (float) – z-component of position
  • u (float) – x-component of velocity
  • v (float) – y-component of velocity
  • w (float) – z-component of velocity
  • mu (float) – gravitational parameter
  • iters (int) – number of iterations of Newton-Raphson process
  • atol (float) – absolute tolerance of Newton-Raphson process
Returns:

alpha = 1/a

Return type:

float

Least squares

Computes the least-squares optimal Keplerian elements for a sequence of cartesian position observations.

orbitdeterminator.kep_determination.least_squares.gauss_LS_mpc(filename, bodyname, obs_arr, r2_root_ind_vec=None, obs_arr_ls=None, gaussiters=0, plot=True)[source]

Minor planets orbit determination high-level function from MPC-formatted ra/dec tracking data. Preliminary orbit determination via Gauss method is performed. Roots of 8-th order Gauss polynomial are computed using np.roots function. Note that if r2_root_ind_vec is not specified by the user, then the first positive root returned by np.roots is used by default.

Args:
filename (string): path to MPC-formatted observation data file bodyname (string): user-defined name of minor planet obs_arr (int vector): line numbers in data file to be processed in Gauss preliminary orbit determination r2_root_ind_vec (1xlen(obs_arr) int array): indices of Gauss polynomial roots. obs_arr (int vector): line numbers in data file to be processed in least-squares fit gaussiters (int): number of refinement iterations to be performed plot (bool): if True, plots data.
Returns:
x (tuple): set of Keplerian orbital elements (a, e, taup, omega, I, omega, T)
orbitdeterminator.kep_determination.least_squares.gauss_LS_sat(filename, bodyname, obs_arr, r2_root_ind_vec=None, obs_arr_ls=None, gaussiters=0, plot=True)[source]

Earth satellites orbit determination high-level function from IOD-formatted ra/dec tracking data. IOD angle subformat 2 is assumed. Preliminary orbit determination via Gauss method is performed. Roots of 8-th order Gauss polynomial are computed using np.roots function. Note that if r2_root_ind_vec is not specified by the user, then the first positive root returned by np.roots is used by default.

Args:
filename (string): path to IOD-formatted observation data file bodyname (string): user-defined name of satellite obs_arr (int vector): line numbers in data file to be processed in Gauss preliminary orbit determination r2_root_ind_vec (1xlen(obs_arr) int array): indices of Gauss polynomial roots. obs_arr (int vector): line numbers in data file to be processed in least-squares fit gaussiters (int): number of refinement iterations to be performed plot (bool): if True, plots data.
Returns:
x (tuple): set of Keplerian orbital elements (a, e, taup, omega, I, omega, T)
orbitdeterminator.kep_determination.least_squares.get_weights(resid)[source]

This function calculates the weights per (x,y,z) by using the inverse of the squared residuals divided by the total sum of the inverse of the squared residuals.

orbitdeterminator.kep_determination.least_squares.radec_res_vec_rov_mpc_w(x, inds, mpc_object_data, mpc_observatories_data, rov, weights)[source]

Compute vector of observed minus computed weighted (O-C) residuals for ra/dec MPC-formatted observations of minor planets (asteroids, comets, etc.), with pre-computed observed radec values vector. Assumes ra/dec observed values vector is contained in rov, and they are stored as rov = [ra1, dec1, ra2, dec2, …].

Args:
x (1x6 float array): set of orbital elements (a, e, taup, omega, I, Omega) inds (int array): line numbers of data in file mpc_object_data (ndarray): observation data mpc_observatories_data (ndarray): MPC observatories data rov (1xlen(inds) float-like array): vector of observed ra/dec values
Returns:
rv (1xlen(inds) array): vector of ra/dec (O-C) residuals.
orbitdeterminator.kep_determination.least_squares.radec_res_vec_rov_sat_w(x, inds, iod_object_data, sat_observatories_data, rov, weights)[source]

Compute vector of observed minus computed (O-C) weighted residuals for ra/dec Earth-orbiting satellite observations with pre-computed observed radec values vector. Assumes ra/dec observed values vector is contained in rov, and they are stored as rov = [ra1, dec1, ra2, dec2, …].

Args:
x (1x6 float array): set of orbital elements (a, e, taup, omega, I, Omega) inds (int array): line numbers of data in file iod_object_data (ndarray): observation data sat_observatories_data (ndarray): satellite tracking stations data rov (1xlen(inds) float-like array): vector of observed ra/dec values
Returns:
rv (1xlen(inds) array): vector of ra/dec weighted (O-C) residuals.

Propagation:

Propagation Model

class orbitdeterminator.propagation.sgp4.SGP4[source]
__init__()[source]

Initializes flag variable to check for FlagCheckError (custom exception).

compute_necessary_kep(kep, b_star=2.1109e-05)[source]

Initializes the necessary class variables using keplerian elements which are needed in the computation of the propagation model.

Parameters:
  • kep (list) – kep elements in order [axis, inclination, ascension, eccentricity, perigee, anomaly]
  • b_star (float) – bstar drag term
Returns:

NIL

compute_necessary_tle(line1, line2)[source]

Initializes the necessary class variables using TLE which are needed in the computation of the propagation model.

Parameters:
  • line1 (str) – line 1 of the TLE
  • line2 (str) – line 2 of the TLE
Returns:

NIL

propagate(t1, t2)[source]

Invokes the function to compute state vectors and organises the final result.

The function first checks if compute_necessary_xxx() is called or not if not then a custom exception is raised stating that call this function first. Then it computes the state vector for the next 8 hours (28800 seconds in 8 hours) at every time epoch (28800 time epcohs) using the sgp4 propagation model. The values of state vector is formatted upto five decimal points and then all the state vectors got appended in a list which stores the final output.

Parameters:
  • t1 (int) – start time epoch
  • t2 (int) – end time epoch
Returns:

vector containing all state vectors

Return type:

numpy.ndarray

propagation_model(tsince)[source]

From the time epoch and information from TLE, applies SGP4 on it.

The function applies the Simplified General Perturbations algorithm SGP4 on the information extracted from the TLE at the given time epoch ‘tsince’ and computes the state vector from it.

Parameters:tsince (int) – time epoch
Returns:position and velocity vector
Return type:tuple
classmethod recover_tle(pos, vel)[source]

Recovers TLE back from state vector.

First of all, only necessary information (which are inclination, right ascension of the ascending node, eccentricity, argument of perigee, mean anomaly, mean motion and bstar) that are needed in the computation of SGP4 propagation model are recovered. It is using a general format of TLE. State vectors are used to find orbital elements which are then inserted into the TLE format at their respective positions. Mean motion and bstar is calculated separately as it is not a part of orbital elements. Format of TLE: x denotes that there is a digit, c denotes a character value, underscore(_) denotes a plus/minus(+/-) sign value and period(.) denotes a decimal point.

Parameters:
  • pos (list) – position vector
  • vel (list) – velocity vector
Returns:

line1 and line2 of TLE

Return type:

list

class orbitdeterminator.propagation.sgp4.FlagCheckError[source]

Raised when compute_necessary_xxx() function is not called.

Cowell Method

Numerical orbit propagator based on RK4. Takes into account J2 and drag perturbations.

orbitdeterminator.propagation.cowell.drag(s)[source]

Returns the drag acceleration for a given state.

Parameters:s (1x6 numpy array) – the state vector [rx,ry,rz,vx,vy,vz]
Returns:the drag acceleration [ax,ay,az]
Return type:1x3 numpy array
orbitdeterminator.propagation.cowell.j2_pert(s)[source]

Returns the J2 acceleration for a given state.

Parameters:s (1x6 numpy array) – the state vector [rx,ry,rz,vx,vy,vz]
Returns:the J2 acceleration [ax,ay,az]
Return type:1x3 numpy array
orbitdeterminator.propagation.cowell.propagate_state(s, t0, tf)[source]

Equivalent to the rk4 function.

orbitdeterminator.propagation.cowell.rk4(s, t0, tf, h=30)[source]

Runge-Kutta 4th Order Numerical Integrator

Args:
s(1x6 numpy array): the state vector [rx,ry,rz,vx,vy,vz] t0(float) : initial time tf(float) : final time h(float) : step-size
Returns:the state at time tf
Return type:1x6 numpy array
orbitdeterminator.propagation.cowell.rkf45(s, t0, tf, h=10, tol=1e-06)[source]

Runge-Kutta Fehlberg 4(5) Numerical Integrator

Args:
s(1x6 numpy array): the state vector [rx,ry,rz,vx,vy,vz] t0(float) : initial time tf(float) : final time h(float) : step-size tol(float) : tolerance of error
Returns:the state at time tf
Return type:1x6 numpy array
orbitdeterminator.propagation.cowell.sdot(s)[source]

Returns the time derivative of a given state.

Parameters:s (1x6 numpy array) – the state vector [rx,ry,rz,vx,vy,vz]
Returns:the time derivative of s [vx,vy,vz,ax,ay,az]
Return type:1x6 numpy array
orbitdeterminator.propagation.cowell.time_period(s, h=30)[source]

Returns the nodal time period of an orbit.

Parameters:
  • s (1x6 numpy array) – the state vector [rx,ry,rz,vx,vy,vz]
  • h (float) – step-size
Returns:

the nodal time period of the orbit

Return type:

float

Simulator

class orbitdeterminator.propagation.simulator.Simulator(params)[source]

A class for the simulator.

__init__(params)[source]

Initializes the simulator.

Parameters:params – A SimParams object containing kep,t0,t,period,speed, and op_writer
Returns:nothing
calc()[source]

Calculates the satellite state at current time and calls itself after a certain amount of time.

simulate()[source]

Starts the calculation thread and waits for keyboard input. Press q or Ctrl-C to quit the simulator cleanly.

stop()[source]

Stops the simulator cleanly.

class orbitdeterminator.propagation.simulator.SimParams[source]

SimParams class. This is just a container for all the parameters required to start the simulation.

kep(1x6 numpy array): the intial osculating keplerian elements epoch(float): the epoch of the above kep period(float): maximum time period between observations t0(float): starting time of the simulation speed(float): speed of the simulation op_writer(OpWriter): output handling object

class orbitdeterminator.propagation.simulator.OpWriter[source]

Base output writer class. Inherit this class and override the methods.

close()[source]

Anything that has to be executed after finishing writing the output. Runs once.

Example: Closing connection to a database

open()[source]

Anything that has to be executed before starting to write output. Runs once.

Example: Establishing connection to database

static write(t, s)[source]

This method is called everytime the calc thread finishes a computation.

Parameters:
  • t – the current time of simulation
  • s – the state vector at t [rx,ry,rz,vx,vy,vz]
class orbitdeterminator.propagation.simulator.print_r[source]

Bases: orbitdeterminator.propagation.simulator.OpWriter

Prints the position vector

class orbitdeterminator.propagation.simulator.save_r(name)[source]

Bases: orbitdeterminator.propagation.simulator.OpWriter

Saves the position vector to a file

__init__(name)[source]

Initialize the class.

Parameters:name (string) – file name

DGSN Simulator

Kalman Filter

sgp4_prop

SGP4 propagator. This is a wrapper around the PyPI SGP4 propagator. However, this does not generate an artificial TLE. So there is no string manipulation involved. Hence this is faster than sgp4_prop_string.

orbitdeterminator.propagation.sgp4_prop.kep_to_sat(kep, epoch, bstar=0.21109E-4, whichconst=wgs72, afspc_mode=False)[source]

Converts a set of keplerian elements into a Satellite object.

Args:
kep(1x6 numpy array): the osculating keplerian elements at epoch epoch(float): the epoch bstar(float): bstar drag coefficient whichconst(float): gravity model. refer pypi sgp4 documentation afspc_mode(boolean): refer pypi sgp4 documentation
Returns:an sgp4 satellite object encapsulating the arguments
Return type:Satellite object
orbitdeterminator.propagation.sgp4_prop.propagate_kep(kep, t0, tf, bstar=2.1109e-05)[source]

Propagates a set of keplerian elements.

Parameters:
  • kep (1x6 numpy array) – osculating keplerian elements at epoch
  • t0 (float) – initial time (epoch)
  • tf (float) – final time
Returns:

the position at tf vel(1x3 numpy array): the velocity at tf

Return type:

pos(1x3 numpy array)

orbitdeterminator.propagation.sgp4_prop.propagate_state(r, v, t0, tf, bstar=2.1109e-05)[source]

Propagates a state vector

Parameters:
  • r (1x3 numpy array) – the position vector at epoch
  • v (1x3 numpy array) – the velocity vector at epoch
  • t0 (float) – initial time (epoch)
  • tf (float) – final time
Returns:

the position at tf vel(1x3 numpy array): the velocity at tf

Return type:

pos(1x3 numpy array)

sgp4_prop_string

SGP4 propagator. This is a wrapper around PyPI SGP4 propagator. It constructs an artificial TLE and passes it to the PyPI module.

orbitdeterminator.propagation.sgp4_prop_string.propagate(kep, init_time, final_time, bstar=2.1109e-05)[source]

Propagates a set of keplerian elements.

Parameters:
  • kep (1x6 numpy array) – osculating keplerian elements at epoch
  • init_time (float) – initial time (epoch)
  • final_time (float) – final time
  • bstar (float) – bstar drag coefficient
Returns:

the position at tf vel(1x3 numpy array): the velocity at tf

Return type:

pos(1x3 numpy array)

Utils:

kep_state

Takes a set of keplerian elements (a, e, i, ω, Ω, v) and transforms it into a state vector (x, y, z, vx, vy, vz) where v is the velocity of the satellite

orbitdeterminator.util.kep_state.kep_state(kep)[source]

Converts the keplerian elements to position and velocity vector

Parameters:kep (numpy array) – a 1x6 matrix which contains the following variables kep(0): semi major axis (km) kep(1): eccentricity (number) kep(2): inclination (degrees) kep(3): argument of perigee (degrees) kep(4): right ascension of the ascending node (degrees) kep(5): true anomaly (degrees)
Returns:1x6 matrix which contains the position and velocity vector r(0),r(1),r(2): position vector (x,y,z) km r(3),r(4),r(5): velocity vector (vx,vy,vz) km/s
Return type:numpy array

read_data

Reads the positional data set from a .csv file

orbitdeterminator.util.read_data.load_data(filename)[source]

Loads the data in numpy array for further processing in tab delimiter format

Parameters:filename (string) – name of the csv file to be parsed
Returns:array of the orbit positions, each point of the orbit is of the format (time, x, y, z)
Return type:numpy array
orbitdeterminator.util.read_data.save_orbits(source, destination)[source]

Saves objects returned from load_data

Parameters:
  • source – path to raw csv files.
  • destination – path where objects need to be saved.

state_kep

Takes a state vector (x, y, z, vx, vy, vz) where v is the velocity of the satellite and transforms it into a set of keplerian elements (a, e, i, ω, Ω, v)

orbitdeterminator.util.state_kep.state_kep(r, v)[source]

Converts state vector to orbital elements.

Parameters:
  • r (numpy array) – position vector
  • v (numpy array) – velocity vector
Returns:

array of the computed keplerian elements kep(0): semimajor axis (kilometers) kep(1): orbital eccentricity (non-dimensional)

(0 <= eccentricity < 1)

kep(2): orbital inclination (degrees) kep(3): argument of perigee (degress) kep(4): right ascension of ascending node (degrees) kep(5): true anomaly (degrees)

Return type:

numpy array

input_transf

Converts cartesian co-ordinates to spherical co-ordinates and vice versa

orbitdeterminator.util.input_transf.cart_to_spher(data)[source]

Takes as an input a data set containing points in cartesian format (time, x, y, z) and returns the computed spherical coordinates (time, azimuth, elevation, r)

Parameters:data (numpy array) – containing the cartesian coordinates in format of (time, x, y, z)
Returns:array of spherical coordinates in format of (time, azimuth, elevation, r)
Return type:numpy array
orbitdeterminator.util.input_transf.spher_to_cart(data)[source]

Takes as an input a data set containing points in spherical format (time, azimuth, elevation, r) and returns the computed cartesian coordinates (time, x, y, z).

Parameters:data (numpy array) – containing the spherical coordinates in format of (time, azimuth, elevation, r)
Returns:array of cartesian coordinates in format of (time, x, y, z)
Return type:numpy array

rkf78

Uses Runge Kutta Fehlberg 7(8) numerical integration method to compute the state vector in a time interval tf

orbitdeterminator.util.rkf78.rkf78(neq, ti, tf, h, tetol, x)[source]

Runge-Kutta-Fehlberg 7[8] method, solve first order system of differential equations

Parameters:
  • neq (int) – number of differential equations
  • ti (float) – initial simulation time
  • tf (float) – final simulation time
  • h (float) – initial guess for integration step size
  • tetol (float) – truncation error tolerance [non-dimensional]
  • x (numpy array) – integration vector at time = ti
Returns:

array of state vector at time tf

Return type:

numpy array

orbitdeterminator.util.rkf78.ypol_a(y)[source]

Computes velocity and acceleration values by using the state vector y and keplerian motion

Parameters:y (numpy array) – state vector (position + velocity)
Returns:derivative of the state vector (velocity + acceleration)
Return type:numpy array

golay_window

orbitdeterminator.util.golay_window.window(error, data)[source]

Calculates the constant c which is needed to determine the savintzky - golay filter window window = len(data) / c ,where c is a constant strongly related to the error contained in the data set

Parameters:
  • error (float) – the a-priori error estimation for each measurment
  • data (numpy array) – the positional data set
Returns:

constant which describes the window that needs to be inputed to the savintzky - golay filter

Return type:

float

anom_conv

Vectorized anomaly conversion scripts

orbitdeterminator.util.anom_conv.ecc_to_mean(E, e)[source]

Converts eccentric anomaly to mean anomaly.

Parameters:
  • E (numpy array) – array of eccentric anomalies (in radians)
  • e (float) – eccentricity
Returns:

array of mean anomalies (in radians)

Return type:

numpy array

orbitdeterminator.util.anom_conv.mean_to_t(M, a)[source]

Converts mean anomaly to time elapsed.

Parameters:
  • M (numpy array) – array of mean anomalies (in radians)
  • a (float) – semi-major axis
Returns:

numpy array of time elapsed

Return type:

numpy array

orbitdeterminator.util.anom_conv.true_to_ecc(theta, e)[source]

Converts true anomaly to eccentric anomaly.

Parameters:
  • theta (numpy array) – array of true anomalies (in radians)
  • e (float) – eccentricity
Returns:

array of eccentric anomalies (in radians)

Return type:

numpy array

new_tle_kep_state

This module computes the state vector from keplerian elements.

orbitdeterminator.util.new_tle_kep_state.kep_to_state(kep)[source]

This function converts from keplerian elements to the position and velocity vector

Parameters:kep (1x6 numpy array) – kep contains the following variables kep[0] = semi-major axis (kms) kep[1] = eccentricity (number) kep[2] = inclination (degrees) kep[3] = argument of perigee (degrees) kep[4] = right ascension of ascending node (degrees) kep[5] = true anomaly (degrees)

Returns: r: 1x6 numpy array which contains the position and velocity vector

r[0],r[1],r[2] = position vector [rx,ry,rz] km r[3],r[4],r[5] = velocity vector [vx,vy,vz] km/s
orbitdeterminator.util.new_tle_kep_state.tle_to_state(tle)[source]

This function converts from TLE elements to position and velocity vector

Parameters:tle (1x6 numpy array) – tle contains the following variables tle[0] = inclination (degrees) tle[1] = right ascension of the ascending node (degrees) tle[2] = eccentricity (number) tle[3] = argument of perigee (degrees) tle[4] = mean anomaly (degrees) tle[5] = mean motion (revs per day)

Returns: r: 1x6 numpy array which contains the position and velocity vector

r[0],r[1],r[2] = position vector [rx,ry,rz] km r[3],r[4],r[5] = velocity vector [vx,vy,vz] km/s

teme_to_ecef

Converts coordinates in TEME frame to ECEF frame.

orbitdeterminator.util.teme_to_ecef.conv_to_ecef(coords)[source]

Converts coordinates in TEME frame to ECEF frame.

Parameters:coords (nx4 numpy array) – list of coordinates in the format [t,x,y,z]
Returns:
list of coordinates in the format
[t, latitude, longitude, altitude]

Note that these coordinates are with respect to the surface of the Earth. Latitude, longitude are in degrees.

Return type:nx4 numpy array