Data Structures | Macros | Enumerations | Functions | Variables
aephem.h File Reference

Header file for AEPHEM library. More...

#include <math.h>
#include <stdio.h>

Go to the source code of this file.

Data Structures

struct  ae_plantbl_t
 For holding data on planetary ephemerides. More...
 
struct  ae_orbit_t
 For holding information on an orbit. More...
 
struct  ae_star_t
 For holding information on a star/planet. More...
 
struct  ae_jpl_handle_t
 For holding information on a JPL ephemeris file. More...
 
struct  ae_physical_term_t
 The coefficients for sine/cosine terms in a body's rotation data. More...
 
struct  ae_physical_t
 For holding information on a body's physical and rotational properties. More...
 
struct  ae_ring_geometry_t
 For storing information on a planet's (i.e., Saturn's) rings. More...
 
struct  ae_stokes_t
 For holding Stokes parameters. More...
 
struct  ae_saturn_temp_model_t
 For modelling Saturn's emission. More...
 
struct  ae_saturn_components_t
 
struct  ae_planet_map_t
 For storing a planetary map. More...
 

Macros

#define AE_DTR   1.7453292519943295769e-2
 Degrees to radians.
 
#define AE_RTD   5.7295779513082320877e1
 Radians to degrees.
 
#define AE_RTS   2.0626480624709635516e5
 Arcseconds per radian.
 
#define AE_STR   4.8481368110953599359e-6
 Radians per arcsecond.
 
#define AE_STD   (1.0 / 3600.0)
 Degrees per arcsecond.
 
#define AE_STDAY   0.00001157407407407407
 Days per second.
 
