orbitalMotion
- orbitalMotion.E2M(Ecc, e)[source]
Maps the eccentric anomaly angle into the corresponding mean elliptic anomaly angle. Both 2D and 1D elliptic orbit are allowed.
- Parameters:
Ecc – eccentric anomaly (rad)
e – eccentricity (0 <= e < 1)
- Returns:
M, mean elliptic anomaly (rad)
- orbitalMotion.E2f(Ecc, e)[source]
Maps eccentric anomaly angles into true anomaly angles This function requires the orbit to be either circular or non-rectilinear elliptic orbit
- Parameters:
Ecc – eccentric anomaly (rad)
e – eccentric (0 <= e < 1)
- Returns:
f, true anomaly (rad)
- orbitalMotion.H2N(H, e)[source]
Maps the hyperbolic anomaly angle H into the corresponding mean hyperbolic anomaly angle N.
- Parameters:
H – hyperbolic anomaly (rad)
e – eccentricity (e > 1)
- Returns:
N, mean hyperbolic anomaly (rad)
- orbitalMotion.H2f(H, e)[source]
Maps hyperbolic anomaly angles into true anomaly angles. This function requires the orbit to be hyperbolic
- Parameters:
H – hyperbolic anomaly (rad)
e – eccentricity (e > 1)
- Returns:
f, true anomaly angle (rad)
- orbitalMotion.M2E(M, e)[source]
Maps the mean elliptic anomaly angle into the corresponding eccentric anomaly angle. Both 2D and 1D elliptic orbit are allowed.
- Parameters:
M – mean elliptic anomaly (rad)
e – eccentricity (0 <= e < 1)
- Returns:
Ecc, eccentric anomaly (rad)
- orbitalMotion.N2H(N, e)[source]
Maps the mean hyperbolic anomaly angle N into the corresponding hyperbolic anomaly angle H.
- Parameters:
N – mean hyperbolic anomaly (rad)
e – eccentricity (e > 1)
- Returns:
H, hyperbolic anomaly (rad)
- orbitalMotion.atmosphericDensity(alt)[source]
This program computes the atmospheric density based on altitude supplied by user. This function uses a curve fit based on atmospheric data from the Standard Atmosphere 1976 Data. This function is valid for altitudes ranging from 100km to 1000km.
Note
This code can only be applied to spacecraft orbiting the Earth
- Parameters:
alt – altitude in km
- Returns:
density at the given altitude in kg/m^3
- orbitalMotion.atmosphericDrag(Cd, A, m, rvec, vvec)[source]
This program computes the atmospheric drag acceleration vector acting on a spacecraft. Note the acceleration vector output is inertial, and is only valid for altitudes up to 1000 km. Afterwards the drag force is zero. Only valid for Earth.
- Parameters:
Cd – drag coefficient of the spacecraft
A – cross-sectional area of the spacecraft in m^2
m – mass of the spacecraft in kg
rvec – Inertial position vector of the spacecraft in km [x;y;z]
vvec – Inertial velocity vector of the spacecraft in km/s [vx;vy;vz]
- Returns:
The inertial acceleration vector due to atmospheric drag in km/sec^2
- orbitalMotion.clElem2eqElem(elements_cl, elements_eq)[source]
conversion from classical orbital elements (a,e,i,Omega,omega,f) to equinoctial orbital elements (a,P1,P2,Q1,Q2,l,L)
- Parameters:
elements_cl – classical elements
- Returns:
elements_eq, equinoctial elements
- orbitalMotion.clMeanOscMap(req, J2, oe, oep, sign)[source]
First-order J2 Mapping Between Mean and Osculating Orbital Elements
Analytical Mechanics of Space Systems, Hanspeter Schaub, John L. Junkins, 4th edition. [m] or [km] should be the same both for req and elements.a
- Parameters:
req – equatorial radius
J2
oe – classical orbit element set
oep
sign – sgn=1:mean to osc, sgn=-1:osc to mean
- orbitalMotion.debyeLength(alt)[source]
This program computes the debyeLength length for a given altitude and is valid for altitudes ranging from 200 km to GEO (35000km). However, all values above 1000 km are HIGHLY speculative at this point.
- Parameters:
alt – altitude in km
- Returns:
debye length given in m
- orbitalMotion.elem2rv(mu, elements)[source]
Translates the orbit elements:
a
semi-major axis
km
e
eccentricity
i
inclination
rad
AN
ascending node
rad
AP
argument of periapses
rad
f
true anomaly angle
rad
to the inertial Cartesian position and velocity vectors. The attracting body is specified through the supplied gravitational constant mu (units of km^3/s^2).
- Parameters:
mu – gravitational parameter
elements – orbital elements
- Returns:
rVec, position vector
- Returns:
vVec, velocity vector
- orbitalMotion.elem2rv_parab(mu, elements)[source]
Translates the orbit elements:
a
semi-major axis
km
e
eccentricity
i
inclination
rad
AN
ascending node
rad
AP
argument of periapses
rad
f
true anomaly angle
rad
to the inertial Cartesian position and velocity vectors. The attracting body is specified through the supplied gravitational constant mu (units of km^3/s^2).
The code can handle the following cases:
circular:
e = 0
a > 0
elliptical-2D:
0 < e < 1
a > 0
elliptical-1D:
e = 1
a > 0
f = Ecc. Anom. here
parabolic:
e = 1
rp = -a
hyperbolic:
e > 1
a < 0
Note
To handle the parabolic case and distinguish it form the rectilinear elliptical case, instead of passing along the semi-major axis a in the “a” input slot, the negative radius at periapses is supplied. Having “a” be negative and e = 1 is a then a unique identified for the code for the parabolic case.
- Parameters:
mu – gravitational parameter
elements – orbital elements
- Returns:
rVec = position vector, vVec = velocity vector
- orbitalMotion.f2E(f, e)[source]
Maps true anomaly angles into eccentric anomaly angles. This function requires the orbit to be either circular or non-rectilinear elliptic orbit.
- Parameters:
f – true anomaly angle (rad)
e – eccentricity (0 <= e < 1)
- Returns:
Ecc, eccentric anomaly (rad)
- orbitalMotion.f2H(f, e)[source]
Maps true anomaly angles into hyperbolic anomaly angles. This function requires the orbit to be hyperbolic
- Parameters:
f – true anomaly angle (rad)
e – eccentricity (e > 1)
- Returns:
H, hyperbolic anomaly (rad)
- orbitalMotion.hill2rv(rc_N, vc_N, rho_H, rhoPrime_H)[source]
Map the deputy position and velocity vector relative to the chief Hill frame to inertial frame.
- Parameters:
rc_N – chief inertial position vector
vc_N – chief inertial velocity vector
rho_H – deputy Hill relative position vector
rhoPrime_H – deputy Hill relative velocity vector
- Returns:
rd_N, vd_N: Deputy inertial position and velocity vectors
- orbitalMotion.hillFrame(rc_N, vc_N)[source]
Compute the Hill frame DCM HN :param rc_N: inertial position vector :param vc_N: inertial velocity vector :return: HN: DCM that maps from the inertial frame N to the Hill (i.e. orbit) frame H
- orbitalMotion.jPerturb(rvec, num, planet)[source]
Computes the J2_EARTH-J6_EARTH zonal gravitational perturbation accelerations.
- Parameters:
rvec – Cartesian Position vector in kilometers [x;y;z].
num – Corresponds to which J components to use, must be an integer between 2 and 6. (note: Additive- 2 corresponds to J2_EARTH while 3 will correspond to J2_EARTH + J3_EARTH)
planet – planet variable, can be CELESTIAL_MERCURY CELESTIAL_VENUS CELESTIAL_EARTH CELESTIAL_MOON CELESTIAL_MARS CELESTIAL_JUPITER CELESTIAL_URANUS CELESTIAL_NEPTUNE
- Returns:
ajtot, The total acceleration vector due to the J perturbations in km/sec^2 [accelx;accely;accelz]
- orbitalMotion.rv2elem(mu, rVec, vVec)[source]
Translates the orbit elements inertial Cartesian position vector rVec and velocity vector vVec into the corresponding classical orbit elements where
a
semi-major axis
km
e
eccentricity
i
inclination
rad
AN
ascending node
rad
AP
argument of periapses
rad
f
true anomaly angle
rad
If the orbit is rectilinear, then this will be the eccentric or hyperbolic anomaly
The attracting body is specified through the supplied gravitational constant mu (units of km^3/s^2).
- Parameters:
mu – gravitational parameter
rVec – position vector
vVec – velocity vector
- Returns:
orbital elements
- orbitalMotion.rv2elem_parab(mu, rVec, vVec)[source]
Translates the orbit elements inertial Cartesian position vector rVec and velocity vector vVec into the corresponding classical orbit elements where
a
semi-major axis
km
e
eccentricity
i
inclination
rad
AN
ascending node
rad
AP
argument of periapses
rad
f
true anomaly angle
rad
If the orbit is rectilinear, then f will be the eccentric or hyperbolic anomaly
The attracting body is specified through the supplied gravitational constant mu (units of km^3/s^2).
The code can handle the following cases:
circular:
e = 0
a > 0
elliptical-2D:
0 < e < 1
a > 0
elliptical-1D:
e = 1
a > 0
parabolic:
e = 1
a = -rp
hyperbolic:
e > 1
a < 0
For the parabolic case the semi-major axis is not defined. In this case -rp (radius at periapses) is returned instead of a. For the circular case, the AN and AP are ill-defined, along with the associated ie and ip unit direction vectors of the perifocal frame. In this circular orbit case, the unit vector ie is set equal to the normalized inertial position vector ir.
- Parameters:
mu – gravitational parameter
rVec – position vector
vVec – velocity vector
- Returns:
orbital elements
Todo: Modify this code to return true longitude of periapsis (non-circular, equatorial), argument of latitude (circular, inclined), and true longitude (circular, equatorial) when appropriate instead of simply zeroing out omega and Omega
- orbitalMotion.rv2hill(rc_N, vc_N, rd_N, vd_N)[source]
Express the deputy position and velocity vector as chief by the chief Hill frame.
- Parameters:
rc_N – chief inertial position vector
vc_N – chief inertial velocity vector
rd_N – chief deputy position vector
vd_N – chief deputy velocity vector
- Returns:
rho_H, rhoPrime_H: Hill frame relative position and velocity vectors
- orbitalMotion.solarRad(A, m, sunvec)[source]
Computes the inertial solar radiation force vectors based on cross-sectional Area and mass of the spacecraft and the position vector of the planet to the sun.
Note
It is assumed that the solar radiation pressure decreases quadratically with distance from sun (in AU)
Solar Radiation Equations obtained from Earth Space and Planets Journal Vol. 51, 1999 pp. 979-986
- Parameters:
A – Cross-sectional area of the spacecraft that is facing the sun in m^2.
m – The mass of the spacecraft in kg.
sunvec – Position vector to the Sun in units of AU. Earth has a distance of 1 AU.
- Returns:
arvec, The inertial acceleration vector due to the effects of Solar Radiation pressure in km/sec^2. The vector components of the output are the same as the vector components of the sunvec input vector.