orbitalMotion

Defines

N_DEBYE_PARAMETERS

Enums

enum CelestialObject_t

Values:

enumerator CELESTIAL_MERCURY
enumerator CELESTIAL_VENUS
enumerator CELESTIAL_EARTH
enumerator CELESTIAL_MOON
enumerator CELESTIAL_MARS
enumerator CELESTIAL_PHOBOS
enumerator CELESTIAL_DEIMOS
enumerator CELESTIAL_JUPITER
enumerator CELESTIAL_SATURN
enumerator CELESTIAL_URANUS
enumerator CELESTIAL_NEPTUNE
enumerator CELESTIAL_PLUTO
enumerator CELESTIAL_SUN
enumerator MAX_CELESTIAL

Functions

double E2f(double E, double e)

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 = eccentricity (0 <= e < 1) Outputs: f = true anomaly (rad)

double E2M(double E, double e)

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)

double f2E(double f, double e)

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)

double f2H(double f, double e)

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)

double H2f(double H, double e)

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)

double H2N(double H, double e)

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)

double M2E(double M, double e)

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)

double N2H(double N, double e)

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)

void elem2rv(double mu, classicElements *elements, double *rVec, double *vVec)

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

void rv2elem(double mu, double *rVec, double *vVec, classicElements *elements)

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 (zero if parabolic) 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 rp - radius at periapses ra - radius at apoapses (zero if parabolic) 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

void clMeanOscMap(double req, double J2, classicElements *elements, classicElements *elements_p, double sgn)

maps classical mean orbit elements to Osculating elements

void clElem2eqElem(classicElements *elements_cl, equinoctialElements *elements_eq)

maps from classical orbit elements to equinoctial elements

void hillFrame(double *rc_N, double *vc_N, double HN[3][3])

Purpose: maps inertial position and velocity vectors in the Hill frame DCM HN Inputs: rc_N: inertial position vector vc_N: inertial velocity vector Outputs: HN: Hill frame DCM relative to inertial frame

void hill2rv(double *rc_N, double *vc_N, double *rho_H, double *rhoPrime_H, double *rd_N, double *vd_N)

Purpose: maps Hill frame deputy states to inertial inertial position and velocity vectors Inputs: rc_N: chief inertial position vector vc_N: chief inertial velocity vector rho_H: deputy Hill position vector rhoPrime_H: deputy Hill velocity vector Outputs: rd_N: deputy inertial position vector vd_N: deputy inertial velocity vector

void rv2hill(double *rc_N, double *vc_N, double *rd_N, double *vd_N, double *rho_H, double *rhoPrime_H)

Purpose: maps inertial frame deputy states to Hill inertial position and velocity vectors Inputs: rc_N: chief inertial position vector vc_N: chief inertial velocity vector rd_N: deputy inertial position vector vd_N: deputy inertial velocity vector Outputs: rho_H: deputy Hill position vector rhoPrime_H: deputy Hill velocity vector

double atmosphericDensity(double alt)

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

double debyeLength(double alt)

Purpose: This program computes the Debye Length 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

void atmosphericDrag(double Cd, double A, double m, double *rvec, double *vvec, double *advec)

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

void jPerturb(double *rvec, int num, double *ajtot, ...)

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]

void solarRad(double A, double m, double *sunvec, double *arvec)

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.

struct equinoctialElements
#include <orbitalMotion.h>

This structure contains the set of Keplerian orbital elements that define the spacecraft translational state. It is operated on by the orbital element routines and the OrbElemConvert module.

equinoctial elment struct definition

Public Members

double a

semi-major axis

double P1

e*sin(omega+Omega)

double P2

e*cos(omega+Omega)

double Q1

tan(i/2)*sin(Omega)

double Q2

tan(i/2)*cos(Omega)

double l

Omega+omega+M.

double L

Omega+omega+f.