libnova  v 0.15.0
Functions
Parabolic Motion

Functions

double LIBNOVA_EXPORT ln_solve_barker (double q, double t)
 Solve Barkers equation.
double LIBNOVA_EXPORT ln_get_par_true_anomaly (double q, double t)
 Calculate the true anomaly.
double LIBNOVA_EXPORT ln_get_par_radius_vector (double q, double t)
 Calculate the radius vector.
void LIBNOVA_EXPORT ln_get_par_geo_rect_posn (struct ln_par_orbit *orbit, double JD, struct ln_rect_posn *posn)
 Calculate an objects rectangular geocentric position.
void LIBNOVA_EXPORT ln_get_par_helio_rect_posn (struct ln_par_orbit *orbit, double JD, struct ln_rect_posn *posn)
 Calculate an objects rectangular heliocentric position.
void LIBNOVA_EXPORT ln_get_par_body_equ_coords (double JD, struct ln_par_orbit *orbit, struct ln_equ_posn *posn)
 Calculate a bodies equatorial coordinates.
double LIBNOVA_EXPORT ln_get_par_body_earth_dist (double JD, struct ln_par_orbit *orbit)
 Calculate the distance between a body and the Earth.
double LIBNOVA_EXPORT ln_get_par_body_solar_dist (double JD, struct ln_par_orbit *orbit)
 Calculate the distance between a body and the Sun.
double LIBNOVA_EXPORT ln_get_par_body_phase_angle (double JD, struct ln_par_orbit *orbit)
 Calculate the phase angle of the body.
double LIBNOVA_EXPORT ln_get_par_body_elong (double JD, struct ln_par_orbit *orbit)
 Calculate the bodies elongation to the Sun.
