orbitalMotion¶
-
orbitalMotion.
D2R
= 0.017453292519943295¶ Gravitational Constants mu = G*m, where m is the planet of the attracting planet. All units are km^3/s^2. Values are obtained from SPICE kernels in http://naif.jpl.nasa.gov/pub/naif/generic_kernels/
-
orbitalMotion.
E2M
(Ecc, e)[source]¶ Function: E2M Purpose: Maps the eccentric anomaly angle into the corresponding
mean elliptic anomaly angle. Both 2D and 1D elliptic orbit are allowed.
- Inputs:
Ecc = eccentric anomaly (rad) e = eccentricity (0 <= e < 1)
- Outputs:
M = mean elliptic anomaly (rad)
-
orbitalMotion.
E2f
(Ecc, e)[source]¶ Function: E2f Purpose: Maps eccentric anomaly angles into true anomaly angles
This function requires the orbit to be either circular or non-rectilinar elliptic orbit
- Inputs:
Ecc = eccentric anomaly (rad) e = eccentric (0 <= e < 1)
- Outputs:
f = true anomaly (rad)
-
orbitalMotion.
H2N
(H, e)[source]¶ Function: H2N Purpose: Maps the hyperbolic anomaly angle H into the corresponding
mean hyperbolic anomaly angle N.
- Inputs:
H = hyperbolic anomaly (rad) e = eccentricity (e > 1)
- Outputs:
N = mean hyperbolic anomaly (rad)
-
orbitalMotion.
H2f
(H, e)[source]¶ Function: H2f Purpose: Maps hyperbolic anomaly angles into true anomaly angles.
This function requires the orbit to be hyperbolic
- Inputs:
H = hyperbolic anomaly (rad) e = eccentricity (e > 1)
- Outputs:
f = true anomaly angle (rad)
-
orbitalMotion.
M2E
(M, e)[source]¶ - Purpose: Maps the mean elliptic anomaly angle into the corresponding
eccentric anomaly angle. Both 2D and 1D elliptic orbit are allowed.
- Inputs:
M = mean elliptic anomaly (rad) e = eccentricity (0 <= e < 1)
- Outputs:
Ecc = eccentric anomaly (rad)
-
orbitalMotion.
MU_PLUTO
= 983.055¶ planet information for major solar system bodies. Units are in km. data taken from http://nssdc.gsfc.nasa.gov/planetary/planets.html
-
orbitalMotion.
N2H
(N, e)[source]¶ Function: N2H Purpose: Maps the mean hyperbolic anomaly angle N into the corresponding
hyperbolic anomaly angle H.
- Inputs:
N = mean hyperbolic anomaly (rad) e = eccentricity (e > 1)
- Outputs:
H = hyperbolic anomaly (rad)
-
orbitalMotion.
atmosphericDensity
(alt)[source]¶ Function: atmosphericDensity
Purpose: 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
Inputs:
alt = altitude in km
Outputs:
density = density at the given altitude in kg/m^3
-
orbitalMotion.
atmosphericDrag
(Cd, A, m, rvec, vvec)[source]¶ Function: atmosphericDrag
Purpose: 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.
Inputs:
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]
Outputs:
advec = The inertial acceleration vector due to atmospheric
drag in km/sec^2
-
orbitalMotion.
debyeLength
(alt)[source]¶ Function: debyeLength
Purpose: 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.
Inputs:
alt = altitude in km
Outputs:
debye = debye length given in m
-
orbitalMotion.
elem2rv
(mu, elements)[source]¶ Function: elem2rv Purpose: 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).
- Inputs:
mu = gravitational parameter elements = orbital elements
- Outputs:
rVec = position vector vVec = velocity vector
-
orbitalMotion.
elem2rv_parab
(mu, elements)[source]¶ Function: elem2rv Purpose: 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.
- Inputs:
mu = gravitational parameter elements = orbital elements
- Outputs:
rVec = position vector vVec = velocity vector
-
orbitalMotion.
f2E
(f, e)[source]¶ Function: f2E Purpose: Maps true anomaly angles into eccentric anomaly angles.
This function requires the orbit to be either circular or non-rectilinar elliptic orbit.
- Inputs:
f = true anomaly angle (rad) e = eccentricity (0 <= e < 1)
- Outputs:
Ecc = eccentric anomaly (rad)
-
orbitalMotion.
f2H
(f, e)[source]¶ Function: f2H Purpose: Maps true anomaly angles into hyperbolic anomaly angles.
This function requires the orbit to be hyperbolic
- Inputs:
f = true anomaly angle (rad) e = eccentricity (e > 1)
- Outputs:
H = hyperbolic anomaly (rad)
-
orbitalMotion.
jPerturb
(rvec, num, planet)[source]¶ Function: jPerturb
Purpose: Computes the J2_EARTH-J6_EARTH zonal graviational perturbation
accelerations.
Inputs:
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)
Outputs:
ajtot = The total acceleration vector due to the J
perturbations in km/sec^2 [accelx;accely;accelz]
-
orbitalMotion.
rv2elem
(mu, rVec, vVec)[source]¶ Function: rv2elem Purpose: 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).
- Inputs:
mu = gravitational parameter rVec = position vector vVec = velocity vector
- Outputs:
elements = orbital elements
-
orbitalMotion.
rv2elem_parab
(mu, rVec, vVec)[source]¶ Function: rv2elem Purpose: 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).
- 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.
- Inputs:
mu = gravitational parameter rVec = position vector vVec = velocity vector
- Outputs:
elements = orbital elements
TODO: (SAO) 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.
solarRad
(A, m, sunvec)[source]¶ Function: solarRad
Purpose: 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
Inputs:
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.
Outputs:
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.