libnova  v 0.16.0
Functions
Elliptic Motion

Functions

double ln_solve_kepler (double e, double M)
 Calculate the eccentric anomaly. More...
 
double ln_get_ell_mean_anomaly (double n, double delta_JD)
 Calculate the mean anomaly. More...
 
double ln_get_ell_true_anomaly (double e, double E)
 Calculate the true anomaly. More...
 
double ln_get_ell_radius_vector (double a, double e, double E)
 Calculate the radius vector. More...
 
double ln_get_ell_smajor_diam (double e, double q)
 Calculate the semi major diameter. More...
 
double ln_get_ell_sminor_diam (double e, double a)
 Calculate the semi minor diameter. More...
 
double ln_get_ell_mean_motion (double a)
 Calculate the mean daily motion (degrees/day). More...
 
void ln_get_ell_geo_rect_posn (struct ln_ell_orbit *orbit, double JD, struct ln_rect_posn *posn)
 Calculate the objects rectangular geocentric position. More...
 
void ln_get_ell_helio_rect_posn (struct ln_ell_orbit *orbit, double JD, struct ln_rect_posn *posn)
 Calculate the objects rectangular heliocentric position. More...
 
double ln_get_ell_orbit_len (struct ln_ell_orbit *orbit)
 Calculate the orbital length in AU. More...
 
double ln_get_ell_orbit_vel (double JD, struct ln_ell_orbit *orbit)
 Calculate orbital velocity in km/s. More...
 
double ln_get_ell_orbit_pvel (struct ln_ell_orbit *orbit)
 Calculate orbital velocity at perihelion in km/s. More...
 
double ln_get_ell_orbit_avel (struct ln_ell_orbit *orbit)
 Calculate the orbital velocity at aphelion in km/s. More...
 
double ln_get_ell_body_phase_angle (double JD, struct ln_ell_orbit *orbit)
 Calculate the phase angle of the body. The angle Sun - body - Earth. More...
 
double ln_get_ell_body_elong (double JD, struct ln_ell_orbit *orbit)
 Calculate the bodies elongation to the Sun.. More...
 
double ln_get_ell_body_solar_dist (double JD, struct ln_ell_orbit *orbit)
 Calculate the distance between a body and the Sun. More...
 
double ln_get_ell_body_earth_dist (double JD, struct ln_ell_orbit *orbit)
 Calculate the distance between a body and the Earth. More...
 
void ln_get_ell_body_equ_coords (double JD, struct ln_ell_orbit *orbit, struct ln_equ_posn *posn)
 Calculate a bodies equatorial coords. More...
 