int LIBNOVA_EXPORT ln_get_par_body_rst (double JD, struct ln_lnlat_posn *observer, struct ln_par_orbit *orbit, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with a parabolic orbit.
int LIBNOVA_EXPORT ln_get_par_body_rst_horizon (double JD, struct ln_lnlat_posn *observer, struct ln_par_orbit *orbit, double horizon, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with a parabolic orbit.
int LIBNOVA_EXPORT ln_get_par_body_next_rst (double JD, struct ln_lnlat_posn *observer, struct ln_par_orbit *orbit, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an parabolic orbit.
int LIBNOVA_EXPORT ln_get_par_body_next_rst_horizon (double JD, struct ln_lnlat_posn *observer, struct ln_par_orbit *orbit, double horizon, struct ln_rst_time *rst)
 Calculate the time of rise, set and transit for a body with an parabolic orbit.
int LIBNOVA_EXPORT ln_get_par_body_next_rst_horizon_future (double JD, struct ln_lnlat_posn *observer, struct ln_par_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 parabolic orbit.

Detailed Description

Functions relating to the Parabolic motion of bodies.

All angles are expressed in degrees.


Function Documentation

double ln_get_par_body_earth_dist ( double  JD,
struct ln_par_orbit orbit 
)

Calculate the distance between a body and the Earth.

Parameters:
JDJulian day.
orbitOrbital parameters
Returns:
Distance in AU

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

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

double ln_get_par_body_elong ( double  JD,
struct ln_par_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..

References ln_par_orbit::JD, ln_get_earth_solar_dist(), ln_get_par_body_solar_dist(), ln_get_par_radius_vector(), ln_rad_to_deg(), ln_range_degrees(), and ln_par_orbit::q.

void ln_get_par_body_equ_coords ( double  JD,
struct ln_par_orbit orbit,
struct ln_equ_posn posn 
)

Calculate a bodies equatorial coordinates.

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

Calculate a bodies equatorial coordinates for the given julian day.

References ln_equ_posn::dec, ln_get_par_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_par_body_next_rst_horizon(), ln_get_par_body_next_rst_horizon_future(), and ln_get_par_body_rst_horizon().

double ln_get_par_body_next_rst ( double  JD,
struct ln_lnlat_posn observer,
struct ln_par_orbit orbit,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with an parabolic 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 parabolic 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_par_body_next_rst_horizon().

double ln_get_par_body_next_rst_horizon ( double  JD,
struct ln_lnlat_posn observer,
struct ln_par_orbit orbit,
double  horizon,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with an parabolic 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 parabolic 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_par_body_equ_coords().

Referenced by ln_get_par_body_next_rst().

double ln_get_par_body_next_rst_horizon_future ( double  JD,
struct ln_lnlat_posn observer,
struct ln_par_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 parabolic orbit.

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 parabolic 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_par_body_equ_coords().

double ln_get_par_body_phase_angle ( double  JD,
struct ln_par_orbit orbit 
)

Calculate the phase angle of the body.

Parameters:
JDJulian day
orbitOrbital parameters
Returns:
Phase angle.

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

References ln_par_orbit::JD, ln_get_earth_solar_dist(), ln_get_par_body_solar_dist(), ln_get_par_radius_vector(), ln_rad_to_deg(), ln_range_degrees(), and ln_par_orbit::q.

double ln_get_par_body_rst ( double  JD,
struct ln_lnlat_posn observer,
struct ln_par_orbit orbit,
struct ln_rst_time rst 
)

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

Parameters:
JDJulian day
observerObservers position
orbitOrbital parameters
rstPointer to store Rise, Set and Transit time in JD
Returns:
0 for success, 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 a parabolic orbit for the given Julian day.

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

References ln_get_par_body_rst_horizon().

double ln_get_par_body_rst_horizon ( double  JD,
struct ln_lnlat_posn observer,
struct ln_par_orbit orbit,
double  horizon,
struct ln_rst_time rst 
)

Calculate the time of rise, set and transit for a body with a parabolic 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.

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

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

References ln_get_par_body_equ_coords().

Referenced by ln_get_par_body_rst().

double ln_get_par_body_solar_dist ( double  JD,
struct ln_par_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.

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

Referenced by ln_get_par_body_elong(), ln_get_par_body_phase_angle(), and ln_get_par_comet_mag().

void ln_get_par_geo_rect_posn ( struct ln_par_orbit orbit,
double  JD,
struct ln_rect_posn posn 
)

Calculate an 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.

References ln_get_earth_helio_coords(), ln_get_par_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_par_body_earth_dist().

void ln_get_par_helio_rect_posn ( struct ln_par_orbit orbit,
double  JD,
struct ln_rect_posn posn 
)

Calculate an 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.

References ln_par_orbit::i, ln_par_orbit::JD, ln_deg_to_rad(), ln_get_par_radius_vector(), ln_get_par_true_anomaly(), ln_par_orbit::omega, ln_par_orbit::q, ln_par_orbit::w, ln_rect_posn::X, ln_rect_posn::Y, and ln_rect_posn::Z.

Referenced by ln_get_par_body_equ_coords(), ln_get_par_body_solar_dist(), and ln_get_par_geo_rect_posn().

double ln_get_par_radius_vector ( double  q,
double  t 
)

Calculate the radius vector.

Parameters:
qPerihelion distance in AU
tTime since perihelion in days
Returns:
Radius vector AU

Calculate the radius vector.

References ln_solve_barker().

Referenced by ln_get_par_body_elong(), ln_get_par_body_phase_angle(), ln_get_par_comet_mag(), and ln_get_par_helio_rect_posn().

double ln_get_par_true_anomaly ( double  q,
double  t 
)

Calculate the true anomaly.

Parameters:
qPerihelion distance in AU
tTime since perihelion
Returns:
True anomaly (degrees)

Calculate the true anomaly.

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

Referenced by ln_get_par_helio_rect_posn().

double ln_solve_barker ( double  q,
double  t 
)

Solve Barkers equation.

Parameters:
qPerihelion distance in AU
tTime since perihelion in days
Returns:
Solution of Barkers equation

Solve Barkers equation. LIAM add more

Referenced by ln_get_par_radius_vector(), and ln_get_par_true_anomaly().