gooding_module Module

Gooding's Kepler and universal elements conversion routines.

Notes

The Gooding universal elements are:

  • alpha - mu/a [km^2/s^2]
  • rp - periapsis radius [km]
  • inc - inclination [rad]
  • raan - right ascension of the ascending node [rad]
  • w - argument of periapsis [rad]
  • tau - time since last periapsis passage [sec]

References

  1. A. W. Odell, R. H. Gooding, "Procedures for solving Kepler's equation" Celestial Mechanics 38 (1986), 307-334.
  2. R. H. Gooding, "On universal elements, and conversion procedures to and from position and velocity" Celestial Mechanics 44 (1988), 283-298.
  3. R. H. Gooding, A. W. Odell. "The hyperbolic Kepler equation (and the elliptic equation revisited)" Celestial Mechanics 44 (1988), 267-282.

Uses

  • module~~gooding_module~~UsesGraph module~gooding_module gooding_module module~kind_module kind_module module~gooding_module->module~kind_module module~numbers_module numbers_module module~gooding_module->module~numbers_module iso_fortran_env iso_fortran_env module~kind_module->iso_fortran_env module~numbers_module->module~kind_module

Used by

  • module~~gooding_module~~UsedByGraph module~gooding_module gooding_module module~fortran_astrodynamics_toolkit fortran_astrodynamics_toolkit module~fortran_astrodynamics_toolkit->module~gooding_module proc~lambert_test lambert_module::lambert_test proc~lambert_test->module~gooding_module

Variables

Type Visibility Attributes Name Initial
real(kind=wp), private, parameter :: ntwo = -two
real(kind=wp), private, parameter :: pineg = -pi
real(kind=wp), private, parameter :: half = 0.5_wp
real(kind=wp), private, parameter :: athird = one/three
real(kind=wp), private, parameter :: asixth = one/six

Functions

public pure function ekepl(em, e1)

Kepler's equation, em = ekepl - (1 - e1)*sin(ekepl), with e1 in range 1 to 0 inclusive, solved accurately (based on ekepl3, but entering e1, not e)

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: em
real(kind=wp), intent(in) :: e1

Return Value real(kind=wp)

public pure function ekepl1(em, e)

Solve kepler's equation, em = ekepl - e*sin(ekepl), with legendre-based starter and halley iterator (function has also been used under the name eafkep)

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: em
real(kind=wp), intent(in) :: e

Return Value real(kind=wp)

public pure function ekepl2(em, e)

Kepler's equation, em = ekepl - e*sin(ekepl) with e in range 0 to 1 inclusive, solved accurately

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: em
real(kind=wp), intent(in) :: e

Return Value real(kind=wp)

public pure function emkepl(e, ee)

Accurate computation of ee - e*sin(ee) when (e, ee) is close to (1, 0)

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: e
real(kind=wp), intent(in) :: ee

Return Value real(kind=wp)

public pure function emkep(e1, ee)

Similar to emkepl, except input is 1-e.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: e1
real(kind=wp), intent(in) :: ee

Return Value real(kind=wp)

public pure function shkepl(el, g1)

Equation el = shkepl + (g1 - 1)*asinh(shkepl), with g1 in range 0 to 1 inclusive, solved accurately.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: el
real(kind=wp), intent(in) :: g1

Return Value real(kind=wp)

public pure function shmkep(g1, s)

Accurate computation of s - (1 - g1)*asinh(s) when (g1, s) is close to (0, 0)

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: g1
real(kind=wp), intent(in) :: s

Return Value real(kind=wp)

private pure function dcbsol(a, b, c) result(x)

Solution to a*x**3 + 3*b*x - 2c = 0, where a and b**3 + a*c**2 are both non-negative (zero generated, in lieu of infinity, if a = b = 0)

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: a
real(kind=wp), intent(in) :: b
real(kind=wp), intent(in) :: c

Return Value real(kind=wp)

private pure function dcubrt(x) result(c)

Cube root computed accurately, by incorporating one Newton-Raphson iteration.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x

Return Value real(kind=wp)


Subroutines

public pure subroutine propagate(mu, rv0, dt, rvf)

Author
Jacob Williams

Basic two-body propagator using the Gooding universal element routines.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: mu

grav. parameter [km^3/s^2]

real(kind=wp), intent(in), dimension(6) :: rv0

initial state [km, km/s]

real(kind=wp), intent(in) :: dt

time step [sec]

real(kind=wp), intent(out), dimension(6) :: rvf

final state [km, km/s]

private pure subroutine els2pv(gm, al, q, om, tau, r, u, vr, vt)

Algorithm for two-dimensional conversion from orbital elements to position and velocity.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: gm

grav. parameter [km^3/s^2]

real(kind=wp), intent(in) :: al

alpha [km^2/s^2]

real(kind=wp), intent(in) :: q

periapsis distance [km]

real(kind=wp), intent(in) :: om

argument of periapsis relative to assumed reference direction [rad]

real(kind=wp), intent(in) :: tau

time from periapsis [sec]

real(kind=wp), intent(out) :: r

radial distance [km]

real(kind=wp), intent(out) :: u

angle from reference direction [rad]

real(kind=wp), intent(out) :: vr

radial velocity [km/2]

real(kind=wp), intent(out) :: vt

transverse velocity >=0 [km/s]

public pure subroutine els3pv(gm, e, pv)

Algorithm for three-dimensional conversion from orbital elements to position and velocity

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: gm

grav. parameter [km^3/sec^2]

real(kind=wp), intent(in), dimension(6) :: e

[al, q, ei, bom, om, tau]

real(kind=wp), intent(out), dimension(6) :: pv

[x, y, z, xdot, ydot, zdot]

private pure subroutine pv2els(gm, r, u, vr, vt, al, q, om, tau)

Algorithm for two-dimensional conversion from position and velocity to orbital elements.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: gm

grav. parameter [km^3/s^2]

real(kind=wp), intent(in) :: r

radial distance [km]

real(kind=wp), intent(in) :: u

angle from assumed reference direction [rad]

real(kind=wp), intent(in) :: vr

radial velocity [km/2]

real(kind=wp), intent(in) :: vt

transverse velocity >=0 [km/s]

real(kind=wp), intent(out) :: al

alpha: gm/a [km^2/s^2]

real(kind=wp), intent(out) :: q

periapsis distance [km]

real(kind=wp), intent(out) :: om

argument of periapsis relative to reference direction [rad]

real(kind=wp), intent(out) :: tau

time from periapsis [sec]

public pure subroutine pv3els(gm, pv, e)

Algorithm for three-dimensional conversion from position and velocity to orbital elements.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: gm

grav. parameter [km^3/s^2]

real(kind=wp), intent(in), dimension(6) :: pv

[x, y, z, xdot, ydot, zdot]

real(kind=wp), intent(out), dimension(6) :: e

[al, q, ei, bom, om, tau]