Utils
Matrices
The rotation matrices defined here are taken from David Vallado’s Fundamentals of Astrodynamics and Applications. They allow to change reference frame, and as such are inverse of classical rotation matrices.
- beyond.utils.matrix.expand(m, rate=None)
Transform a 3x3 rotation matrix into a 6x6 rotation matrix
- Parameters:
m (numpy.ndarray) – 3x3 matrix transforming a position vector from frame1 to frame2
rate (numpy.array) – 1D 3 elements vector rate of frame2 expressed in frame1
- Returns:
6x6 rotation matrix
- Return type:
numpy.ndarray
Example:
>>> m = np.array([[0, -1, 0], [-1, 0, 0], [0, 0, 1]]) >>> print(expand(m)) [[ 0. -1. 0. 0. 0. 0.] [-1. 0. 0. 0. 0. 0.] [ 0. 0. 1. 0. 0. 0.] [ 0. 0. 0. 0. -1. 0.] [ 0. 0. 0. -1. 0. 0.] [ 0. 0. 0. 0. 0. 1.]] >>> print(expand(m, rate=[1,2,3])) [[ 0. -1. 0. 0. 0. 0.] [-1. 0. 0. 0. 0. 0.] [ 0. 0. 1. 0. 0. 0.] [ 3. 0. -1. 0. -1. 0.] [ 0. -3. 2. -1. 0. 0.] [ 2. -1. 0. 0. 0. 1.]]
- beyond.utils.matrix.rot1(theta)
- Parameters:
theta (float) – Angle in radians
- Returns:
Rotation matrix of angle theta around the X-axis
- beyond.utils.matrix.rot2(theta)
- Parameters:
theta (float) – Angle in radians
- Returns:
Rotation matrix of angle theta around the Y-axis
- beyond.utils.matrix.rot3(theta)
- Parameters:
theta (float) – Angle in radians
- Returns:
Rotation matrix of angle theta around the Z-axis
Interpolations
- class beyond.utils.interp.Interp(xs, ys, method, order=None)
Linear or Lagrangian interpolation
The constructed object is callable
Example
f = Interp(xs, ys, "lagrange", 8) y1 = f(x1) # interpolated value at x1
- __call__(x)
Compute interpolation at x
- Parameters:
x (float) – value to interpolate to
- Return type:
numpy.array
- __init__(xs, ys, method, order=None)
- Parameters:
xs (numpy.array) – 1-D array of real values. Shall be monotonically increasing.
ys (numpy.array) – 1-D or 2-D array of real values
method (str) – Selection of the interpolation method Either ‘linear’ or ‘lagrange’.
order (int) – Order of the interpolation. Has no effect on a linear interpolation.
- class beyond.utils.interp.DatedInterp(dates, ys, method, order=None)
Bases:
Interp
Interpolation for time series
- __call__(date)
Compute interpolation at a given date
- Parameters:
date (Date) – Date to interpolate to
- Return type:
numpy.array
- __init__(dates, ys, method, order=None)
- Parameters:
dates (list of Date) – 1-D array of dates. Shall be monotonically increasing.
ys (numpy.array) – N-D array of real values
method (str) – Selection of the interpolation method Either ‘linear’ or ‘lagrange’.
order (int) – Order of the interpolation. Has no effect on a linear interpolation.
Node
This module helps to find the shortest path between two element in a hierarchy or in a graph.
- class beyond.utils.node.Node(name)
Bases:
object
Class representing a node in a graph, relations may be circular.
A = Node('A') B = Node('B') C = Node('C') D = Node('D') E = Node('E') F = Node('F') A + B + C + D + E + F + A F + C # A # / \ # B F # | / | # C E # \ / # D A.path('E') # [A, F, E] A.steps('E') # [(A, F), (F, E)] E.path('B') # [E, F, A, B] or [E, D, C, B]
- name
Name of the node
- neighbors
List of all direct neighbors in the graph. OrderedDict is only used as OrderedSet, so only the keys of the dict matter
- path(goal)
Get the shortest way between two nodes of the graph
- Parameters:
goal (str) – Name of the targeted node
- Returns:
list of Node
- routes
Route mapping. What direction to follow in order to reach a particular target
- steps(goal)
Get the list of individual relations leading to the targeted node
- Parameters:
goal (str) – Name of the targeted node
- Returns:
list of tuple of Node
Beta
The beta angle is the angle between the plane of the orbit and the direction of a distant body (see wikipedia).
If the distant body is the Sun, it will be useful to characterize the spacecraft illumination during its orbit. If the distant body is another spacecraft, it will be useful to compute the visibility between the two spacecrafts.
- beyond.utils.beta.beta(orb, ref='Sun')
Compute beta angle
- beyond.utils.beta.beta_limit(orb)
Compute minimal beta angle for a constant visibility on another body.
Below this threshold the spacecraft will experience eclipses during a portion of its orbit. Above, it will be fully illuminated during all its orbit.
- Parameters:
orb (Orbit) – Orbit of the primary spacecraft, expressed in a reference frame whose centre is the obscuring body.
- Returns:
Angle in radians
- Return type:
float
LTAN
Utilities to compute the Local Time at Ascending Node (LTAN)
Both True and Mean LTAN are available, the difference between them being the equation of time
- beyond.utils.ltan.ltan2raan(date, ltan, type='mean')
Conversion to Longitude
- Parameters:
date (Date) – Date of the conversion
ltan (float) – LTAN in seconds
type (str) – either “mean” or “true”
- Returns:
RAAN in radians in EME2000
- Return type:
float
Constellation
Utilities to compute the parameters of a constellation.
At the moment, only the Walker Star and Walker Delta are available (see wikipedia)
- class beyond.utils.constellation.WalkerDelta(total, planes, spacing, raan0=0)
Bases:
WalkerStar
Definition of the Walkek Delta constellation
Example: Galileo is a Walker Delta 24/3/1 constellation so to generate this, one has to call
WalkerDelta(24, 3, 1)
- nu(i_plane, i_sat)
- Parameters:
i_plane (int) – index of the plane
i_sat (int) – index of the satellite
- Returns:
True anomaly in radians
- Return type:
float
- raan(i_plane)
- Parameters:
i_plane (int) – index of the plane
- Returns:
Right Ascension of Ascending Node in radians
- Return type:
float
- class beyond.utils.constellation.WalkerStar(total, planes, spacing, raan0=0)
Bases:
object
Definition of the WalkerStar constellation
Example: Iridium is a Walker Star 66/6/2 constellation so to generate this, one has to call
WalkerStar(66, 6, 2)
- __init__(total, planes, spacing, raan0=0)
- Parameters:
total (int) – Total number of satellites
planes (int) – Number of planes
spacing (int) – relative spacing between satellites of adjacent planes
raan0 (float) – RAAN of the first plane (in radians)
This call order is compliant with Walker notation total/planes/spacing.
- nu(i_plane, i_sat)
- Parameters:
i_plane (int) – index of the plane
i_sat (int) – index of the satellite
- Returns:
True anomaly in radians
- Return type:
float
- property per_plane
Number of satellites per orbital plane
- raan(i_plane)
- Parameters:
i_plane (int) – index of the plane
- Returns:
Right Ascension of Ascending Node in radians
- Return type:
float
LEO
Utilities for Low Earth Orbit design
- beyond.utils.leo.frozen(a, i)
Compute the e and ω for a frozen orbit
- Parameters:
a (float) – Semi major axis in meters
i (float) – Inclination in radians
- Returns:
eccentricity (e) and argument of perigee (ω, in radians)
- Return type:
Tuple[float, float]
- beyond.utils.leo.sso(*, a=None, e=None, i=None)
Compute keplerian elements for a Sun-Synchronous Orbit
Given two elements among a, e and i, compute the remaining one.
Example
e = sso(a=a, i=i) i = sso(a=a, e=e) a = sso(e=e, i=i)
- Parameters:
a (float) – Semi major axis, in meters
e (float) – Eccentricity
i (float) – Inclination in radians
- beyond.utils.leo.sso_frozen(a)
Iterate to find a SSO frozen orbit
- Parameters:
a (float) – Semi major axis in meters
- Returns:
eccentricity (e), inclination (i, in radians) and argument of perigee (ω, in radians)
- Return type:
Tuple[float, float, float]
Lambert’s problem
Lambert’s problem solvers
- beyond.utils.lambert.lambert(orb0, orb1, prograde=True)
Solve Lamber’s problem, with the solution provided by Howard D. Curtis in chapter 5.3 of his book, “Orbital Mechanics for Engineering Students” (ed. 2014)
- Parameters:
- Return
- TupleInitial and Target
Orbits
patched with solution’s velocities
- TupleInitial and Target
The initial orbit reference frame is the one used for the computation. So if one desires to compute an interplanetary opportunity, orb0 should be expressed in the “Sun” or “SolarSystemBarycenter” reference frame.
Warning
This is only compatible with elliptical orbits
Interplanetary
Module to compute planetary captures and flybys
- class beyond.utils.interplanetary.BPlane(B, theta, S, T, R, e, h)
B-plane characteristics
- Parameters:
B – 3d B vector
theta – angle
S – vector collinear to v_infinity
T – vector along the ecliptic
R – S x T
e – 3d eccentricity vector
h – angular momentum vector
- beyond.utils.interplanetary.bplane(orb)
Compute the B vector details
- beyond.utils.interplanetary.flyby(v_inf_in, v_inf_out, μ)
Compute B vector characteristics for a flyby
- Parameters:
v_inf_in (numpy.array, shape 3) – arrival excess velocity
v_inf_out (numpy.array, shape 3) – exit excess velocity
µ (float) – standard gravitational parameter of the body flown-by
- Returns:
- tuple with length 2, the first element being the B vector
(a numpy array), the second being the periapsis radius
- Return type:
tuple
Rendez-vous and Proximity Operations Helper
- class beyond.utils.rpohelper.RpoHelper(propagator)
This class provides computation helpers for relative motion positioning and maneuvers, to be used with the
rpo
propagatorsSee the Docking script in the documentation for an example of utilisation.
- coelliptic(date, radial, tangential)
Create an Orbit instance at a given radial and tangential distance with the guarantee that it’s coelliptic
- coelliptic_velocity(radial)
- Parameters:
radial (float) – radial distance between the target and the chaser
- Returns:
Necessary tangential velocity for a coelliptic orbit
- Return type:
float
- eccentric_boost(tangential, date, continuous=False)
Perform an eccentric boost to move tangentially, when the radial distance is zero.
- Parameters:
tangential (float) – Tangential distance to cover
date (Date) – Begining of the eccentric boost
continuous (bool)
- Returns:
List[Man]
- hohmann(radial, date, continuous=False)
Perform a Hohmann transfer
- Parameters:
radial (float) – Radial distance to cover
date (Date) – Begining of the Hohmann transfer
continuous (bool)
- Returns:
List[Man]
- hohmann_distance(radial, continuous=False)
Compute the tangential distance traveled during a hohmann transfer This is interesting to anticipate and place the arrival at a desired position.
- Parameters:
radial (float) – radial distance from the target to the chaser
continuous (bool) – The Hohmann transfer will be done with a ContinuousMan object
- Returns:
tangential distance
- Return type:
float
- property period
Period of the target orbit