#define AE_COSSUN   -0.014543897651582657
 Cos(90d 50').
 
#define AE_COSZEN   -9.8900378587411476e-3
 Cos(90d 34').
 
#define AE_J2000   2451545.0
 1.5 Jan 2000.
 
#define AE_B1950   2433282.423
 0.923 Jan 1950, Besselian.
 
#define AE_J1900   2415020.0
 0 Jan 1900. More...
 
#define AE_MJD_START   2400000.5
 The Julian Date when the Modified Julian Date begins.
 
#define AE_AU   1.49597870691e8
 Astronomical unit in km.
 
#define AE_CLIGHT   2.99792458e5
 Speed of light in km/s.
 
#define AE_S_PER_D   86400.0
 Number of seconds per day.
 
#define AE_D_PER_S   (1.0 / AE_S_PER_D)
 Number of days per second. More...
 
#define AE_CLIGHT_AUD   (86400.0 * AE_CLIGHT / AE_AU)
 The speed of light in AU / day.
 
#define AE_R_EARTH   6378137.0
 Radius of earth in m. More...
 
#define AE_R_EARTH_AU   (0.001 * AE_R_EARTH / AE_AU)
 Radius of the earth in AU. More...
 
#define AE_E_M_RAT   81.300585
 Earth/moon mass ratio. More...
 
#define AE_FLAT   298.257222101
 The reciprocal of flattening for the earth. More...
 
#define AE_T_CMB   2.725
 The temperature of the cosmic background radiation, in degrees kelvin.
 
#define AE_JPL_HEAD_SIZE   (5 * sizeof(double) + 41 * sizeof(int32_t))
 Size of the JPL binary header.
 
#define AE_JPL_CONST_LEN   6
 Constant name lengths in JPL ephemerides.
 
#define AE_N_SS_BODIES   16
 Number of entries in ae_ss_bodies_t.
 
#define AE_N_SS_BODIES_JPL   13
 Number of bodies given in the JPL ephemeris.
 
#define AE_PLANTBL_N_HARMONIC   18
 Number of arguments for ae_plantbl_t. More...
 
#define ae_get_d(a)   ((double)((int)(a)))
 Extract the degrees/hours from a decimal angle/time. More...
 
#define ae_get_m(a)   ((double)fabs(((int)(((a) - ae_get_d(a)) * 60))))
 Extract the minutes from a decimal angle/time. More...
 
#define ae_get_s(a)   (fabs(((a) - ae_get_d(a) - ((double)((int)(((a) - ae_get_d(a)) * 60))) / 60) * 3600))
 Extract the seconds from a decimal angle/time. More...
 
#define ae_dms_fmt   "%dd %d' %.2f\""
 Defines a format string for printing D:M:S. More...
 
#define ae_dms_sfmt   "%4dd %02d' %05.2f\""
 Defines a format string for printing D:M:S with fixed spacing. More...
 
#define ae_dms_arg(a)   (int)ae_get_d(a), (int)ae_get_m(a), ae_get_s(a)
 Defines the arguments for printing D:M:S. More...
 
#define ae_hms_fmt   "%dh %dm %.2fs"
 Defines a format string for printing H:M:S. More...
 
#define ae_hms_sfmt   "%3dh %02dm %05.2fs"
 Defines a format string for printing H:M:S with fixed spacing. More...
 
#define ae_hms_arg(a)   (int)ae_get_d(a / 15.0), (int)ae_get_m(a / 15.0), ae_get_s(a / 15.0)
 Defines the arguments for printing D:M:S. More...
 

Enumerations

enum  ae_retcode_t {
  AE_RET_BAD_PATH = 1, AE_RET_UNEXPECTED_EOF, AE_RET_READ_ERROR, AE_RET_BAD_JPL_HEADER,
  AE_RET_BAD_JPL_HEADER_KSIZE, AE_RET_BAD_JPL_HEADER_GROUP, AE_RET_BAD_JPL_HEADER_NUM, AE_RET_BAD_JPL_HEADER_RANGE,
  AE_RET_BAD_JPL_NOT_ASCII, AE_RET_BAD_JPL_NOT_BIN, AE_RET_BAD_JPL_BLOCK, AE_RET_BAD_JPL_CORRUPT,
  AE_RET_BAD_JPL_DATE, AE_RET_BAD_JPL_FLOAT, AE_RET_BAD_JPL_CHEBY_TIME, AE_RET_BAD_OBJ_INDEX,
  AE_RET_OBJ_NOT_FOUND, AE_RET_CFITSIO_NOT_ENABLED, AE_RET_CFITSIO_NO_KEYWORD, AE_RET_CFITSIO_BAD_KEYWORD,
  AE_RET_WCS_UNKNOWN_PROJ, AE_RET_CFITSIO_ERROR, AE_RET_BAD_NATIVE_POLE, AE_RET_NO_INVERSE,
  AE_RET_UNDER_CONSTRUCTION
}
 Error codes that functions might return. More...
 
enum  ae_ss_bodies_t {
  AE_SS_MERCURY = 0, AE_SS_VENUS, AE_SS_EMBARY, AE_SS_MARS,
  AE_SS_JUPITER, AE_SS_SATURN, AE_SS_URANUS, AE_SS_NEPTUNE,
  AE_SS_PLUTO, AE_SS_MOON_EMB, AE_SS_SUN, AE_SS_NUTATION,
  AE_SS_LIBRATION, AE_SS_EARTH, AE_SS_MOON, AE_SS_SSBARY
}
 Enumerate Solar System bodies. More...
 
enum  ae_precess_direction_t { AE_FROM_J2000 = -1, AE_TO_J2000 = 1 }
 Define directions for precession and nutation. More...
 
enum  ae_aberration_direction_t { AE_ADD_ABERRATION = -1, AE_REMOVE_ABERRATION = 1 }
 Define directions for aberration. More...
 
enum  ae_physical_use_d_or_t { AE_PHYSICAL_END = 0, AE_PHYSICAL_D, AE_PHYSICAL_T }
 Define whether a term uses T or d for physical ephemeride terms. More...
 

Functions

void ae_altaz_to_radec (double last, double glat, double alt, double az, double *ra, double *dec)
 Convert alt/az to ra/dec. More...
 
void ae_annual_aberration (double v_earth[], double p[], int direction)
 Correct for the annual aberration. More...
 
double ae_cal_to_jd (long year, int month, double day)
 Calculate Julian date from a Gregorian calendar date. More...
 
int ae_cat_to_constel_index (const char *cat_name, char *full_name, int len)
 Given a catalogue star name, find a star's constellation. More...
 
int ae_coord_to_constel_index (double ra, double dec, double epoch)
 Given a sky coordinate, find the constellation it is in. More...
 
double ae_ctime_to_jd (double t)
 Convert a C time to a Julian date. More...
 
double ae_ctime_to_last (double t, double delta_t, double tlong, double nutl, double eps)
 Calculate the local apparent sidereal time from a C time. More...
 
double aes_ctime_to_last (double t, double tlong)
 Calculate the local sidereal time from a C time. More...
 
double ae_delta_t (double jd_ut1)
 Get the value of delta T == Ephemeris Time - Universal Time. More...
 
void ae_delta_q (double q0[], double q1[], double *dra, double *ddec)
 Convert change in rectangular coordinatates to change in ra/dec. More...
 
double ae_disc_semiminor (double a, double c, double lat)
 Calculate the semi-minor axis of an observed disc. More...
 
double ae_disc_solid_angle (double a, double b, double dist)
 Calculate the solid angle of an observed disc. More...
 
double aes_disc_solid_angle (double jd_ut1, const struct ae_orbit_t *o_orb, const struct ae_orbit_t *q_orb, const struct ae_physical_t *phys)
 Calculate the solid angle of an observed disc. More...
 
void ae_diurnal_aberration (double last, double tlat, double trho, int direction, double *ra, double *dec)
 Correct for diurnal aberration. More...
 
void ae_diurnal_parallax (double last, double tlat, double trho, double dist, double *ra, double *dec)
 Correct for diurnal parallax. More...
 
double ae_dut1 (double jd_utc)
 Given a UTC, return DUT1. More...
 
double ae_epsilon (double jd_tt)
 Calculate the obliquity of the ecliptic. More...
 
void ae_fk4_to_fk5 (struct ae_star_t *star)
 Convert FK4 B1950.0 catalogue coordinates to FK5 J2000.0 coordinates. More...
 
void ae_geocentric (double jd_tt, double q[], double e[], double v_e[], double *ra, double *dec, double *dist)
 Given heliocentric coordinates, reduce to geocentric coordinates. More...
 
void ae_geocentric_moon_from_orbit (double jd_tt, const struct ae_orbit_t *o_orb, double *ra, double *dec, double *dist)
 Calculate the geocentric coordinates of the moon from orbital elements. More...
 
void ae_geocentric_from_cat (double jd_tt, double e[], double v_e[], const struct ae_star_t *star, double *ra, double *dec)
 Reduce a star's catalogue coordinates to its apparent place. More...
 
int ae_geocentric_from_jpl (struct ae_jpl_handle_t *jh, double jd_tt, int obj_num, double *ra, double *dec, double *dist)
 Determine geocentric coordinates using JPL ephemerides. More...
 
void ae_geocentric_from_orbit (double jd_tt, const struct ae_orbit_t *o_orb, const struct ae_orbit_t *q_orb, double *ra, double *dec, double *dist)
 Determine geocentric coordinates using orbital elements. More...
 
void ae_geocentric_sun_from_orbit (double jd_tt, const struct ae_orbit_t *o_orb, double *ra, double *dec, double *dist)
 Calculate the geocentric coordinates of the sun from orbital elements. More...
 
void ae_geocentric_lat (double glat, double height, double *tlat, double *trho)
 Calculate the geocentric latitude and the distance to the earth's centre. More...
 
int ae_read_orbit_from_cat (const char *path, const char *name, struct ae_orbit_t *orb)
 Parse orbital elements from a catalogue file. More...
 
int ae_read_star_from_cat (const char *path, const char *name, struct ae_star_t *star)
 Parse star positional data from a catalogue file. More...
 
void ae_gmoon (double jd, double rect[], double pol[])
 Compute geocentric moon position. More...
 
double ae_gmst (double jd_ut1, double jd_tt)
 Get the Greenwich mean sidereal time. More...
 
double ae_gast (double jd_ut1, double jd_tt, double nutl, double eps)
 Get Greenwich apparent sidereal time. More...
 
double ae_lmst (double jd_ut1, double jd_tt, double tlong)
 Get the local mean sidereal time. More...
 
double ae_last (double jd_ut1, double jd_tt, double tlong, double nutl, double eps)
 Get the local apparent sidereal time. More...
 
double aes_last (double jd_ut1, double tlong)
 Get the local apparent sidereal time. More...
 
void ae_gplan (double jd_tt, struct ae_plantbl_t *plan, double pobj[])
 Use a table to find a planet's polar, heliocentric position. More...
 
double ae_g1plan (double jd_tt, struct ae_plantbl_t *plan)
 Accumulate the sum of trigonometric series in one variable. More...
 
void ae_g2plan (double jd_tt, struct ae_plantbl_t *plan, double pobj[], double *lp_equinox)
 Accumulate the sum of trigonometric series in two variables. More...
 
void ae_g3plan (double jd_tt, struct ae_plantbl_t *plan, double pobj[], int objnum)
 Accumulate the sum of trigonometric series in three variables. More...
 
void ae_jd_to_cal (double jd, int *year, int *month, double *day)
 Calculate month, day, and year from a Julian date. More...
 
double ae_jpl_get_const_val (const struct ae_jpl_handle_t *jh, const char *const_name)
 Get the value of a constant defined in the JPL ephemeris header. More...
 
int ae_jpl_init (const char *path, struct ae_jpl_handle_t *jh)
 Initialise to use a JPL ephemeris file. More...
 
int ae_jpl_init_ascii (const char *path, const char *header_path, int search_dates, struct ae_jpl_handle_t *jh)
 Read a JPL ephemeris ASCII header. More...
 
int ae_jpl_init_bin (const char *path, struct ae_jpl_handle_t *jh)
 Read a JPL Ephemeris header in binary format. More...
 
int ae_jpl_get_coords (struct ae_jpl_handle_t *jh, double jd_tt, int obj_num, double r[], double v[], char is_planetary)
 Get object positions from a JPL ephemeris file. More...
 
void ae_jpl_close (struct ae_jpl_handle_t *jh)
 Finish with a JPL handle. More...
 
void ae_kepler (double jd_tt, const struct ae_orbit_t *orb, double q[])
 Find a Keplerian orbit. More...
 
double ae_light_t (double p[], double q[], int do_retardation)
 Get light-travel time between two objects. More...
 
int ae_light_t_jpl (struct ae_jpl_handle_t *jh, double jd_tt, double e[], int obj_num, double q[], double v_q[], char is_planetary)
 Get the apparent position of an object, corrected for light-travel time. More...
 
void ae_light_t_orbit (double jd_tt, double e[], const struct ae_orbit_t *orb, double q[])
 Get the apparent position of an object, corrected for light-travel time. More...
 
double ae_mjd (double jd)
 Given a Julian Date, get the Modified Julian Date. More...
 
double ae_mod_2pi (double x)
 Reduce x modulo 2 pi. More...
 
double ae_mod_360 (double x)
 Reduce x modulo 360 degrees. More...
 
double ae_mod_180 (double x)
 Reduce x modulo 360 degrees, but centre at 0 degrees. More...
 
void ae_nutation_lon_ob (double jd_tt, double *nutl, double *nuto)
 Calculate the nutation in longitude and oblation. More...
 
void ae_nutate (double nutl, double nuto, double eps, double p[], int direction)
 Correct for nutation. More...
 
void aes_nutate (double jd_tt, double p[], int direction)
 Correct for nutation. More...
 
void ae_phys_pole (double jd_ut1, double jd_tt, const struct ae_physical_t *phys, double n[], double *w)
 Calculate the axis of a body's north pole and its prime meridean. More...
 
int ae_is_retrograde (const struct ae_physical_t *phys)
 Determine whether a body's rotation is retrograde. More...
 
void ae_polar_to_rect (double ra, double dec, double radius, double rect[])
 Convert ra/dec/radius to a rectangular vector. More...
 
void ae_precess (double jd_tt, double r[], int direction)
 Precess a coordinate. More...
 
void ae_radec_to_altaz (double last, double glat, double ra, double dec, double *alt, double *az)
 Convert ra/dec to alt/az. More...
 
void ae_radec_to_gal (double ra, double dec, double *l, double *b, int fk5)
 Convert ra/dec to l/b (galactic coordinates). More...
 
void ae_gal_to_radec (double ra, double dec, double *l, double *b, int fk5)
 Convert l/b (galactic coordinates) to ra/dec. More...
 
void ae_rect_to_polar (const double rect[], double *ra, double *dec, double *radius)
 Convert an equatorial rectangular unit vector to ra/dec/radius. More...
 
void ae_rect_to_polar_with_jd (double pp[], double jd, char ofdate, double polar[])
 Convert equatorial rectangular coordinates to polar coordinates. More...
 
double ae_refrac_visible (double alt, void *param)
 Atmospheric refraction correction (in visible bands). More...
 
double ae_refrac_ulich (double alt, void *param)
 Atmospheric refraction correction (in millimetre bands). More...
 
void ae_relativity (double p[], double q[], double o[])
 Correct for light deflection due to solar gravitation. More...
 
const char * ae_retcode_string (int retcode)
 For returning human readable return codes. More...
 
double ae_tdb (double jd)
 Find Barycentric Dynamical Time from Terrestrial Dynamical Time. More...
 
void ae_subobs_point (double j[], double n[], double w, double f, int retrograde, double *lat, double *lon)
 Calculate the sub-observer point. More...
 
void aes_subobs_point (double jd_ut1, const struct ae_orbit_t *o_orb, const struct ae_orbit_t *q_orb, const struct ae_physical_t *phys, double *lat, double *lon, double *dist)
 Calculate the sub-observer point. More...
 
double ae_flattening (const struct ae_physical_t *phys)
 Calculate a body's flattening. More...
 
void ae_topocentric (double jd_tt, double jd_ut1, double tlat, double glat, double tlong, double trho, double(*refrac)(double, void *), void *refrac_param, double dist, double *ra, double *dec)
 Apply diurnal aberrations and calculate topocentric altitude and azimuth. More...
 
void aes_topocentric (double jd_ut1, double glat, double tlong, double dist, double *ra, double *dec)
 Apply diurnal aberrations and calculate topocentric altitude and azimuth. More...
 
void ae_v_orbit (double jd_tt, const struct ae_orbit_t *orb, double v[])
 Calculate the velocity of an object as it moves in its orbit. More...
 
double ae_zatan2 (double x, double y)
 Quadrant correct inverse circular tangent. More...
 
void ae_map_saturn (double res, const struct ae_physical_t *b, const struct ae_ring_geometry_t *r, double inclination, const struct ae_saturn_temp_model_t *sm, struct ae_planet_map_t *map, struct ae_saturn_components_t *comp)
 
int ae_map_to_gnuplot (const struct ae_planet_map_t *map, char stokes, const char *path)
 
void ae_map_t_eff (const struct ae_planet_map_t *map, const struct ae_physical_t *b, double subobs_lat, double background, struct ae_stokes_t *t_eff)
 
void ae_saturn_t_eff (const struct ae_ring_geometry_t *r, const struct ae_saturn_temp_model_t *sm, const struct ae_saturn_components_t *components, struct ae_stokes_t *t_eff)
 
void ae_set_stokes (struct ae_stokes_t *s, double i, double q, double u, double v)
 
void ae_copy_stokes (const struct ae_stokes_t *origin, struct ae_stokes_t *dest)
 
void ae_saturn_components_alloc (int n_ring, struct ae_saturn_components_t *comp)
 
void ae_saturn_components_free (struct ae_saturn_components_t *comp)
 

Variables

struct ae_orbit_t ae_orb_mercury
 Default orbital parameters for Mercury.
 
struct ae_orbit_t ae_orb_venus
 Default orbital parameters for Venus.
 
struct ae_orbit_t ae_orb_earth
 Default orbital parameters for Earth.
 
struct ae_orbit_t ae_orb_mars
 Default orbital parameters for Mars.
 
struct ae_orbit_t ae_orb_jupiter
 Default orbital parameters for Jupiter.
 
struct ae_orbit_t ae_orb_saturn
 Default orbital parameters for Saturn.
 
struct ae_orbit_t ae_orb_uranus
 Default orbital parameters for Uranus.
 
struct ae_orbit_t ae_orb_neptune
 Default orbital parameters for Neptune.
 
struct ae_orbit_t ae_orb_pluto
 Default orbital parameters for Pluto.
 
struct ae_orbit_tae_orb_planet []
 Pointers to the planetary orbit objects. More...
 
struct ae_plantbl_t ae_mer404
 Ephemeris for Mercury.
 
struct ae_plantbl_t ae_ven404
 Ephemeris for Venus.
 
struct ae_plantbl_t ae_ear404
 Ephemeris for the earth-moon barycentre.
 
struct ae_plantbl_t ae_mar404
 Ephemeris for Mars.
 
struct ae_plantbl_t ae_jup404
 Ephemeris for Jupiter.
 
struct ae_plantbl_t ae_sat404
 Ephemeris for Saturn.
 
struct ae_plantbl_t ae_ura404
 Ephemeris for Uranus.
 
struct ae_plantbl_t ae_nep404
 Ephemeris for Neptune.
 
struct ae_plantbl_t ae_plu404
 Ephemeris for Pluto.
 
struct ae_plantbl_t ae_mlat404
 Ephemeris for the moon. More...
 
struct ae_plantbl_t ae_mlr404
 Ephemeris data for the moon. More...
 
struct ae_physical_t ae_phys_sun
 Physical parameters for the Sun.
 
struct ae_physical_t ae_phys_mercury
 Physical parameters for Mercury.
 
struct ae_physical_t ae_phys_venus
 Physical parameters for Venus.
 
struct ae_physical_t ae_phys_earth
 Physical parameters for Earth.
 
struct ae_physical_t ae_phys_mars
 Physical parameters for Mars.
 
struct ae_physical_t ae_phys_jupiter
 Physical parameters for Jupiter.
 
struct ae_physical_t ae_phys_saturn
 Physical parameters for Saturn.
 
struct ae_ring_geometry_t ae_saturn_rings
 Physical parameters for Saturn's rings.
 
struct ae_physical_t ae_phys_uranus
 Physical parameters for Uranus.
 
struct ae_physical_t ae_phys_neptune
 Physical parameters for Neptune.
 
struct ae_physical_t ae_phys_pluto
 Physical parameters for Pluto.
 
struct ae_physical_tae_phys_planet []
 Pointers to the planetary physical information objects. More...
 
struct ae_physical_term_t ae_phys_no_term []
 An empty series for ephemeride terms with no extra terms.
 
const char * ae_constel_name []
 A list of constellation names. More...
 
const char * ae_ss_name []
 The names of the solar system bodies. More...
 

Detailed Description

Header file for AEPHEM library.

Macro Definition Documentation

#define AE_D_PER_S   (1.0 / AE_S_PER_D)

Number of days per second.

Referenced by ae_ctime_to_last(), aes_ctime_to_last(), aes_last(), and aes_subobs_point().

#define ae_dms_arg (   a)    (int)ae_get_d(a), (int)ae_get_m(a), ae_get_s(a)

Defines the arguments for printing D:M:S.

Use in conjunction with ae_dms_fmt.

#define ae_dms_fmt   "%dd %d' %.2f\""

Defines a format string for printing D:M:S.

Use in conjuction with ae_dms_arg.

#define ae_dms_sfmt   "%4dd %02d' %05.2f\""

Defines a format string for printing D:M:S with fixed spacing.

Use in conjuction with ae_dms_arg.

#define AE_E_M_RAT   81.300585

Earth/moon mass ratio.

Referenced by ae_jpl_init_ascii(), and ae_kepler_embofs().

#define AE_FLAT   298.257222101

The reciprocal of flattening for the earth.

Taken from the World Geodetic System 1984 (WGS84) definitions.

Referenced by ae_geocentric_lat().

#define ae_get_d (   a)    ((double)((int)(a)))

Extract the degrees/hours from a decimal angle/time.

Parameters
aThe decimal angle.
Returns
The degrees/hours.
#define ae_get_m (   a)    ((double)fabs(((int)(((a) - ae_get_d(a)) * 60))))

Extract the minutes from a decimal angle/time.

Parameters
aThe decimal angle.
Returns
The minutes.
#define ae_get_s (   a)    (fabs(((a) - ae_get_d(a) - ((double)((int)(((a) - ae_get_d(a)) * 60))) / 60) * 3600))

Extract the seconds from a decimal angle/time.

Parameters
aThe decimal angle.
Returns
The seconds.
#define ae_hms_arg (   a)    (int)ae_get_d(a / 15.0), (int)ae_get_m(a / 15.0), ae_get_s(a / 15.0)

Defines the arguments for printing D:M:S.

Use in conjunction with ae_dms_fmt.

#define ae_hms_fmt   "%dh %dm %.2fs"

Defines a format string for printing H:M:S.

Use in conjuction with ae_hms_arg.

#define ae_hms_sfmt   "%3dh %02dm %05.2fs"

Defines a format string for printing H:M:S with fixed spacing.

Use in conjuction with ae_hms_arg.

#define AE_J1900   2415020.0

0 Jan 1900.

Referenced by ae_read_star_from_cat().

#define AE_PLANTBL_N_HARMONIC   18

Number of arguments for ae_plantbl_t.

Referenced by ae_g1plan(), ae_g2plan(), and ae_g3plan().

#define AE_R_EARTH   6378137.0

Radius of earth in m.

Referenced by ae_diurnal_parallax(), and ae_geocentric_lat().

#define AE_R_EARTH_AU   (0.001 * AE_R_EARTH / AE_AU)

Radius of the earth in AU.

Radius of earth in AU.

Referenced by ae_geocentric_moon_from_orbit().

Enumeration Type Documentation

Define directions for aberration.

Enumerator
AE_ADD_ABERRATION 

Add the effect of aberration.

AE_REMOVE_ABERRATION 

Remove the effect of aberration.

Define whether a term uses T or d for physical ephemeride terms.

Enumerator
AE_PHYSICAL_END 

For marking the end of the series.

AE_PHYSICAL_D 

Use d = # days since J2000.

AE_PHYSICAL_T 

Use T = # Julian centuries since J2000.

Define directions for precession and nutation.

Enumerator
AE_FROM_J2000 

Transform from J2000 to JD.

AE_TO_J2000 

Transform from JD to J2000.

Error codes that functions might return.

Functions will generally return negative values, e.g., -AE_RET_BAD_PATH == -1, not AE_RET_BAD_PATH == 1.

Enumerator
AE_RET_BAD_PATH 

A path could not be accessed for I/O.

AE_RET_UNEXPECTED_EOF 

An unexpected EOF was found in a file.

AE_RET_READ_ERROR 

An unknown read error occured in a file.

AE_RET_BAD_JPL_HEADER 

The JPL header was corrupted.

AE_RET_BAD_JPL_HEADER_KSIZE 

Bad JPL header; ksize != 2 * n_coeff.

AE_RET_BAD_JPL_HEADER_GROUP 

Bad JPL header; could not find a group.

AE_RET_BAD_JPL_HEADER_NUM 

Bad JPL header; differing numbers of constants and values.

AE_RET_BAD_JPL_HEADER_RANGE 

Bad JPL header; could not read date range.

AE_RET_BAD_JPL_NOT_ASCII 

Was expecting a JPL ASCII file.

AE_RET_BAD_JPL_NOT_BIN 

Was expecting a JPL binary file.

AE_RET_BAD_JPL_BLOCK 

Bad JPL block.

AE_RET_BAD_JPL_CORRUPT 

Data file seems corrupt.

AE_RET_BAD_JPL_DATE 

Julian date outside JPL ephemeris range.

AE_RET_BAD_JPL_FLOAT 

Unable to parse float in JPL file.

AE_RET_BAD_JPL_CHEBY_TIME 

Bad chebyshev time when gettin JPL data.

AE_RET_BAD_OBJ_INDEX 

Unknown object number.

AE_RET_OBJ_NOT_FOUND 

Object not found in a catalogue file.

AE_RET_CFITSIO_NOT_ENABLED 

A function requiring cfitsio was called, but cfitsio is not available.

AE_RET_CFITSIO_NO_KEYWORD 

A required FITS keyword did not exist.

AE_RET_CFITSIO_BAD_KEYWORD 

A required FITS keyword was corrupt.

AE_RET_WCS_UNKNOWN_PROJ 

An unknown projection type was detected.

AE_RET_CFITSIO_ERROR 

An error occurred in a call to libcfitsio.

AE_RET_BAD_NATIVE_POLE 

Could not compute projection's native pole.

AE_RET_NO_INVERSE 

A matrix was not invertible.

AE_RET_UNDER_CONSTRUCTION 

The function does not work as it is still being developed!

Enumerate Solar System bodies.

The first AE_N_SS_BODIES_JPL are in the the JPL ephemeris order. The rest are derived from these.

Enumerator
AE_SS_MERCURY 

Mercury, the Winged Messenger.

AE_SS_VENUS 

Venus, the Bringer of Peace.

AE_SS_EMBARY 

Earth-moon barycentre.

AE_SS_MARS 

Mars, the Bringer of War.

AE_SS_JUPITER 

Jove, the Bringer of Jollity.

AE_SS_SATURN 

Saturn, the Bringer of Old Age.

AE_SS_URANUS 

Uranus, the Magician.

AE_SS_NEPTUNE 

Neptune, the Mystic.

AE_SS_PLUTO 

Pluto, now a dwarf planet.

AE_SS_MOON_EMB 

Our moon, relative to the Earth-moon barycentre.

AE_SS_SUN 

The sun.

AE_SS_NUTATION 

Nutation.

AE_SS_LIBRATION 

Libration.

AE_SS_EARTH 

The earth.

AE_SS_MOON 

The moon.

AE_SS_SSBARY 

The Solar System Barycenter.

Function Documentation

void ae_altaz_to_radec ( double  last,
double  glat,
double  alt,
double  az,
double *  ra,
double *  dec 
)

Convert alt/az to ra/dec.

Parameters
lastThe local aparent sidereal time, in degrees.
glatThe geodetic latitude, in degrees.
altThe altitude, in degrees.
azThe azimuth, in degrees.
raFor returning the right ascension, in degrees.
decFor returning the declination, in degrees.

References AE_DTR, ae_mod_360(), AE_RTD, and ae_zatan2().

Referenced by ae_topocentric().

void ae_annual_aberration ( double  v_earth[],
double  p[],
int  direction 
)

Correct for the annual aberration.

See AA pages B17, B37, C24.

Parameters
v_earthThe heliocentric rectangular velocity of Earth in AU per day.
pA unit vector pointing from the earth to the object; the corrected position is returned in this parameter.
directionTo add aberration to p, pass AE_ADD_ABERRATION; to remove it, pass AE_REMOVE_ABERRATION; if neither is passed, then AE_REMOVE_ABERRATION is assumed. AE_ADD_ABERRATION is appropriate, for example, if one wishes to compute the apparent position of a star from a catalogue position, while AE_REMOVE_ABERRATION is appropriate, for example, if one wants to compute the true position of an observed star.

References AE_ADD_ABERRATION, AE_CLIGHT_AUD, and ae_stokes_t::i.

Referenced by ae_geocentric(), and ae_geocentric_from_cat().

double ae_cal_to_jd ( long  year,
int  month,
double  day 
)

Calculate Julian date from a Gregorian calendar date.

The Julian date is double precision floating point with the origin used by astronomers.

There is no year 0. Enter B.C. years as negative; i.e., 2 B.C. = -2.

The approximate range of dates handled is 4713 B.C. to 54,078 A.D. This should be adequate for most applications.

B.C. dates are calculated by extending the Gregorian sequence of leap years and century years into the past. This seems the only sensible definition, but it may not be the official one.

Note that the astronomical Julian day starts at noon on the previous calendar day. Thus at midnight in the morning of the present calendar day the Julian date ends in 0.5; it rolls over to tomorrow at noon today.

The month-finding algorithm is attributed to Meeus.

Parameters
yearThe Gregorian year.
monthThe month of the year, with January = 1, February = 2, etc.
dayThe fractional day of the month, starting at 1.
Returns
The Julian day number.

Referenced by ae_ctime_to_jd().

int ae_cat_to_constel_index ( const char *  cat_name,
char *  full_name,
int  len 
)

Given a catalogue star name, find a star's constellation.

Optionally get a full star-constellation name.

Parameters
cat_nameThe star name from the catalogue. It should be of the form {2-char Greek letter abbrev.}[-#]{3-char constellation name}({star name}). For example: `beCet(Diphda)', `al-2Lib(Zubenelgen)'.
full_nameFor returning an expanded star name; set to NULL if you do not need this.
lenThe length of full_name; ignored if full_name is NULL.
Returns
The index of ae_constel_name representing the constellation that this star lives in; or -1 if it cannot be found.

References AE_CONSTEL_N_CON, AE_CONSTEL_N_GREEK, ae_constel_name, ae_greek_letter_name, and ae_skip_white().

int ae_coord_to_constel_index ( double  ra,
double  dec,
double  epoch 
)

Given a sky coordinate, find the constellation it is in.

Parameters
raThe right ascension.
decThe declination.
epochThe Julian date for the ra/dec coordinates.
Returns
The index of ae_constel_name representing the constellation that this coordinate lives in; or -1 if it cannot be found.

References ae_constel_boundary, AE_CONSTEL_N_BOUNDARIES, AE_FROM_J2000, ae_polar_to_rect(), ae_precess(), AE_RTD, and AE_TO_J2000.

double ae_ctime_to_jd ( double  t)

Convert a C time to a Julian date.

The C time is the time since the Epoch (00:00:00 UTC, January 1, 1970), measured in seconds.

Parameters
tThe C time at which to calculate the Julian date.
Returns
The Julian date.

References ae_cal_to_jd(), and AE_STDAY.

Referenced by ae_ctime_to_last(), and aes_ctime_to_last().

double ae_ctime_to_last ( double  t,
double  delta_t,
double  tlong,
double  nutl,
double  eps 
)

Calculate the local apparent sidereal time from a C time.

The C time is the time since the Epoch (00:00:00 UTC, January 1, 1970), measured in seconds.

Parameters
tThe UNIX time in UT1 measure at which to calculate the LST.
delta_tTT - UT1 (see ae_delta_t()), in seconds.
tlongThe longitude of the observer, in degrees.
nutlThe nutation in longitude in seconds of arc (see ae_nutation_lon_ob()).
epsThe obliquity of the ecliptic in seconds of arc (see ae_epsilon()).
Returns
The LAST in degrees.

References ae_ctime_to_jd(), AE_D_PER_S, ae_dut1(), and ae_last().

void ae_delta_q ( double  q0[],
double  q1[],
double *  dra,
double *  ddec 
)

Convert change in rectangular coordinatates to change in ra/dec.

For changes greater than about 0.1 degree, the coordinates are converted directly to ra and dec and the results subtracted. For small changes, the change is calculated to first order by differentiating $tan(ra) = y/x$ to obtain

\[\frac{dra}{\cos^2(ra)} = \frac{dy}{x} - \frac{y dx}{x^2},\]

where

\[\cos^2(ra) = \frac{1}{1 + (y/x)^2}.\]

The change in declination $arcsin(z/R)$ is

\[d\arcsin(u) = \frac{du}{\sqrt{1-u^2}},\]

where $u = z/R.$

Parameters
q0The initial object - earth vector.
q1The vector after motion or aberration.
draThe change in right ascension.
ddecThe change in declination.

References ae_zatan2().

double ae_delta_t ( double  jd_ut1)

Get the value of delta T == Ephemeris Time - Universal Time.

This routine uses the table ae_delta_t_tab for historic values and near-future predictions.

The program adjusts for a value of secular tidal acceleration ndot. It is -25.8 arcsec per century squared for JPL's DE403 ephemeris. ELP2000 and DE200 use the value -23.8946.

For dates earlier than the tabulated range, the program calculates approximate formulae of Stephenson and Morrison or K. M. Borkowski. These approximations have an estimated error of 15 minutes at 1500 B.C. They are not adjusted for small improvements in the current estimate of ndot because the formulas were derived from studies of ancient eclipses and other historical information, whose interpretation depends only partly on ndot.

A quadratic extrapolation formula, that agrees in value and slope with current data, predicts future values of deltaT.

References:

  • Stephenson, F. R., and L. V. Morrison, "Long-term changes in the rotation of the Earth: 700 B.C. to A.D. 1980," Philosophical Transactions of the Royal Society of London Series A 313, 47-70 (1984)
  • Borkowski, K. M., "ELP2000-85 and the Dynamical Time - Universal Time relation," Astronomy and Astrophysics 205, L8-L10 (1988)
    • Borkowski's formula is derived from eclipses going back to 2137 BC and uses lunar position based on tidal coefficient of -23.9 arcsec/cy^2.
  • Chapront-Touze, Michelle, and Jean Chapront, Lunar Tables and Programs from 4000 B.C. to A.D. 8000, Willmann-Bell 1991
    • Their table agrees with the one here, but the entries are rounded to the nearest whole second.
  • Stephenson, F. R., and M. A. Houlden, Atlas of Historical Eclipse Maps, Cambridge U. Press (1986)
Parameters
jd_ut1The Julian date.
Returns
The ET - UT (delta T) in seconds.

References ae_delta_t_end_year, ae_delta_t_start_year, ae_delta_t_tab, and AE_J2000.

Referenced by aes_last(), aes_subobs_point(), and aes_topocentric().

double ae_disc_semiminor ( double  a,
double  c,
double  lat 
)

Calculate the semi-minor axis of an observed disc.

When an oblate spheroid is viewed, it appears as an ellipse (to a sufficiently distant observer). This routine calculates the semi-minor axis (i.e., polar radius) of the projected ellipse.

Parameters
aThe major axis (i.e., equatorial radius) of the planet.
cThe minor axis (i.e., polar radius) of the planet.
latThe sub-observer latitude, in degrees.
Returns
The semi-minor axis of the projected ellipse. Note that the semi-major axis of the ellipse is the same as the input value.

References AE_DTR.

Referenced by aes_disc_solid_angle().

double ae_disc_solid_angle ( double  a,
double  b,
double  dist 
)

Calculate the solid angle of an observed disc.

It is assumed that the disc is elliptical. The semi-minor axis of the projected ellipse can be obtained from ae_disc_semiminor(). Note that this calculation assumes that the size of the disc is small compared to the distance of the observer.

Parameters
aThe semi-major axis of the disc, in kilometres.
bThe semi-minor axis of the disc, in kilometres.
distThe distance from the observer to the planet, in AU.
Returns
The solid angle of the disc, in steradians.

References AE_AU.

Referenced by aes_disc_solid_angle().

void ae_diurnal_aberration ( double  last,
double  tlat,
double  trho,
int  direction,
double *  ra,
double *  dec 
)

Correct for diurnal aberration.

This formula is less rigorous than the method used for annual aberration (see ae_annual_aberration()). However, the correction is small.

Parameters
lastThe local apparent sidereal time, in degrees.
tlatThe geocentric latitude of the observer, in degrees.
trhoThe distance from the centre of the earth to the observer, in earth radii.
directionTo add aberration to p, pass AE_ADD_ABERRATION; to remove it, pass AE_REMOVE_ABERRATION; if neither is passed, then AE_REMOVE_ABERRATION is assumed. AE_ADD_ABERRATION is appropriate, for example, if one wishes to compute the apparent position of a star from a catalogue position, while AE_REMOVE_ABERRATION is appropriate, for example, if one wants to compute the true position of an observed star.
raThe right ascension of the object, in degrees; this routine modifies this parameter to correct for diurnal aberration.
decThe declination of the object, in degrees; this routine modifies this parameter to correct for diurnal aberration.

References AE_ADD_ABERRATION, AE_DTR, and AE_RTD.

Referenced by ae_topocentric().

void ae_diurnal_parallax ( double  last,
double  tlat,
double  trho,
double  dist,
double *  ra,
double *  dec 
)

Correct for diurnal parallax.

See AA page D3.

This function does not bother to calculate anything unless the equatorial horizontal parallax is at least 0.005".

Parameters
lastThe local apparent sidereal time, in degrees.
tlatThe geocentric latitude of the observer, in degrees.
trhoThe distance from the centre of the earth to the observer, in earth radii.
distThe earth-object distance, in AU.
raThe right ascension of the object, in degrees; this routine modifies this parameter to correct for diurnal parallax.
decThe declination of the object, in degrees; this routine modifies this parameter to correct for diurnal parallax.

References AE_AU, AE_DTR, AE_R_EARTH, AE_RTD, and ae_zatan2().

Referenced by ae_topocentric().

double ae_dut1 ( double  jd_utc)

Given a UTC, return DUT1.

DUT1 is the offset between UTC and UT1, i.e., $ DUT1 \equiv UT1 - UTC $. Functions that find sidereal time expect UT1, so this function will be required to convert UTC if reasonable precision (better than roughly 10 seconds of arc).

Values are taken from a look-up table. Past values start on 19 May 1976 are are tabulated daily until the release date of this version of aephem. After this, predicted values are used through about one year.

If DUT1 is requested for a date before the tabulation starts, 0 will be returned. After the last predicted value, a value will be estimated using an interpolation formula.

Past values, predicted values and the interpolation formula were taken from the U.S. Naval Observatory IERS Bulletin A. This is available online at http://maia.usno.navy.mil/.

Parameters
jd_utcThe Julian date in UTC.
Returns
DUT1, in seconds. If jd_utc is before the first tabulated value, 0 will be returned.

References ae_dut1_first_mjd, ae_dut1_last_mjd, ae_dut1_predict_mjd, ae_dut1_predict_offset, ae_dut1_predict_scale, ae_dut1_tab, and ae_mjd().

Referenced by ae_ctime_to_last(), and aes_ctime_to_last().

double ae_epsilon ( double  jd_tt)

Calculate the obliquity of the ecliptic.

J. L. Simon, P. Bretagnon, J. Chapront, M. Chapront-Touze', G. Francou, and J. Laskar, "Numerical Expressions for precession formulae and mean elements for the Moon and the planets," Astronomy and Astrophysics 282, 663-683 (1994).

IAU Coefficients are from: J. H. Lieske, T. Lederle, W. Fricke, and B. Morando, "Expressions for the Precession Quantities Based upon the IAU (1976) System of Astronomical Constants," Astronomy and Astrophysics 58, 1-16 (1977).

Before or after 200 years from J2000, the formula used is from: J. Laskar, "Secular terms of classical planetary theories using the results of general theory," Astronomy and Astrophysics 157, 59070 (1986).

See precess.c and page B18 of the Astronomical Almanac.

Parameters
jd_ttThe Julian date in TT.
Returns
The obliquity, in seconds of arc.

Referenced by ae_geocentric_moon_from_orbit(), ae_gmoon(), ae_kepler(), ae_precess(), ae_rect_to_polar_with_jd(), ae_topocentric(), aes_last(), and aes_nutate().

void ae_fk4_to_fk5 ( struct ae_star_t el)

Convert FK4 B1950.0 catalogue coordinates to FK5 J2000.0 coordinates.

See AA page B58.

Parameters
elThe catalogue object to update to J2000.0.

References AE_DTR, ae_fk4_to_fk5_abb, ae_fk4_to_fk5_abb_d, ae_fk4_to_fk5_mat, AE_J2000, AE_RTD, ae_zatan2(), ae_star_t::dec, ae_star_t::epoch, ae_star_t::mudec, ae_star_t::mura, ae_star_t::px, ae_star_t::ra, and ae_star_t::v.

Referenced by ae_geocentric_from_cat().

double ae_flattening ( const struct ae_physical_t phys)

Calculate a body's flattening.

Takes the equatorial radius, $a$, and polar radius, $c$, and calculates the flattening of the oblate sphere:

\[f \equiv \frac{a - b}{a}\]

Parameters
physThe physical parameters of the body.
Returns
The flattening.

References ae_physical_t::r_eq, and ae_physical_t::r_pole.

Referenced by aes_subobs_point().

double ae_g1plan ( double  jd_tt,
struct ae_plantbl_t plan 
)

Accumulate the sum of trigonometric series in one variable.

This is a generic routine to accumulate sum of trigonometric series in one variable of the same list of arguments.

Parameters
jd_ttThe Julian date in TT.
planA table of planetary orbit data.
Returns
The accumlated sum.

References ae_gplan_cc, ae_gplan_mean_element(), ae_gplan_ss, ae_gplan_sscc(), AE_J2000, AE_PLANTBL_N_HARMONIC, ae_plantbl_t::arg_tbl, ae_plantbl_t::lon_tbl, ae_plantbl_t::max_harmonic, ae_plantbl_t::timescale, and ae_plantbl_t::trunclvl.

Referenced by ae_gmoon().

void ae_g2plan ( double  jd_tt,
struct ae_plantbl_t plan,
double  pobj[],
double *  lp_equinox 
)

Accumulate the sum of trigonometric series in two variables.

This is a generic routine to accumulate sum of trigonometric series in two variables (e.g., longitude, radius) of the same list of arguments.

Parameters
jd_ttThe Julian date in TT.
planA table of planetary orbit data.
pobjFor returning the polar position of the planet.
lp_equinoxFor returning the lp_equinox (?); set to NULL if it is not not needed

References ae_gplan_cc, ae_gplan_mean_element(), ae_gplan_ss, ae_gplan_sscc(), AE_J2000, AE_PLANTBL_N_HARMONIC, AE_STR, ae_plantbl_t::arg_tbl, ae_plantbl_t::lon_tbl, ae_plantbl_t::max_harmonic, ae_plantbl_t::maxargs, ae_plantbl_t::rad_tbl, ae_plantbl_t::timescale, and ae_plantbl_t::trunclvl.

Referenced by ae_gmoon().

void ae_g3plan ( double  jd_tt,
struct ae_plantbl_t plan,
double  pobj[],
int  objnum 
)

Accumulate the sum of trigonometric series in three variables.

This is a generic program to accumulate sum of trigonometric series in three variables (e.g., longitude, latitude, radius) of the same list of arguments.

Parameters
jd_ttThe Julian date in TT.
planA table of planetary orbit data.
pobjFor returning the polar position of the planet.
objnumThe planet number.

References ae_gplan_cc, ae_gplan_mean_element(), ae_gplan_ss, ae_gplan_sscc(), AE_J2000, AE_PLANTBL_N_HARMONIC, AE_STR, ae_plantbl_t::arg_tbl, ae_plantbl_t::distance, ae_plantbl_t::lat_tbl, ae_plantbl_t::lon_tbl, ae_plantbl_t::max_harmonic, ae_plantbl_t::maxargs, ae_plantbl_t::rad_tbl, ae_plantbl_t::timescale, and ae_plantbl_t::trunclvl.

Referenced by ae_kepler().

void ae_gal_to_radec ( double  l,
double  b,
double *  ra,
double *  dec,
int  fk5 
)

Convert l/b (galactic coordinates) to ra/dec.

Parameters
lThe galactic longitude, in degrees.
bThe galactic latitude, in degrees.
raFor returning the right ascension, in degrees.
decFor returning the declination, in degrees.
fk5If non-zero, assume FK5 system (i.e., J2000 compatible); if zero, use the old FK4 (B1950 compatible) galactic pole.

References AE_DTR, ae_mod_360(), and AE_RTD.

double ae_gast ( double  jd_ut1,
double  jd_tt,
double  nutl,
double  eps 
)

Get Greenwich apparent sidereal time.

The Greenwich apparent sidereal time is the Greenwich mean sidereal time corrected for nutation. See ae_gmst() for reference.

Parameters
jd_ut1The Julian date, in UT1.
jd_ttThe Julian date, in TT. (See ae_delta_t().)
nutlThe nutation in longitude in seconds of arc (see ae_nutation_lon_ob()).
epsThe obliquity of the ecliptic in seconds of arc (see ae_epsilon()).
Returns
Greenwich apparent sidereal time in degrees.

References ae_gmst(), ae_mod_360(), AE_STD, and AE_STR.

Referenced by ae_last().

void ae_geocentric ( double  jd_tt,
double  q[],
double  e[],
double  v_e[],
double *  ra,
double *  dec,
double *  dist 
)

Given heliocentric coordinates, reduce to geocentric coordinates.

This routine reduces heliocentric coordinates to equitorial coordinates. Unlike ae_geocentric_from_orbit() or ae_geocentric_from_jpl(), it does not automatically compute the heliocentric coordinates. This means that light travel-time is not automatically accounted for! The input position of the planet is what is used.

Corrections done during the reduction are: light deflection, annual aberration, precession and nutation.

This routine should be thought of as a macro. One can call the functions that perform the corrections listed above individually. This routine simply wraps them into one call.

Parameters
jd_ttThe Julian date in TT.
qThe heliocentric rectangular coordinates of the object being observed, in AU.
eThe heliocentric rectangular coordinates of the earth, in AU.
v_eThe velocity of the earth, in AU per day.
raFor returning the right ascension, in degrees.
decFor returning the declination, in degrees.
distFor returning the geocentric distance between earth and the object, in AU. Set to NULL if you do not require this value.

References AE_ADD_ABERRATION, ae_annual_aberration(), AE_FROM_J2000, ae_precess(), ae_rect_to_polar(), ae_relativity(), and aes_nutate().

Referenced by ae_geocentric_from_jpl(), ae_geocentric_from_orbit(), ae_geocentric_sun_from_orbit(), and aes_subobs_point().

void ae_geocentric_from_cat ( double  jd_tt,
double  e[],
double  v_e[],
const struct ae_star_t star,
double *  ra,
double *  dec 
)

Reduce a star's catalogue coordinates to its apparent place.

This routine puts the coordinates in J2000.0 (if not already there), and performs precession, proper motion, parallax, light deflection, aberration and nutation.

This routine should be thought of as a macro. One can call the functions that perform the corrections listed above individually. This routine simply wraps them into one call.

Parameters
jd_ttThe Julian date in TT.
eThe heliocentric rectangular coordinates of the earth, in AU. These can be obtained from ae_kepler() or ae_jpl_get_coords().
v_eThe velocity of the earth, in AU per day. This can be obtained from ae_v_orbit() or ae_jpl_get_coords().
starThe catalogue object. If it is in the FK4 B1950.0 coordinates, this routine will convert it to FK5 J2000.0 coordinates.
raFor returning the apparent (geocentric) right ascension, in degrees.
decFor returning the apparent (geocentric) declination, in degrees.

References AE_ADD_ABERRATION, ae_annual_aberration(), AE_B1950, AE_DTR, ae_fk4_to_fk5(), AE_FROM_J2000, ae_precess(), ae_rect_to_polar(), ae_relativity(), AE_STR, AE_TO_J2000, aes_nutate(), ae_star_t::dec, ae_star_t::epoch, ae_star_t::mudec, ae_star_t::mura, ae_star_t::px, ae_star_t::ra, and ae_star_t::v.

int ae_geocentric_from_jpl ( struct ae_jpl_handle_t jh,
double  jd_tt,
int  obj_num,
double *  ra,
double *  dec,
double *  dist 
)

Determine geocentric coordinates using JPL ephemerides.

This routine reduces the heliocentric equitorial coordinates, obtained by calling ae_jpl_get_coords(), to the geocentric right ascension and declination. The position is corrected for light travel time, light deflection, annual aberration, precession and nutation.

It is assumed that the JPL planetary ephemeris file is being used, since it automatically polls the file for Earth's position and velocity. If you are using a different ephemeris file, you cannot use this routine.

This routine should be thought of as a macro. One can call ae_jpl_get_coords() and the functions that perform the corrections listed above individually. This routine simply wraps them into one call.

To obtain geocentric coordinates from an orbital elements struct, use ae_geocentric_from_orbit().

Parameters
jhA handle initialised by ae_jpl_init().
jd_ttThe Julian date in TT.
obj_numThe object for which to get coordinates. If a planetary ephemeris file is being used, one can use the predefined ae_ss_bodies_t constants for clarity.
raFor returning the geocentric right ascension, in degrees.
decFor returning the geocentric declination, in degrees.
distFor returning the geocentric distance between earth and the object, in AU. Set to NULL if you do not require this value.
Returns
0 on success; on failure, a negative error code from ae_retcode_t.

References ae_geocentric(), ae_jpl_get_coords(), ae_light_t_jpl(), and AE_SS_EARTH.

void ae_geocentric_from_orbit ( double  jd_tt,
const struct ae_orbit_t o_orb,
const struct ae_orbit_t q_orb,
double *  ra,
double *  dec,
double *  dist 
)

Determine geocentric coordinates using orbital elements.

This routine reduces the heliocentric equitorial coordinates, obtained by calling ae_kepler(), to the geocentric right ascension and declination. The position is corrected for light travel time, light deflection, annual aberration, precession and nutation.

This routine should be thought of as a macro. One can call ae_kepler() and the functions that perform the corrections listed above individually. This routine simply wraps them into one call.

To obtain geocentric coordinates from a JPL ephemeris file, use ae_geocentric_from_jpl().

Parameters
jd_ttThe Julian date in TT.
o_orbOrbital elements of the observer.
q_orbThe orbital elements of the object.
raFor returning the geocentric right ascension, in degrees.
decFor returning the geocentric declination, in degrees.
distFor returning the geocentric distance between earth and the object, in AU. Set to NULL if you do not require this value.

References ae_geocentric(), ae_kepler(), ae_light_t_orbit(), and ae_v_orbit().

void ae_geocentric_lat ( double  glat,
double  height,
double *  tlat,
double *  trho 
)

Calculate the geocentric latitude and the distance to the earth's centre.

This function uses the reciprocal of flattening and the earth's radius (or, more precisely, semi-major axis) from from the WGS84 definition. (See AE_FLAT and AE_R_EARTH.)

Parameters
glatThe geodetic latitude of the observer, in degrees.
heightThe height of the observer above sea level, in metres.
tlatFor returning the geocentric latitude, in degrees.
trhoFor returning the distance to the earth's centre, in earth radii.

References AE_DTR, AE_FLAT, AE_R_EARTH, and AE_RTD.

Referenced by aes_topocentric().

void ae_geocentric_moon_from_orbit ( double  jd_tt,
const struct ae_orbit_t o_orb,
double *  ra,
double *  dec,
double *  dist 
)

Calculate the geocentric coordinates of the moon from orbital elements.

Corrects for precession, nutation, annual aberration. Neglects light travel time.

Note that if using JPL ephemeris, one can use ae_geocentric_from_jpl() to do this calculation since the moon is an object in those ephemerides.

The Astronomical Almanac (Section D, Daily Polynomial Coefficients) seems to omit the annual aberration at one point, even though the reference ephemeris is inertial . . .

Parameters
jd_ttThe Julian date in TT.
o_orbOrbital elements of the observer.
raFor returning the geocentric right ascension, in degrees.
decFor returning the geocentric declination, in degrees.
distFor returning the geocentric distance between earth and the object, in AU. Set to NULL if you do not require this value.

References AE_DTR, ae_epsilon(), AE_FROM_J2000, ae_gmoon(), ae_kepler(), ae_nutate(), ae_nutation_lon_ob(), AE_R_EARTH_AU, ae_rect_to_polar(), and AE_STR.

void ae_geocentric_sun_from_orbit ( double  jd_tt,
const struct ae_orbit_t o_orb,
double *  ra,
double *  dec,
double *  dist 
)

Calculate the geocentric coordinates of the sun from orbital elements.

Corrects for precession, nutation, annual aberration.

Note that if using JPL ephemeris, one can use ae_geocentric_from_jpl() to do this calculation since the sun is an object in those ephemerides.

Parameters
jd_ttThe Julian date in TT.
o_orbOrbital elements of the observer.
raFor returning the geocentric right ascension, in degrees.
decFor returning the geocentric declination, in degrees.
distFor returning the geocentric distance between earth and the object, in AU. Set to NULL if you do not require this value.

References ae_geocentric(), ae_kepler(), and ae_v_orbit().

void ae_gmoon ( double  jd_tt,
double  rect[],
double  pol[] 
)

Compute geocentric moon position.

Parameters
jd_ttThe Julian date in TT.
rectFor returning the position in rectangular coordinates.
polFor returning the position in polar coordinates.

References ae_epsilon(), ae_g1plan(), ae_g2plan(), ae_mlat404, ae_mlr404, AE_STR, and ae_plantbl_t::distance.

Referenced by ae_geocentric_moon_from_orbit(), and ae_kepler_embofs().

double ae_gmst ( double  jd_ut1,
double  jd_tt 
)

Get the Greenwich mean sidereal time.

Get the mean sidereal time at Greenwich (i.e., longitude 0). The mean sidereal time does not include nutation corrections. Coefficients are from the IAU and can be found in:

  • George H. Kaplan, "The IAU Resolutions on Astronomical Reference Systems, Time Scales, and Earth Rotation Models," United States Naval Observatory Circular No. 179, 2005.
Parameters
jd_ut1The Julian date, in UT1.
jd_ttThe Julian date, in TT. (See ae_delta_t().)
Returns
The Greenwich mean sidereal time, in degrees.

References AE_J2000, ae_mod_360(), and AE_STD.

Referenced by ae_gast(), and ae_lmst().

void ae_gplan ( double  jd_tt,
struct ae_plantbl_t plan,
double  pobj[] 
)

Use a table to find a planet's polar, heliocentric position.

Parameters
jd_ttThe Julian date in TT.
planThe planet's orbital table.
pobjFor returning the polar position of the planet.

References ae_gplan_cc, ae_gplan_freq, ae_gplan_mods3600, AE_GPLAN_N_HARMONIC, ae_gplan_phase, ae_gplan_ss, ae_gplan_sscc(), AE_J2000, AE_STR, ae_plantbl_t::arg_tbl, ae_plantbl_t::distance, ae_plantbl_t::lat_tbl, ae_plantbl_t::lon_tbl, ae_plantbl_t::max_harmonic, ae_plantbl_t::maxargs, ae_plantbl_t::rad_tbl, and ae_plantbl_t::timescale.

Referenced by ae_kepler().

int ae_is_retrograde ( const struct ae_physical_t phys)

Determine whether a body's rotation is retrograde.

This function merely looks at the sign of the variable w_d in the struct ae_physical_t: negative values are retrograde.

Parameters
physThe body to examine.
Returns
1 if the rotation is retrograde; 0 if the rotation is prograde.

References ae_physical_t::w_d.

Referenced by aes_subobs_point().

void ae_jd_to_cal ( double  jd,
int *  year,
int *  month,
double *  day 
)

Calculate month, day, and year from a Julian date.

Parameters
jdA Julian day number.
yearFor returning the year of jd.
monthFor returning the month of jd, in the range [1:12].
dayFor returning the fractional day of jd, with the first day of the month being 1.
void ae_jpl_close ( struct ae_jpl_handle_t jh)

Finish with a JPL handle.

Closes any open file pointers and frees any allocated arrays.

Parameters
jhA JPL handle.

References ae_jpl_close_wrap().

double ae_jpl_get_const_val ( const struct ae_jpl_handle_t jh,
const char *  const_name 
)

Get the value of a constant defined in the JPL ephemeris header.

Parameters
jhAn initialised JPL handle.
const_nameThe name of the constant, as it appears in the header.
Returns
The constant value; 0 if the constant name does not exist.

References AE_JPL_CONST_LEN, ae_jpl_handle_t::const_name, ae_jpl_handle_t::const_val, and ae_jpl_handle_t::n_const.

Referenced by ae_jpl_init_ascii().

int ae_jpl_get_coords ( struct ae_jpl_handle_t jh,
double  jd_tt,
int  obj_num,
double  r[],
double  v[],
char  is_planetary 
)

Get object positions from a JPL ephemeris file.

Parameters
jhA handle initialised by ae_jpl_init().
jd_ttThe Julian date in TT.
obj_numThe object for which to get coordinates. If a planetary ephemeris file is being used, one can use the predefined ae_ss_bodies_t constants for clarity.
rFor returning the rectangular position coordinates, in AU; for nutations and librations, units are seconds of arc.
vFor returning the rectangular velocity, in AU per day (or seconds of arc per day for nutations or librations) Set to NULL if you do not require this information (though setting to NULL does not make this routine faster).
is_planetaryIf 0, assume these are not planetary ephemerides and barycentric coordinates are returned. Otherwise, it will be assumed that a Solar System ephemeris file is being used and that the object numbers are as in ae_ss_bodies_t; heliocentric coordinates will be returned.
Returns
0 on success; on failure, a negative code from ae_retcode_t.

References AE_AU, ae_jpl_cheby(), ae_jpl_get_block(), AE_N_SS_BODIES_JPL, AE_RET_BAD_JPL_DATE, AE_RET_BAD_OBJ_INDEX, AE_RTS, AE_SS_EARTH, AE_SS_EMBARY, AE_SS_LIBRATION, AE_SS_MOON, AE_SS_MOON_EMB, AE_SS_NUTATION, AE_SS_SSBARY, AE_SS_SUN, ae_jpl_handle_t::block, ae_jpl_handle_t::em_ratio, ae_jpl_handle_t::end, ae_jpl_handle_t::start, ae_jpl_handle_t::step, ae_jpl_handle_t::sun_pv, and ae_jpl_handle_t::sun_pv_t.

Referenced by ae_geocentric_from_jpl(), and ae_light_t_jpl().

int ae_jpl_init ( const char *  path,
struct ae_jpl_handle_t jh 
)

Initialise to use a JPL ephemeris file.

This function opens the file, determines whether it is binary or ASCII, and parses the header information. This information is stored in the handle jh and is passed to subsequent JPL function calls.

If the ephemeris is in ASCII format and has a separate header file, use ae_jpl_init_ascii() instead.

If you know the format of the ephemeris file, you can use ae_jpl_init_bin() or ae_jpl_init_ascii() directly.

Parameters
pathThe location of the ephemeris file.
jhFor returning an initialised handle.
Returns
0 on success; otherwise, a negative code from ae_retcode_t; this assuming that it is a binary file, which is the last format attempted.

References ae_jpl_init_ascii(), ae_jpl_init_bin(), ae_jpl_handle_t::fp, and ae_jpl_handle_t::is_bin.

int ae_jpl_init_ascii ( const char *  path,
const char *  header_path,
int  search_dates,
struct ae_jpl_handle_t jh 
)

Read a JPL ephemeris ASCII header.

Parameters
pathThe location of the ephemeris file.
header_pathSet to NULL, unless the header file is separate (as it might be an ASCII ephermeris).
search_datesSince JPL often splits its ASCII ephermerides into multiple files, the date range specified in the header might not correspond to the file being used. To search through the file to get the dates, set this to 1 (recommended). To trust the header, set to 0.
jhFor returning an initialised handle, on success.
Returns
0 on success; on error, a negative error code from ae_retcode_t.

References AE_AU, AE_E_M_RAT, ae_jpl_ascii_find_group(), ae_jpl_close_wrap(), AE_JPL_CONST_LEN, AE_JPL_DATA_GROUP, ae_jpl_get_const_val(), ae_jpl_init_common_end(), ae_jpl_init_common_start(), AE_JPL_MAXLINE, ae_jpl_read_ascii_float(), AE_RET_BAD_JPL_CORRUPT, AE_RET_BAD_JPL_HEADER, AE_RET_BAD_JPL_HEADER_KSIZE, AE_RET_BAD_JPL_HEADER_NUM, AE_RET_BAD_JPL_HEADER_RANGE, AE_RET_BAD_PATH, AE_RET_UNEXPECTED_EOF, ae_jpl_handle_t::ascii_start, ae_jpl_handle_t::au, ae_jpl_handle_t::const_name, ae_jpl_handle_t::const_val, ae_jpl_handle_t::em_ratio, ae_jpl_handle_t::end, ae_jpl_handle_t::ephemeris_version, ae_jpl_handle_t::fp, ae_jpl_handle_t::ipt, ae_jpl_handle_t::is_bin, ae_jpl_handle_t::n_coeff, ae_jpl_handle_t::n_const, ae_jpl_handle_t::rec_size, ae_jpl_handle_t::start, and ae_jpl_handle_t::step.

Referenced by ae_jpl_init().

int ae_jpl_init_bin ( const char *  path,
struct ae_jpl_handle_t jh 
)
void ae_kepler ( double  jd_tt,
const struct ae_orbit_t orb,
double  q[] 
)

Find a Keplerian orbit.

This routine solves for a Keplerian orbit, given orbital parameters and the time.

The routine detects several cases of given orbital elements. If a program for perturbations is pointed to, it is called to calculate all the elements. If there is no program, then the mean longitude is calculated from the mean anomaly and daily motion. If the daily motion is not given, it is calculated by Kepler's law. If the eccentricity is given to be 1.0, it means that meandistance is really the perihelion distance, as in a comet specification, and the orbit is parabolic.

Reference:

  • Taff, L.G., "Celestial Mechanics, A Computational Guide for the Practitioner." Wiley, 1985.
Parameters
jd_ttThe Julian date in TT.
orbThe orbital elements.
qFor returning the heliocentric equatorial rectangular coordinates of the object, in AU.

References ae_orbit_t::a, AE_DTR, ae_epsilon(), ae_g3plan(), ae_gplan(), AE_J2000, ae_kepler_embofs(), ae_mod_2pi(), ae_mod_360(), ae_precess(), AE_STR, AE_TO_J2000, ae_zatan2(), ae_orbit_t::dm, ae_orbit_t::ecc, ae_orbit_t::epoch, ae_orbit_t::equinox, ae_orbit_t::i, ae_orbit_t::L, ae_orbit_t::M, ae_orbit_t::name, ae_orbit_t::ptable, ae_orbit_t::W, and ae_orbit_t::w.

Referenced by ae_geocentric_from_orbit(), ae_geocentric_moon_from_orbit(), ae_geocentric_sun_from_orbit(), ae_light_t_orbit(), ae_v_orbit(), and aes_subobs_point().

double ae_last ( double  jd_ut1,
double  jd_tt,
double  tlong,
double  nutl,
double  eps 
)

Get the local apparent sidereal time.

The local apparent sidereal time is simply the Greenwich apparent sidereal time plus the local longitude. See ae_gast().

Parameters
jd_ut1The Julian date in UT1.
jd_ttThe Julian date in TT. (See ae_gast().)
tlongThe longitude of the observer in degrees.
nutlThe nutation in longitude in seconds of arc (see ae_nutation_lon_ob()).
epsThe obliquity of the ecliptic in seconds of arc (see ae_epsilon()).
Returns
Local sidereal time in degrees.

References ae_gast(), and ae_mod_360().

Referenced by ae_ctime_to_last(), ae_topocentric(), and aes_last().

double ae_light_t ( double  p[],
double  q[],
int  do_retardation 
)

Get light-travel time between two objects.

Given the position of two objects in heliocentric equatorial rectangular coordinates, returns the light-travel time between them. To correct for light-travel time, this should be done iteratively a couple of times. The routines ae_light_t_jpl() and ae_light_t_orbit() do this automatically.

This routine can approximate the effect due to the sun's gravitational retardation (AA, p. 36). However, it is usually a small effect. It will blow up if one of the objects is the sun.

Parameters
pThe first object, in heliocentric rectangular coordinates, in AU.
qThe second object, in heliocentric rectangular coordinates, in AU.
do_retardationIf 0, neglect gravitational retardation from the sun; otherwise, compute gravitational retardation from the sun.
Returns
The light-travel time, in days.

References AE_CLIGHT_AUD.

Referenced by ae_light_t_jpl(), ae_light_t_orbit(), and aes_subobs_point().

int ae_light_t_jpl ( struct ae_jpl_handle_t jh,
double  jd_tt,
double  o[],
int  obj_num,
double  q[],
double  v_q[],
char  is_planetary 
)

Get the apparent position of an object, corrected for light-travel time.

The routine does three iterations of light-time correction. Gravitational retardation from the sun is neglected. (However, it can be added if necessary by using ae_light_t() directly.)

This function uses JPL ephemerides. For orbital elements, use ae_light_t_orbit().

Parameters
jhA handle initialised by ae_jpl_init().
jd_ttThe Julian date in TT.
oThe heliocentric rectangular position of the earth, in AU.
obj_numThe object for which to get coordinates. If a planetary ephemeris file is being used, one can use the predefined ae_ss_bodies_t constants for clarity.
qFor returning the heliocentric rectangular position of the object in AU.
v_qFor returning the velocity of the object in AU. Set to NULL if you do not require the velocity.
is_planetaryIf 0, assume these are not planetary ephemerides and barycentric coordinates are returned. Otherwise, it will be assumed that a Solar System ephemeris file is being used and that the object numbers are as in ae_ss_bodies_t; heliocentric coordinates will be returned.
Returns
0 on success; on failure, a negative error code from ae_retcode_t.

References ae_jpl_get_coords(), and ae_light_t().

Referenced by ae_geocentric_from_jpl().

void ae_light_t_orbit ( double  jd_tt,
double  o[],
const struct ae_orbit_t orb,
double  q[] 
)

Get the apparent position of an object, corrected for light-travel time.

The routine does three iterations of light-time correction. Gravitational retardation from the sun is neglected. (However, it can be added if necessary by using ae_light_t() directly.)

This routine uses orbital elements. For JPL ephemerides, use ae_light_t_jpl().

Parameters
jd_ttThe Julian date in TT.
oThe heliocentric rectangular position of the observer, in AU.
orbThe orbital elements of the object.
qFor returning the heliocentric rectangular position of the object in AU.

References ae_kepler(), and ae_light_t().

Referenced by ae_geocentric_from_orbit().

double ae_lmst ( double  jd_ut1,
double  jd_tt,
double  tlong 
)

Get the local mean sidereal time.

The local mean sidereal time is simply the Greenwich mean sidereal time plus the local longitude. See ae_gmst().

Parameters
jd_ut1The Julian date in UT1.
jd_ttThe Julian date in TT. (See ae_delta_t().)
tlongThe longitude of the observer in degrees.
Returns
Local sidereal time in degrees.

References ae_gmst(), and ae_mod_360().

double ae_mjd ( double  jd)

Given a Julian Date, get the Modified Julian Date.

The difference between the two is simply a constant, AE_MJD_START.

Returns
The Modified Julian Date.

References AE_MJD_START.

Referenced by ae_dut1().

double ae_mod_180 ( double  x)

Reduce x modulo 360 degrees, but centre at 0 degrees.

Parameters
xThe number to be modulo'd.
Returns
A number in the range (-180: 180].

References ae_mod_360().

Referenced by ae_subobs_point().

double ae_mod_2pi ( double  x)

Reduce x modulo 2 pi.

Parameters
xThe number to be modulo'd.
Returns
The number in the range [0: 2 pi].

Referenced by ae_kepler().

double ae_mod_360 ( double  x)

Reduce x modulo 360 degrees.

Parameters
xThe number to be modulo'd.
Returns
A number in the range [0: 360.0].

Referenced by ae_altaz_to_radec(), ae_gal_to_radec(), ae_gast(), ae_gmst(), ae_kepler(), ae_last(), ae_lmst(), ae_mod_180(), ae_phys_pole(), ae_radec_to_gal(), and ae_subobs_point().

void ae_nutate ( double  nutl,
double  nuto,
double  eps,
double  p[],
int  direction 
)

Correct for nutation.

From AA page B20.

Parameters
nutlThe nutation in longitude in seconds of arc (see ae_nutation_lon_ob()).
nutoThe nutation in oblation in seconds of arc (see ae_nutation_lon_ob()).
epsThe obliquity of the ecliptic in seconds of arc (see ae_epsilon()).
pThe equitorial rectangular position vector for mean ecliptic and equinox of date. This is modified by the routine for nutation.
directionAE_TO_J2000 to nutate from jd_tt to J2000; AE_FROM_J2000 to nutate from J2000 to jd_tt. (If neither is passed, AE_TO_J2000 is assumed.) To nutate from jd_1 to jd_2, first go from jd_1 to J2000, and then from J2000 to jd_2.

References AE_FROM_J2000, and AE_STR.

Referenced by ae_geocentric_moon_from_orbit(), and aes_nutate().

void ae_nutation_lon_ob ( double  jd_tt,
double *  nutl,
double *  nuto 
)

Calculate the nutation in longitude and oblation.

References:

  • "Summary of 1980 IAU Theory of Nutation (Final Report of the IAU Working Group on Nutation)", P. K. Seidelmann et al., in Transactions of the IAU Vol. XVIII A, Reports on Astronomy, P. A. Wayman, ed.; D. Reidel Pub. Co., 1982.
  • "Nutation and the Earth's Rotation", I.A.U. Symposium No. 78, May, 1977, page 256. I.A.U., 1980.
  • Woolard, E.W., "A redevelopment of the theory of nutation", The Astronomical Journal, 58, 1-3 (1953).

This program implements all of the 1980 IAU nutation series. Results checked at 100 points against the 1986 AA; all agreed.

Parameters
jd_ttThe Julian date in TT.
nutlFor returning the nutation in longitude, in seconds of arc.
nutoFor returning the nutation in oblation, in seconds of arc.

References ae_mod_3600, ae_nt_param, ae_nut_cc, ae_nut_ss, ae_nut_sscc(), and AE_STR.

Referenced by ae_geocentric_moon_from_orbit(), ae_topocentric(), aes_last(), and aes_nutate().

void ae_phys_pole ( double  jd_ut1,
double  jd_tt,
const struct ae_physical_t phys,
double  n[],
double *  w 
)

Calculate the axis of a body's north pole and its prime meridean.

Get the precessed, nutated direction of the body's north pole; also return the prime meridean.

Parameters
jd_ut1The Julian date in UT1. For reasonable accuracy, pass a date that has been corrected for light travel time, to get the position of the north pole and meridean as seen by the observer.
jd_ttThe Julian date in TT. This should not be corrected for light-travel time.
physThe physical parameters of the body.
nFor returning the unit vector of the north pole, in equatorial units of date, in degrees.
wFor returning the prime meridean, in degrees.

References ae_physical_term_t::a, AE_DTR, AE_FROM_J2000, AE_J2000, ae_mod_360(), AE_PHYSICAL_D, AE_PHYSICAL_END, ae_polar_to_rect(), ae_precess(), aes_nutate(), ae_physical_term_t::b, ae_physical_term_t::c, ae_physical_term_t::d, ae_physical_t::dec_cos_term, ae_physical_t::pole_dec, ae_physical_t::pole_dec_t, ae_physical_t::pole_ra, ae_physical_t::pole_ra_t, ae_physical_t::ra_sin_term, ae_physical_term_t::time_var, ae_physical_t::w, ae_physical_t::w_d, ae_physical_t::w_d_sq, and ae_physical_t::w_sin_term.

Referenced by aes_subobs_point().

void ae_polar_to_rect ( double  ra,
double  dec,
double  radius,
double  rect[] 
)

Convert ra/dec/radius to a rectangular vector.

Parameters
raThe right ascension, in degrees.
decThe declination, in degrees.
radiusThe radius; the units of radius will determine the units of the output.
rectFor returning the rectangular vector.

References AE_DTR.

Referenced by ae_coord_to_constel_index(), ae_phys_pole(), ae_subobs_point(), and aes_subobs_point().

void ae_precess ( double  jd_tt,
double  r[],
int  direction 
)

Precess a coordinate.

The precession is calcuated using ae_precess_coef, ae_precess_node_coef and ae_precess_incl_coef.

Parameters
jd_ttThe Julian date in TT.
rA rectangular coordinate vector to be precessed. This routine returns the precessed coordinate in this parameter.
directionAE_TO_J2000 to precess from jd_tt to J2000; AE_FROM_J2000 to precess from J2000 to jd_tt. (If neither is passed, AE_TO_J2000 is assumed.) To precess from jd_1 to jd_2, first go from jd_1 to J2000, and then from J2000 to jd_2.

References ae_epsilon(), AE_FROM_J2000, AE_J2000, ae_precess_coef, ae_precess_incl_coef, ae_precess_node_coef, AE_STR, and AE_TO_J2000.

Referenced by ae_coord_to_constel_index(), ae_geocentric(), ae_geocentric_from_cat(), ae_kepler(), ae_kepler_embofs(), ae_phys_pole(), and ae_rect_to_polar_with_jd().

void ae_radec_to_altaz ( double  last,
double  glat,
double  ra,
double  dec,
double *  alt,
double *  az 
)

Convert ra/dec to alt/az.

Parameters
lastThe local aparent sidereal time, in degrees.
glatThe geodetic latitude, in degrees.
raThe right ascension, in degrees.
decThe declination, in degrees.
altFor returning the altitude, in degrees.
azFor returning the azimuth, in degrees.

References AE_DTR, AE_RTD, and ae_zatan2().

Referenced by ae_topocentric().

void ae_radec_to_gal ( double  ra,
double  dec,
double *  l,
double *  b,
int  fk5 
)

Convert ra/dec to l/b (galactic coordinates).

Parameters
raThe right ascension, in degrees.
decThe declination, in degrees.
lFor returning the galactic longitude, in degrees.
bFor returning the galactic latitude, in degrees.
fk5If non-zero, assume FK5 system (i.e., J2000 compatible); if zero, use the old FK4 (B1950 compatible) galactic pole.

References AE_DTR, ae_mod_360(), and AE_RTD.

int ae_read_orbit_from_cat ( const char *  path,
const char *  name,
struct ae_orbit_t orb 
)

Parse orbital elements from a catalogue file.

The catalogue file should be a columnated ASCII file. Lines starting with a hash ('#') are treated as comments. The columns are:

  • epoch of elements;
  • inclination, in degrees;
  • longitude of the ascending node, in degrees;
  • argument of the perihelion, in degrees;
  • mean distance, in AU;
  • daily motion, in AU per day;
  • eccentricity;
  • mean anomaly, in degrees;
  • epoch of equinox and ecliptic;
  • name.
Parameters
pathThe location of the catalogue file.
nameThe name of the object, as it appears in the catalogue file.
orbFor returning the orbital elements. If the routine fails, the information in el will be undefined.
Returns
0 on successful load; on failure, a negative error code from ae_retcode_t.

References ae_orbit_t::a, AE_RET_BAD_PATH, AE_RET_OBJ_NOT_FOUND, ae_orbit_t::dm, ae_orbit_t::ecc, ae_orbit_t::epoch, ae_orbit_t::equinox, ae_orbit_t::i, ae_orbit_t::L, ae_orbit_t::M, ae_orbit_t::name, ae_orbit_t::plat, ae_orbit_t::ptable, ae_orbit_t::r, ae_orbit_t::W, and ae_orbit_t::w.

int ae_read_star_from_cat ( const char *  path,
const char *  name,
struct ae_star_t star 
)

Parse star positional data from a catalogue file.

The catalogue file should be a columnated ASCII file. Lines starting with a hash ('#') are treated as comments. The columns are:

  • epoch
  • ra hours
  • ra minutes
  • ra second
  • dec degrees
  • dec minutes
  • dec seconds
  • mu ra (in seconds/century)
  • mu dec (in seconds of arc/century)
  • v (in km/s)
  • px (in seconds of arc)
  • mag
  • name
Parameters
pathThe location of the catalogue file.
nameThe name of the star, as it appears in the catalogue file.
starFor returning the star position data. If the routine fails, the information in el will be undefined.
Returns
0 on successful load; on failure, a negative error code from ae_retcode_t.

References AE_B1950, AE_J1900, AE_J2000, AE_RET_BAD_PATH, AE_RET_OBJ_NOT_FOUND, ae_star_t::dec, ae_star_t::epoch, ae_star_t::mag, ae_star_t::mudec, ae_star_t::mura, ae_star_t::name, ae_star_t::px, ae_star_t::ra, and ae_star_t::v.

void ae_rect_to_polar ( const double  rect[],
double *  ra,
double *  dec,
double *  radius 
)

Convert an equatorial rectangular unit vector to ra/dec/radius.

Parameters
rectAn equatorial rectangular vector.
raFor returning the right ascension, in degrees.
decFor returning the declination, in degrees.
radiusFor returning the radius. Set to NULL if you do not need this value.

References AE_RTD, and ae_zatan2().

Referenced by ae_geocentric(), ae_geocentric_from_cat(), ae_geocentric_moon_from_orbit(), and ae_subobs_point().

void ae_rect_to_polar_with_jd ( double  pp[],
double  jd,
char  ofdate,
double  polar[] 
)

Convert equatorial rectangular coordinates to polar coordinates.

This routine is similar to ae_rect_to_polar(), except that it also performs //! any necessary precession.

Parameters
ppThe equatorial rectangular coordinates.
jdThe Julian date.
ofdateIf 1, precess from J2000.0 to jd.
polarFor returning the polar coordinates.

References ae_epsilon(), AE_FROM_J2000, ae_precess(), AE_STR, and ae_zatan2().

double ae_refrac_ulich ( double  alt,
void *  param 
)

Atmospheric refraction correction (in millimetre bands).

This function is suitable for passing to ae_topocentric().

This function uses the equations from Ulich 1981. It does not cite which millimetre wavelength ranges it is expected to be good for. The ALMA memo 366 does analysis on another model and finds it good to 2% up to frequencies of 1000 GHz. (The worst is at 500 GHz.)

Ulich claims accuracy to 2" above altitudes of 3 degrees.

References: -B. L. Ulich. "Millimeter Wave Radio Telescopes: Gain and Pointing Charactersitics". International Journal of Infrared and Millimeter Waves, 2, 2 (1981).

Parameters
altThe altitude of the observation, in degrees.
paramFor passing the atmospheric pressure, in hPa, the temperature, in degrees centigrade and the surface humidity, as a percentage. These should be stored as doubles, such that: -pressure is in ((double *) param)[0], -temperature is in ((double *) param)[1]. -and humidity is in ((double *) param)[2].
Returns
The correction in degrees to be added to true altitude to obtain apparent altitude.

References AE_DTR, and AE_STD.

double ae_refrac_visible ( double  alt,
void *  param 
)

Atmospheric refraction correction (in visible bands).

This function is suitable for passing to ae_topocentric().

For high altitude angle, AA page B61 is used. The accuracy is `usually about 0.1 arcsecond'.

The formula for low altitude is from the Almanac for Computers. It gives the correction for observed altitude, so has to be inverted numerically to get the observed from the true. The accuracy is about 0.2' for -20C < T < +40C and 970mb < P < 1050mb.

Parameters
altThe altitude of the observation, in degrees.
paramFor passing the atmospheric pressure, in millibar, and the temperature, in degrees centigrade. These should be stored as doubles, such that the pressure is in ((double *) param)[0] and the temperature in ((double *) param)[1].
Returns
The correction in degrees to be added to true altitude to obtain apparent altitude.

References AE_DTR.

void ae_relativity ( double  p[],
double  q[],
double  o[] 
)

Correct for light deflection due to solar gravitation.

Note that this blows up if the object is very near (or is) the sun. Therefore, if the sun-object or sun-object distance is zero, no correction is made.

See AA page B37.

Parameters
pThe unit vector from observer to an object; this routine returns the corrected vector in this parameter.
qThe heliocentric ecliptic rectangular coordinates of the object.
oThe heliocentric ecliptic rectangular coordinates of the observer.

Referenced by ae_geocentric(), and ae_geocentric_from_cat().

const char* ae_retcode_string ( int  retcode)

For returning human readable return codes.

Parameters
retcodeA return code. The absolute value will be used, so either a negative or positive number may be passed.
Returns
A string describing the return code.

References ae_retcode_name.

void ae_subobs_point ( double  j[],
double  n[],
double  w,
double  f,
int  retrograde,
double *  lat,
double *  lon 
)

Calculate the sub-observer point.

The sub-observer point is the latitude and longitude on an observed body at the observer is directly over-head.

Note that all input values should be precessed and nutated to the epoch of observation.

WARNING: the longitude computed by this function is, for some reason, only accurate to tens of arc minutes when compared to HORIZONS (and I don't know why. Latitude is accurate to a fraction of an arc minute.

Parameters
jThe unit vector pointing from the centre of the observed body to the observer.
nThe unit vector of the north pole of the observed body, in equatorial coordinates of date.
wThe prime meridean of the observed body, in degreees.
fThe flattening of the observed body.
retrogradePass non-zero value if the rotation is retrograde, i.e., the body rotates in the opposite sense of the earth's rotation. Pass 0 for prograde rotation.
latFor returning the sub-observer latitude, in degrees.
lonFor returning the sub-observer longitude, in degrees.

References AE_DTR, ae_mod_180(), ae_mod_360(), ae_polar_to_rect(), ae_rect_to_polar(), AE_RTD, ae_zatan2(), ae_physical_t::pole_dec, and ae_physical_t::pole_ra.

Referenced by aes_subobs_point().

double ae_tdb ( double  jd)

Find Barycentric Dynamical Time from Terrestrial Dynamical Time.

See AA page B5.

Parameters
jdA Julian date, in TDT.
Returns
The corresponding time in TDB.

References AE_J2000, and AE_STR.

void ae_topocentric ( double  jd_tt,
double  jd_ut1,
double  tlat,
double  glat,
double  tlong,
double  trho,
double(*)(double, void *)  refrac,
void *  refrac_param,
double  dist,
double *  ra,
double *  dec 
)

Apply diurnal aberrations and calculate topocentric altitude and azimuth.

Apply diurnal aberrations and calculate topocentric altitude and azimuth, given the geocentric apparent right ascension and declination. From AA page B60 and D3.

The geocentric coordinates can be obtained from ae_geocentric_from_orbit() or a similar routine.

Parameters
jd_ttThe Julian date in TT.
jd_ut1The Julian date in UT1 measure.
tlatThe geocentric latitude of the observer, in degrees.
glatThe geodetic latitude of the observer, in degrees.
tlongThe longitude of the observer, in degrees.
trhoThe distance from the centre of the earth to the observer, in earth radii.
refracA routine that calculates refraction (pass NULL to perform no refraction correction). Its first argument is the altitude, in degrees. Its second argument is a void pointer for passing other parameters, which will be particular to the refraction algorithm. These parameters are passed via refrac_param.
refrac_paramFor passing auxiliary parameters to the refrac function. Pass NULL if refrac is not being used.
distThe geocentric distance to the object being observed, in AU. Pass 0 for extra-solar-system objects.
raThe right ascension of the object, in degrees; this routine returns the topocentric right ascension in this parameter.
decThe declination of the object, in degrees; this routine returns the topocentric declination in this parameter.

References AE_ADD_ABERRATION, ae_altaz_to_radec(), ae_diurnal_aberration(), ae_diurnal_parallax(), ae_epsilon(), ae_last(), ae_nutation_lon_ob(), and ae_radec_to_altaz().

Referenced by aes_topocentric().

void ae_v_orbit ( double  jd_tt,
const struct ae_orbit_t orb,
double  v[] 
)

Calculate the velocity of an object as it moves in its orbit.

The difference in velocity between this and assuming a circular orbit is only about 1 part in 10**4. Note that this gives heliocentric, not barycentric, velocity.

Parameters
jd_ttThe Julian date in TT.
orbOrbital elements for the object.
vFor returning the heliocentric, rectangular velocity, in AU per day.

References ae_kepler().

Referenced by ae_geocentric_from_orbit(), ae_geocentric_sun_from_orbit(), and aes_subobs_point().

double ae_zatan2 ( double  x,
double  y 
)

Quadrant correct inverse circular tangent.

Returns radian angle between 0 and +2PI whose tangent is y/x.

Cephes Math Library Release 2.0: April, 1987 Copyright 1984, 1987 by Stephen L. Moshier. Direct inquiries to 30 Frost Street, Cambridge, MA 02140 Certain routines from the Library, including this one, may be used and distributed freely provided this notice is retained and source code is included with all distributions.

Parameters
xThe denominator.
yThe numerator.
Returns
The quadrant correct inverse circular tangent.

Referenced by ae_altaz_to_radec(), ae_delta_q(), ae_diurnal_parallax(), ae_fk4_to_fk5(), ae_kepler(), ae_radec_to_altaz(), ae_rect_to_polar(), ae_rect_to_polar_with_jd(), and ae_subobs_point().

double aes_ctime_to_last ( double  t,
double  tlong 
)

Calculate the local sidereal time from a C time.

This is a simplification of ae_ctime_to_last(). It automatically does the calculations of delta_t, nutl and eps. Use ae_ctime_to_last() for more efficiency, especially with repeated calls.

Parameters
tThe UNIX time at which to calculate the LST.
tlongThe longitude of the observer, in degrees.
Returns
The LAST in degrees.

References ae_ctime_to_jd(), AE_D_PER_S, ae_dut1(), and aes_last().

double aes_disc_solid_angle ( double  jd_ut1,
const struct ae_orbit_t o_orb,
const struct ae_orbit_t q_orb,
const struct ae_physical_t phys 
)

Calculate the solid angle of an observed disc.

This simplified routine is potentially less efficient, since the distance to the planet is calculated as well as the sub-earth latitude to find the semi-minor axis of the projected disc—ae_disc_solid_angle() and ae_disc_semiminor() can make use of precomputed values.

Parameters
jd_ut1The Julian date in UT1.
o_orbThe orbital elements of the observer.
q_orbThe orbital elements of the object being observed.
physThe physical elements of the object being observed.
Returns
The solid angle of the disc, in steradians.

References ae_disc_semiminor(), ae_disc_solid_angle(), aes_subobs_point(), ae_physical_t::r_eq, and ae_physical_t::r_pole.

double aes_last ( double  jd_ut1,
double  tlong 
)

Get the local apparent sidereal time.

This is a simplified version of ae_last(), which automatically calculates jd_tt, nutl and eps. It might be less efficient to use this routine if speed is required.

Parameters
jd_ut1The Julian date in UT1.
tlongThe longitude of the observer in degrees.
Returns
Local sidereal time in degrees.

References AE_D_PER_S, ae_delta_t(), ae_epsilon(), ae_last(), and ae_nutation_lon_ob().

Referenced by aes_ctime_to_last().

void aes_nutate ( double  jd_tt,
double  p[],
int  direction 
)

Correct for nutation.

This is a simplified version of ae_nutate(), which automatically calculates nutl, nuto and eps. It might be less efficient to use this routine if speed is required.

Parameters
jd_ttThe Julian date in TT.
pThe equitorial rectangular position vector for mean ecliptic and equinox of date. This is modified by the routine for nutation.
directionAE_TO_J2000 to nutate from jd_tt to J2000; AE_FROM_J2000 to nutate from J2000 to jd_tt. (If neither is passed, AE_TO_J2000 is assumed.) To nutate from jd_1 to jd_2, first go from jd_1 to J2000, and then from J2000 to jd_2.

References ae_epsilon(), ae_nutate(), and ae_nutation_lon_ob().

Referenced by ae_geocentric(), ae_geocentric_from_cat(), and ae_phys_pole().

void aes_subobs_point ( double  jd_ut1,
const struct ae_orbit_t o_orb,
const struct ae_orbit_t q_orb,
const struct ae_physical_t phys,
double *  lat,
double *  lon,
double *  dist 
)

Calculate the sub-observer point.

This is a simplified version of ae_subobs_point(). It automatically calculates the requisite unit vectors, given the orbital and physical elements of the bodies.

The sub-observer point is the latitude and longitude on an observed body at the observer is directly over-head.

This routine includes corrections for light-travel time, annual aberration, precession and nutation.

WARNING: the longitude computed by this function is, for some reason, only accurate to tens of arc minutes when compared to HORIZONS (and I don't know why. Latitude is accurate to a fraction of an arc minute.

Parameters
jd_ut1The Julian date in UT1.
o_orbThe orbital elements of the observer.
q_orbThe orbital elements of the object being observed.
physThe physical elements of the object being observed.
latFor returning the sub-observer latitude, in degrees.
lonFor returning the sub-observer longitude, in degress.
distFor returning the distance to the object, in AU. Set to NULL if you do not require this value. Note that this returns the distance between the centres of the two bodies.

References AE_D_PER_S, ae_delta_t(), ae_flattening(), ae_geocentric(), ae_is_retrograde(), ae_kepler(), ae_light_t(), ae_phys_pole(), ae_polar_to_rect(), ae_subobs_point(), ae_v_orbit(), and ae_physical_t::w.

Referenced by aes_disc_solid_angle().

void aes_topocentric ( double  jd_ut1,
double  glat,
double  tlong,
double  dist,
double *  ra,
double *  dec 
)

Apply diurnal aberrations and calculate topocentric altitude and azimuth.

This is a simplified version of ae_topocentric(), which automatically calculates jd_tt, tlat and trho. It assumes that the altitude of the observer is at sea level. It also uses no atmospheric model. It might be less efficient to use this routine if speed is required.

Parameters
jd_ut1The Julian date in UT1 measure.
glatThe geodetic latitude of the observer, in degrees.
tlongThe longitude of the observer, in degrees.
distThe geocentric distance to the object being observed, in AU. Pass 0 for extra-solar-system objects.
raThe right ascension of the object, in degrees; this routine returns the topocentric right ascension in this parameter.
decThe declination of the object, in degrees; this routine returns the topocentric declination in this parameter.

References ae_delta_t(), ae_geocentric_lat(), AE_S_PER_D, and ae_topocentric().

Variable Documentation

const char* ae_constel_name[]

A list of constellation names.

Cancri is also abbreviated `Cnc'.

Referenced by ae_cat_to_constel_index().

struct ae_plantbl_t ae_mlat404

Ephemeris for the moon.

These are with respect to the mean ecliptic and equinox of date.

Referenced by ae_gmoon().

struct ae_plantbl_t ae_mlr404

Ephemeris data for the moon.

These are with respect to the mean equinox and ecliptic of date.

Referenced by ae_gmoon().

struct ae_orbit_t* ae_orb_planet[]

Pointers to the planetary orbit objects.

The indices correspond to ae_ss_bodies_t.

struct ae_physical_t* ae_phys_planet[]

Pointers to the planetary physical information objects.

The indices correspond to ae_ss_bodies_t.

const char* ae_ss_name[]

The names of the solar system bodies.

The indices correspond to ae_ss_bodies_t.