int ln_get_ell_body_rst (double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an elliptic orbit. More...
 
int ln_get_ell_body_rst_horizon (double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, double horizon, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an elliptic orbit. More...
 
int ln_get_ell_body_next_rst (double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an elliptic orbit. More...
 
int ln_get_ell_body_next_rst_horizon (double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, double horizon, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an elliptic orbit. More...
 
int ln_get_ell_body_next_rst_horizon_future (double JD, struct ln_lnlat_posn *observer, struct ln_ell_orbit *orbit, double horizon, int day_limit, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an elliptic orbit.
 
double ln_get_ell_last_perihelion (double epoch_JD, double M, double n)
 Calculate the julian day of the last perihelion. More...
 

Detailed Description

Functions relating to the elliptic motion of bodies.

All angles are expressed in degrees.

Function Documentation

double ln_get_ell_body_earth_dist ( double  JD,
struct ln_ell_orbit orbit 
)

Calculate the distance between a body and the Earth.

Parameters
JDJulian day.
orbitOrbital parameters
Returns
Distance in AU

Calculate the distance between an body and the Earth for the given julian day.

Examples:
asteroid.c, and comet.c.

References ln_get_ell_geo_rect_posn(), ln_rect_posn::X, ln_rect_posn::Y, and ln_rect_posn::Z.

Referenced by ln_get_ell_body_phase_angle().

double ln_get_ell_body_elong ( double  JD,
struct ln_ell_orbit orbit 
)

Calculate the bodies elongation to the Sun..

Parameters
JDJulian day
orbitOrbital parameters
Returns
Elongation to the Sun.

Calculate the bodies elongation to the Sun..

Examples:
asteroid.c.

References ln_ell_orbit::a, ln_ell_orbit::e, ln_ell_orbit::JD, ln_get_earth_solar_dist(), ln_get_ell_body_solar_dist(), ln_get_ell_mean_anomaly(), ln_get_ell_mean_motion(), ln_get_ell_radius_vector(), ln_rad_to_deg(), ln_range_degrees(), ln_solve_kepler(), and ln_ell_orbit::n.

void ln_get_ell_body_equ_coords ( double  JD,
struct ln_ell_orbit orbit,
struct ln_equ_posn posn 
)

Calculate a bodies equatorial coords.

Parameters
JDJulian Day.
orbitOrbital parameters.
posnPointer to hold asteroid position.

Calculate a bodies equatorial coordinates for the given julian day.

Examples:
asteroid.c, and comet.c.

References ln_equ_posn::dec, ln_get_ell_helio_rect_posn(), ln_get_solar_geo_coords(), ln_rad_to_deg(), ln_range_degrees(), ln_equ_posn::ra, ln_rect_posn::X, ln_rect_posn::Y, and ln_rect_posn::Z.

Referenced by ln_get_ell_body_next_rst_horizon(), ln_get_ell_body_next_rst_horizon_future(), and ln_get_ell_body_rst_horizon().

double ln_get_ell_body_next_rst ( double  JD,
struct ln_lnlat_posn observer,
struct ln_ell_orbit orbit,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with an elliptic orbit.

Parameters
JDJulian day
observerObservers position
orbitOrbital parameters
rstPointer to store Rise, Set and Transit time in JD
Returns
0 for success, else 1 for circumpolar (above the horizon), -1 for circumpolar (bellow the horizon)

Calculate the time of next rise, set and transit (crosses the local meridian at upper culmination) time of a body with an elliptic orbit for the given Julian day.

This function guarantee, that rise, set and transit will be in <JD, JD+1> range.

Note: this functions returns 1 if the body is circumpolar, that is it remains the whole day above the horizon. Returns -1 when it remains the whole day below the horizon.

References ln_get_ell_body_next_rst_horizon().

double ln_get_ell_body_next_rst_horizon ( double  JD,
struct ln_lnlat_posn observer,
struct ln_ell_orbit orbit,
double  horizon,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with an elliptic orbit.

Parameters
JDJulian day
observerObservers position
orbitOrbital parameters
horizonHorizon height
rstPointer to store Rise, Set and Transit time in JD
Returns
0 for success, else 1 for circumpolar (above the horizon), -1 for circumpolar (bellow the horizon)

Calculate the time of next rise, set and transit (crosses the local meridian at upper culmination) time of a body with an elliptic orbit for the given Julian day.

This function guarantee, that rise, set and transit will be in <JD, JD+1> range.

Note: this functions returns 1 if the body is circumpolar, that is it remains the whole day above the horizon. Returns -1 when it remains the whole day below the horizon.

Parameters
JDJulian day
observerObservers position
orbitOrbital parameters
horizonHorizon height
day_limitMaximal number of days that will be searched for next rise and set
rstPointer to store Rise, Set and Transit time in JD
Returns
0 for success, else 1 for circumpolar (above the horizon), -1 for circumpolar (bellow the horizon)

Calculate the time of next rise, set and transit (crosses the local meridian at upper culmination) time of a body with an elliptic orbit for the given Julian day.

This function guarantee, that rise, set and transit will be in <JD, JD + day_limit> range.

Note: this functions returns 1 if the body is circumpolar, that is it remains the whole day above the horizon. Returns -1 when it remains the whole day below the horizon.

References ln_get_ell_body_equ_coords().

Referenced by ln_get_ell_body_next_rst().

double ln_get_ell_body_phase_angle ( double  JD,
struct ln_ell_orbit orbit 
)

Calculate the phase angle of the body. The angle Sun - body - Earth.

Parameters
JDJulian day
orbitOrbital parameters
Returns
Phase angle.

Calculate the phase angle of the body. The angle Sun - body - Earth.

Examples:
asteroid.c.

References ln_ell_orbit::a, ln_ell_orbit::e, ln_ell_orbit::JD, ln_deg_to_rad(), ln_get_ell_body_earth_dist(), ln_get_ell_body_solar_dist(), ln_get_ell_mean_anomaly(), ln_get_ell_mean_motion(), ln_get_ell_radius_vector(), ln_range_degrees(), ln_solve_kepler(), and ln_ell_orbit::n.

double ln_get_ell_body_rst ( double  JD,
struct ln_lnlat_posn observer,
struct ln_ell_orbit orbit,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with an elliptic orbit.

Parameters
JDJulian day
observerObservers position
orbitOrbital parameters
rstPointer to store Rise, Set and Transit time in JD
Returns
0 for success, else 1 for circumpolar (above the horizon), -1 for circumpolar (bellow the horizon)

Calculate the time the rise, set and transit (crosses the local meridian at upper culmination) time of a body with an elliptic orbit for the given Julian day.

Note: this functions returns 1 if the body is circumpolar, that is it remains the whole day above the horizon. Returns -1 when it remains the whole day below the horizon.

Examples:
asteroid.c, and comet.c.

References ln_get_ell_body_rst_horizon().

double ln_get_ell_body_rst_horizon ( double  JD,
struct ln_lnlat_posn observer,
struct ln_ell_orbit orbit,
double  horizon,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with an elliptic orbit.

Parameters
JDJulian day
observerObservers position
orbitOrbital parameters
horizonHorizon height
rstPointer to store Rise, Set and Transit time in JD
Returns
0 for success, else 1 for circumpolar (above the horizon), -1 for circumpolar (bellow the horizon)

Calculate the time the rise, set and transit (crosses the local meridian at upper culmination) time of a body with an elliptic orbit for the given Julian day.

Note: this functions returns 1 if the body is circumpolar, that is it remains the whole day above the horizon. Returns -1 when it remains the whole day below the horizon.

References ln_get_ell_body_equ_coords().

Referenced by ln_get_ell_body_rst().

double ln_get_ell_body_solar_dist ( double  JD,
struct ln_ell_orbit orbit 
)

Calculate the distance between a body and the Sun.

Parameters
JDJulian Day.
orbitOrbital parameters
Returns
The distance in AU between the Sun and the body.

Calculate the distance between a body and the Sun.

Examples:
asteroid.c, and comet.c.

References ln_get_ell_helio_rect_posn(), ln_rect_posn::X, ln_rect_posn::Y, and ln_rect_posn::Z.

Referenced by ln_get_ell_body_elong(), ln_get_ell_body_phase_angle(), ln_get_ell_comet_mag(), and ln_get_ell_orbit_vel().

void ln_get_ell_geo_rect_posn ( struct ln_ell_orbit orbit,
double  JD,
struct ln_rect_posn posn 
)

Calculate the objects rectangular geocentric position.

Parameters
orbitOrbital parameters of object.
JDJulian day
posnPosition pointer to store objects position

Calculate the objects rectangular geocentric position given it's orbital elements for the given julian day.

Examples:
asteroid.c, and comet.c.

References ln_get_earth_helio_coords(), ln_get_ell_helio_rect_posn(), ln_get_rect_from_helio(), ln_rect_posn::X, ln_rect_posn::Y, and ln_rect_posn::Z.

Referenced by ln_get_ell_body_earth_dist().

void ln_get_ell_helio_rect_posn ( struct ln_ell_orbit orbit,
double  JD,
struct ln_rect_posn posn 
)

Calculate the objects rectangular heliocentric position.

Parameters
orbitOrbital parameters of object.
JDJulian day
posnPosition pointer to store objects position

Calculate the objects rectangular heliocentric position given it's orbital elements for the given julian day.

Examples:
asteroid.c, and comet.c.

References ln_ell_orbit::a, ln_helio_posn::B, ln_ell_orbit::e, ln_ell_orbit::i, ln_ell_orbit::JD, ln_deg_to_rad(), ln_get_ell_mean_anomaly(), ln_get_ell_mean_motion(), ln_get_ell_radius_vector(), ln_get_ell_true_anomaly(), ln_solve_kepler(), ln_ell_orbit::n, ln_ell_orbit::omega, ln_helio_posn::R, ln_ell_orbit::w, ln_rect_posn::X, ln_rect_posn::Y, and ln_rect_posn::Z.

Referenced by ln_get_ell_body_equ_coords(), ln_get_ell_body_solar_dist(), and ln_get_ell_geo_rect_posn().

double ln_get_ell_last_perihelion ( double  epoch_JD,
double  M,
double  n 
)

Calculate the julian day of the last perihelion.

Parameters
epoch_JDJulian day of epoch
MMean anomaly
ndaily motion in degrees

Calculate the julian day of the last perihelion.

Examples:
asteroid.c.
double ln_get_ell_mean_anomaly ( double  n,
double  delta_JD 
)

Calculate the mean anomaly.

Parameters
nMean motion (degrees/day)
delta_JDTime since perihelion
Returns
Mean anomaly (degrees)

Calculate the mean anomaly.

Examples:
comet.c.

Referenced by ln_get_ell_body_elong(), ln_get_ell_body_phase_angle(), ln_get_ell_comet_mag(), and ln_get_ell_helio_rect_posn().

double ln_get_ell_mean_motion ( double  a)

Calculate the mean daily motion (degrees/day).

Parameters
aSemi major diameter in AU
Returns
Mean daily motion (degrees/day)

Calculate the mean daily motion (degrees/day).

Examples:
comet.c.

Referenced by ln_get_ell_body_elong(), ln_get_ell_body_phase_angle(), ln_get_ell_comet_mag(), and ln_get_ell_helio_rect_posn().

double ln_get_ell_orbit_avel ( struct ln_ell_orbit orbit)

Calculate the orbital velocity at aphelion in km/s.

Parameters
orbitOrbital parameters
Returns
Orbital velocity in km/s.

Calculate the orbital velocity at aphelion in km/s.

Examples:
asteroid.c, and comet.c.

References ln_ell_orbit::a, and ln_ell_orbit::e.

double ln_get_ell_orbit_len ( struct ln_ell_orbit orbit)

Calculate the orbital length in AU.

Parameters
orbitOrbital parameters
Returns
Orbital length in AU

Calculate the orbital length in AU.

Accuracy:

  • 0.001% for e < 0.88
  • 0.01% for e < 0.95
  • 1% for e = 0.9997
  • 3% for e = 1
Examples:
asteroid.c, and comet.c.

References ln_ell_orbit::a, ln_ell_orbit::e, and ln_get_ell_sminor_diam().

double ln_get_ell_orbit_pvel ( struct ln_ell_orbit orbit)

Calculate orbital velocity at perihelion in km/s.

Parameters
orbitOrbital parameters
Returns
Orbital velocity in km/s.

Calculate orbital velocity at perihelion in km/s.

Examples:
asteroid.c, and comet.c.

References ln_ell_orbit::a, and ln_ell_orbit::e.

double ln_get_ell_orbit_vel ( double  JD,
struct ln_ell_orbit orbit 
)

Calculate orbital velocity in km/s.

Parameters
JDJulian day.
orbitOrbital parameters
Returns
Orbital velocity in km/s.

Calculate orbital velocity in km/s for the given julian day.

Examples:
asteroid.c, and comet.c.

References ln_ell_orbit::a, and ln_get_ell_body_solar_dist().

double ln_get_ell_radius_vector ( double  a,
double  e,
double  E 
)

Calculate the radius vector.

Parameters
aSemi-Major axis in AU
eOrbital eccentricity
EEccentric anomaly
Returns
Radius vector AU

Calculate the radius vector.

Examples:
comet.c.

References ln_deg_to_rad().

Referenced by ln_get_ell_body_elong(), ln_get_ell_body_phase_angle(), ln_get_ell_comet_mag(), and ln_get_ell_helio_rect_posn().

double ln_get_ell_smajor_diam ( double  e,
double  q 
)

Calculate the semi major diameter.

Parameters
eEccentricity
qPerihelion distance in AU
Returns
Semi-major diameter in AU

Calculate the semi major diameter.

double ln_get_ell_sminor_diam ( double  e,
double  a 
)

Calculate the semi minor diameter.

Parameters
eEccentricity
aSemi-Major diameter in AU
Returns
Semi-minor diameter in AU

Calculate the semi minor diameter.

Referenced by ln_get_ell_orbit_len().

double ln_get_ell_true_anomaly ( double  e,
double  E 
)

Calculate the true anomaly.

Parameters
eOrbital eccentricity
EEccentric anomaly
Returns
True anomaly (degrees)

Calculate the true anomaly.

Examples:
comet.c.

References ln_deg_to_rad(), ln_rad_to_deg(), and ln_range_degrees().

Referenced by ln_get_ell_helio_rect_posn().

double ln_solve_kepler ( double  E,
double  M 
)

Calculate the eccentric anomaly.

Parameters
EOrbital eccentricity
MMean anomaly
Returns
Eccentric anomaly

Calculate the eccentric anomaly. This method was devised by Roger Sinnott. (Sky and Telescope, Vol 70, pg 159)

Examples:
comet.c.

References ln_deg_to_rad(), and ln_rad_to_deg().

Referenced by ln_get_ell_body_elong(), ln_get_ell_body_phase_angle(), ln_get_ell_comet_mag(), and ln_get_ell_helio_rect_posn().