![]() |
SuperNOVAS v1.5
The NOVAS C library, made better
|
SuperNOVAS types, definitions, and function prototypes. More...
Data Structures | |
struct | cat_entry |
Basic astrometric data for any sidereal object located outside the solar system. More... | |
struct | in_space |
data for an observer's location on Earth orbit More... | |
struct | novas_delaunay_args |
Fundamental Delaunay arguments of the Sun and Moon, from Simon section 3.4(b.3). More... | |
struct | novas_frame |
A set of parameters that uniquely define the place and time of observation. More... | |
struct | novas_matrix |
A 3x3 matrix for coordinate transformations. More... | |
struct | novas_observable |
Spherical and spectral coordinate set. More... | |
struct | novas_orbital |
Keplerian orbital elements for NOVAS_ORBITAL_OBJECT type. More... | |
struct | novas_orbital_system |
Specification of an orbital system, in which orbital elements are defined. More... | |
struct | novas_planet_bundle |
Position and velocity data for a set of major planets (which may include the Sun and the Moon also). More... | |
struct | novas_timespec |
A structure, which defines a precise instant of time that can be extpressed in any of the astronomical timescales. More... | |
struct | novas_track |
The spherical and spectral tracking position of a source, and its first and second time derivatives. More... | |
struct | novas_transform |
A transformation between two astronomical coordinate systems for the same observer location and time. More... | |
struct | object |
Celestial object of interest. More... | |
struct | observer |
Observer location. More... | |
struct | on_surface |
Data for an observer's location on the surface of the Earth, and optional local weather data for refraction calculations only. More... | |
struct | ra_of_cio |
[deprecated] More... | |
struct | sky_pos |
Celestial object's place on the sky; contains the output from place() More... | |
Macros | |
#define | ASEC2RAD (DEG2RAD / 3600.0) |
[rad/arcsec] 1 arcsecond in radians | |
#define | ASEC360 (360 * 60 * 60) |
[arcsec] Number of arcseconds in 360 degrees. | |
#define | BARYC NOVAS_BARYCENTER |
Deprecated. | |
#define | CAT_ENTRY_INIT |
Initializer for a NOVAS cat_entry structure. | |
#define | DE405_AU 1.4959787069098932e+11 |
[m] Astronomical unit (AU). | |
#define | DEG2RAD (M_PI / 180.0) |
[rad/deg] 1 degree in radians | |
#define | HELIOC NOVAS_HELIOCENTER |
Deprecated. | |
#define | IN_SPACE_INIT |
Initializer for a NOVAS in_space structure. | |
#define | M_PI 3.14159265358979323846 |
Definition of π in case it's not defined in math.h. | |
#define | NOVAS_ARCMIN (NOVAS_DEGREE / 60.0) |
[rad] An arc minute expressed in radians. | |
#define | NOVAS_ARCSEC (NOVAS_ARCMIN / 60.0) |
[rad] An arc second expressed in radians. | |
#define | NOVAS_AU 1.495978707e+11 |
[m] Astronomical unit (AU). | |
#define | NOVAS_AU_KM ( 1e-3 * NOVAS_AU ) |
[km] Astronomical Unit (AU) in kilometers. | |
#define | NOVAS_AU_SEC ( NOVAS_AU / NOVAS_C ) |
[s] Light-time for one astronomical unit (AU) in seconds. | |
#define | NOVAS_C 299792458.0 |
[m/s] Speed of light in meters/second is a defining physical constant. | |
#define | NOVAS_C_AU_PER_DAY ( NOVAS_DAY / AU_SEC ) |
[AU/day] Speed of light in units of AU/day. | |
#define | NOVAS_DAY 86400.0 |
[s] The length of a synodic day, that is 24 hours exactly. | |
#define | NOVAS_DEFAULT_DISTANCE (1e9 * NOVAS_PARSEC) |
[m] Default distance at which sidereal sources are assumed when not specified otherwise Historically NOVAS placed sources with no parallax at 1 Gpc distance, and so we follow. | |
#define | NOVAS_DEFAULT_WAVELENGTH 0.55 |
[μm] Default wavelength, e.g. | |
#define | NOVAS_DEGREE (M_PI / 180.0) |
[rad] A degree expressed in radians. | |
#define | NOVAS_DELAUNAY_ARGS_INIT |
Empty initializer for novas_delaunay_args. | |
#define | NOVAS_EARTH_ANGVEL 7.2921150e-5 |
[rad/s] Rotational angular velocity of Earth from IERS Conventions (2003). | |
#define | NOVAS_EARTH_FLATTENING NOVAS_GRS80_FLATTENING |
[m] Earth ellipsoid flattening (ITRF / GRS80_MODEL) | |
#define | NOVAS_EARTH_INIT |
object initializer for the planet Earth | |
#define | NOVAS_EARTH_RADIUS NOVAS_GRS80_RADIUS |
[m] Equatorial radius of Earth (ITRF / GRS80 model) | |
#define | NOVAS_EMB_INIT |
object initializer for the the Earth-Moon Barycenter (EMB) | |
#define | NOVAS_EQUATOR_TYPES (NOVAS_GCRS_EQUATOR + 1) |
The number of equator types defined in enum novas_equator_type . | |
#define | NOVAS_FRAME_INIT |
Empty initializer for novas_frame. | |
#define | NOVAS_G_EARTH 3.98600435507e14 |
[m3/s2] Geocentric gravitational constant (GMearth) from DE440, see Park et al., AJ, 161, 105 (2021) | |
#define | NOVAS_G_SUN 1.32712440041279419e20 |
[m3/s2] Heliocentric gravitational constant (GMsun) from DE440, see Park et al., AJ, 161, 105 (2021) | |
#define | NOVAS_GPS_TO_TAI 19.0 |
[s] TAI - GPS time offset | |
#define | NOVAS_GRS80_FLATTENING (1.0 / 298.257222101) |
[m] WGS84 Earth flattening | |
#define | NOVAS_GRS80_RADIUS 6378137.0 |
[m] Equatorial radius of the WGS84 reference ellipsoid. | |
#define | NOVAS_HOURANGLE (M_PI / 12.0) |
[rad] An hour of angle expressed in radians. | |
#define | NOVAS_IERS_EARTH_FLATTENING (1.0 / 298.25642) |
[m] Earth ellipsoid flattening from IERS Conventions (2003). | |
#define | NOVAS_IERS_EARTH_RADIUS 6378136.6 |
[m] Equatorial radius of Earth in meters from IERS Conventions (2003). | |
#define | NOVAS_JD_B1900 15019.81352 |
[day] Julian date at B1900 precession(), transform_cat() | |
#define | NOVAS_JD_B1950 2433282.42345905 |
[day] Julian date at B1950 precession(), transform_cat() | |
#define | NOVAS_JD_HIP 2448349.0625 |
[day] Julian date for J1991.25, which the Hipparcos catalog is referred to. | |
#define | NOVAS_JD_J2000 2451545.0 |
[day] Julian date at J2000 | |
#define | NOVAS_JD_MJD0 2400000.5 |
[day] Julian date at which the Modified Julian Date (MJD) is zero | |
#define | NOVAS_JD_START_GREGORIAN 2299160.5 |
[day] The Julian day the Gregorian calendar was introduced in 15 October 1582. | |
#define | NOVAS_JULIAN_YEAR_DAYS 365.25 |
[day] The length of a tropical year (at J2000) in days. | |
#define | NOVAS_JUPITER_INIT |
object initializer for the planet Jupiter | |
#define | NOVAS_KM 1000.0 |
[m] A kilometer (km) in meters. | |
#define | NOVAS_KMS (NOVAS_KM) |
[m] One km/s in m/s. | |
#define | NOVAS_LIGHT_YEAR ( NOVAS_C * NOVAS_TROPICAL_YEAR_DAYS * NOVAS_DAY ) |
[m] A light-year in meters. | |
#define | NOVAS_MAJOR_VERSION 3 |
Major version of NOVAS on which this library is based. | |
#define | NOVAS_MARS_INIT |
object initializer for the planet Mars | |
#define | NOVAS_MATRIX_IDENTITY |
novas_matric initializer for an indentity matrix | |
#define | NOVAS_MATRIX_INIT |
Empty initializer for novas_matrix. | |
#define | NOVAS_MERCURY_INIT |
object initializer for the planet Venus | |
#define | NOVAS_MINOR_VERSION 1 |
Minor version of NOVAS on which this library is based. | |
#define | NOVAS_MOON_INIT |
object initializer for the Moon | |
#define | NOVAS_NEPTUNE_INIT |
object initializer for the planet Neptune | |
#define | NOVAS_OBJECT_INIT |
Empty object initializer. | |
#define | NOVAS_OBJECT_TYPES (NOVAS_ORBITAL_OBJECT + 1) |
The number of object types distinguished by NOVAS. | |
#define | NOVAS_OBSERVABLE_INIT |
Empty initializer for novas_observable. | |
#define | NOVAS_OBSERVER_PLACES (NOVAS_SOLAR_SYSTEM_OBSERVER + 1) |
The number of observer place types supported. | |
#define | NOVAS_ORBIT_INIT |
Initializer for novas_orbital for heliocentric orbits using GCRS ecliptic pqrametrization. | |
#define | NOVAS_ORBITAL_SYSTEM_INIT |
Default orbital system initializer for heliocentric GCRS ecliptic orbits. | |
#define | NOVAS_ORIGIN_TYPES (NOVAS_HELIOCENTER + 1) |
The number of different ICSR origins available in NOVAS. | |
#define | NOVAS_PARSEC ( NOVAS_AU / NOVAS_ARCSEC ) |
[m] A parsec in meters. | |
#define | NOVAS_PLANET_BUNDLE_INIT |
Empty initializer for novas_planet_bundle. | |
#define | NOVAS_PLANET_INIT(num, name) |
object initializer macro for major planets, the Sun, Moon, and barycenters. | |
#define | NOVAS_PLANETS (NOVAS_PLUTO_BARYCENTER + 1) |
The number of major planets defined in NOVAS. | |
#define | NOVAS_PLUTO_BARYCENTER_INIT |
object initializer for the Pluto system barycenter | |
#define | NOVAS_PLUTO_INIT |
object initializer for the minor planet Pluto | |
#define | NOVAS_REFERENCE_SYSTEMS (NOVAS_ITRS + 1) |
The number of basic coordinate reference systems in NOVAS. | |
#define | NOVAS_REFRACTION_MODELS (NOVAS_WAVE_REFRACTION + 1) |
The number of built-in refraction models available in SuperNOVAS. | |
#define | NOVAS_SATURN_INIT |
object initializer for the planet Saturn | |
#define | NOVAS_SOLAR_RADIUS 696340000.0 |
[m] Solar radius (photosphere) | |
#define | NOVAS_SSB_INIT |
object initializer for the Solar System Barycenter (SSB) | |
#define | NOVAS_SUN_INIT |
object initializer for the Sun | |
#define | NOVAS_SYSTEM_B1950 "B1950" |
The B1950 coordiante system as a string. | |
#define | NOVAS_SYSTEM_FK4 "FK4" |
The 4th catalog (FK4) coordinate system as a string. | |
#define | NOVAS_SYSTEM_FK5 "FK5" |
The 5th catalog (FK5) coordinate system as a string. | |
#define | NOVAS_SYSTEM_HIP "HIP" |
The Hipparcos dataset coordinate system as a string. | |
#define | NOVAS_SYSTEM_ICRS "ICRS" |
The ICRS system as a string. | |
#define | NOVAS_SYSTEM_J2000 "J2000" |
The J2000 coordinate syste, as a string. | |
#define | NOVAS_TAI_TO_TT 32.184 |
[s] TT - TAI time offset | |
#define | NOVAS_TIMESCALES (NOVAS_UT1 + 1) |
The number of asronomical time scales supported. | |
#define | NOVAS_TIMESPEC_INIT |
Empty initializer for novas_timespec. | |
#define | NOVAS_TRACK_INIT |
Empty initializer for novas_track. | |
#define | NOVAS_TRANSFORM_INIT |
Empty initializer for NOVAS_TRANSFORM. | |
#define | NOVAS_TRANSFORM_TYPES (ICRS_TO_J2000 + 1) |
The number of coordinate transfor types in NOVAS. | |
#define | NOVAS_TROPICAL_YEAR_DAYS 365.2421897 |
[day] The length of a tropical year (at J2000) in days. | |
#define | NOVAS_URANUS_INIT |
object initializer for the planet Uranus | |
#define | NOVAS_VENUS_INIT |
object initializer for the planet Mercury | |
#define | NOVAS_VERSION_STRING str_2(NOVAS_MAJOR_VERSION) "." str_2(NOVAS_MINOR_VERSION) |
The version string of the upstream NOVAS library on which this library is based. | |
#define | NOVAS_WGS84_FLATTENING (1.0 / 298.257223563) |
[m] WGS84 Earth flattening | |
#define | NOVAS_WGS84_RADIUS 6378137.0 |
[m] Equatorial radius of the WGS84 reference ellipsoid. | |
#define | NOVAS_WOBBLE_DIRECTIONS (WOBBLE_PEF_TO_ITRS + 1) |
Number of values in enum novas_wobble_direction. | |
#define | OBSERVER_INIT |
Empty initializer for observer. | |
#define | ON_SURFACE_INIT |
Initializer for a NOVAS on_surface data structure. | |
#define | ON_SURFACE_LOC(lon, lat, alt) |
Initializer for a NOVAS on_surface data structure at a specified geodetic location. | |
#define | RAD2DEG (1.0 / DEG2RAD) |
[deg/rad] 1 radian in degrees | |
#define | SIZE_OF_CAT_NAME 6 |
Maximum bytes in catalog IDs including string termination. | |
#define | SIZE_OF_OBJ_NAME 50 |
Maximum bytes in object names including string termination. | |
#define | SKY_POS_INIT |
Initializer for a NOVAS sky_pos structure. | |
#define | SUPERNOVAS_MAJOR_VERSION 1 |
API major version. | |
#define | SUPERNOVAS_MINOR_VERSION 5 |
API minor version. | |
#define | SUPERNOVAS_PATCHLEVEL 0 |
Integer sub version of the release. | |
#define | SUPERNOVAS_RELEASE_STRING "-rc6" |
Additional release information in version, e.g. "-1", or "-rc1", or empty string "" for releases. | |
#define | SUPERNOVAS_VERSION_STRING |
The version string for this library. | |
#define | TWOPI (2.0 * M_PI) |
2π | |
Typedefs | |
typedef int(* | novas_nutation_provider) (double jd_tt_high, double jd_tt_low, double *restrict dpsi, double *restrict deps) |
Function type definition for the IAU 2000 nutation series calculation. | |
typedef double(* | RefractionModel) (double jd_tt, const on_surface *loc, enum novas_refraction_type type, double el) |
A function that returns a refraction correction for a given date/time of observation at the given site on earth, and for a given astrometric source elevation. | |
Enumerations | |
enum | novas_accuracy { NOVAS_FULL_ACCURACY = 0 , NOVAS_REDUCED_ACCURACY } |
Constants to control the precision of NOVAS nutation calculations. More... | |
enum | novas_calendar_type { NOVAS_ROMAN_CALENDAR = -1 , NOVAS_ASTRONOMICAL_CALENDAR , NOVAS_GREGORIAN_CALENDAR } |
Constants to disambiguate which type of calendar yo use for interpreting calendar dates. More... | |
enum | novas_cio_location_type { CIO_VS_GCRS = 1 , CIO_VS_EQUINOX } |
enum | novas_date_format { NOVAS_YMD = 0 , NOVAS_DMY , NOVAS_MDY } |
The general order of date components for parsing. More... | |
enum | novas_debug_mode { NOVAS_DEBUG_OFF = 0 , NOVAS_DEBUG_ON , NOVAS_DEBUG_EXTRA } |
Settings for 'novas_debug()'. More... | |
enum | novas_dynamical_type { NOVAS_DYNAMICAL_MOD = 0 , NOVAS_DYNAMICAL_TOD , NOVAS_DYNAMICAL_CIRS } |
Constants that determine the type of dynamical system. More... | |
enum | novas_earth_rotation_measure { EROT_ERA = 0 , EROT_GST } |
enum | novas_equator_type { NOVAS_MEAN_EQUATOR = 0 , NOVAS_TRUE_EQUATOR , NOVAS_GCRS_EQUATOR } |
Constants that determine the type of equator to be used for the coordinate system. More... | |
enum | novas_equatorial_class { NOVAS_REFERENCE_CLASS = 0 , NOVAS_DYNAMICAL_CLASS } |
enum | novas_equinox_type { NOVAS_MEAN_EQUINOX = 0 , NOVAS_TRUE_EQUINOX } |
The type of equinox used in the pre IAU 2006 (old) methodology. More... | |
enum | novas_frametie_direction { J2000_TO_ICRS = -1 , ICRS_TO_J2000 } |
Direction constant to use for frame_tie(), to determine the direction of transformation between J2000 and ICRS coordinates. More... | |
enum | novas_nutation_direction { NUTATE_TRUE_TO_MEAN = -1 , NUTATE_MEAN_TO_TRUE } |
Direction constant for nutation(), between mean and true equatorial coordinates. More... | |
enum | novas_object_type { NOVAS_PLANET = 0 , NOVAS_EPHEM_OBJECT , NOVAS_CATALOG_OBJECT , NOVAS_ORBITAL_OBJECT } |
The type of astronomical objects distinguied by the NOVAS library. More... | |
enum | novas_observer_place { NOVAS_OBSERVER_AT_GEOCENTER = 0 , NOVAS_OBSERVER_ON_EARTH , NOVAS_OBSERVER_IN_EARTH_ORBIT , NOVAS_AIRBORNE_OBSERVER , NOVAS_SOLAR_SYSTEM_OBSERVER } |
Types of places on and around Earth that may serve a a reference position for the observation. More... | |
enum | novas_origin { NOVAS_BARYCENTER = 0 , NOVAS_HELIOCENTER } |
The origin of the ICRS system for referencing positions and velocities for solar-system bodies. More... | |
enum | novas_planet { NOVAS_SSB = 0 , NOVAS_MERCURY , NOVAS_VENUS , NOVAS_EARTH , NOVAS_MARS , NOVAS_JUPITER , NOVAS_SATURN , NOVAS_URANUS , NOVAS_NEPTUNE , NOVAS_PLUTO , NOVAS_SUN , NOVAS_MOON , NOVAS_EMB , NOVAS_PLUTO_BARYCENTER } |
Enumeration for the 'major planet' numbers in NOVAS to use as the solar-system body number whenever the object type is NOVAS_PLANET. More... | |
enum | novas_pole_offset_type { POLE_OFFSETS_DPSI_DEPS = 1 , POLE_OFFSETS_X_Y } |
enum | novas_reference_ellipsoid { NOVAS_GRS80_ELLIPSOID = 0 , NOVAS_WGS84_ELLIPSOID , NOVAS_IERS_1989_ELLIPSOID , NOVAS_IERS_2003_ELLIPSOID } |
Type of Earth reference ellipsoid. More... | |
enum | novas_reference_plane { NOVAS_ECLIPTIC_PLANE = 0 , NOVAS_EQUATORIAL_PLANE } |
The plane in which values, such as orbital parameters are referenced. More... | |
enum | novas_reference_system { NOVAS_GCRS = 0 , NOVAS_TOD , NOVAS_CIRS , NOVAS_ICRS , NOVAS_J2000 , NOVAS_MOD , NOVAS_TIRS , NOVAS_ITRS } |
The basic types of positional coordinate reference systems supported by NOVAS. More... | |
enum | novas_refraction_model { NOVAS_NO_ATMOSPHERE = 0 , NOVAS_STANDARD_ATMOSPHERE , NOVAS_WEATHER_AT_LOCATION , NOVAS_RADIO_REFRACTION , NOVAS_WAVE_REFRACTION } |
Constants that determine whether what model (if any) to use for implicit refraction calculations. More... | |
enum | novas_refraction_type { NOVAS_REFRACT_OBSERVED = -1 , NOVAS_REFRACT_ASTROMETRIC } |
The type of elevation value for which to calculate a refraction. More... | |
enum | novas_separator_type { NOVAS_SEP_COLONS = 0 , NOVAS_SEP_SPACES , NOVAS_SEP_UNITS , NOVAS_SEP_UNITS_AND_SPACES } |
Separator type to use for broken-down time/angle string representations in HMS/DMS formats. More... | |
enum | novas_timescale { NOVAS_TCB = 0 , NOVAS_TDB , NOVAS_TCG , NOVAS_TT , NOVAS_TAI , NOVAS_GPS , NOVAS_UTC , NOVAS_UT1 } |
Constants to reference various astrnomical timescales used. More... | |
enum | novas_transform_type { PROPER_MOTION = 1 , PRECESSION , CHANGE_EPOCH , CHANGE_J2000_TO_ICRS , CHANGE_ICRS_TO_J2000 } |
The types of coordinate transformations available for tranform_cat(). More... | |
enum | novas_wobble_direction { WOBBLE_ITRS_TO_TIRS = 0 , WOBBLE_TIRS_TO_ITRS , WOBBLE_ITRS_TO_PEF , WOBBLE_PEF_TO_ITRS } |
Direction constants for polar wobble corrections via the wobble() function. More... | |
Functions | |
int | aberration (const double *pos, const double *vobs, double lighttime, double *out) |
Corrects position vector for aberration of light. | |
double | accum_prec (double t) |
Returns the general precession in longitude (Simon et al. | |
short | app_planet (double jd_tt, const object *restrict ss_body, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis) |
reference_system: TOD observer location: geocenter position_type: apparent | |
short | app_star (double jd_tt, const cat_entry *restrict star, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec) |
reference_system: TOD observer location: geocenter position_type: apparent | |
double | app_to_cirs_ra (double jd_tt, enum novas_accuracy accuracy, double ra) |
Converts an apparent right ascension coordinate (measured from the true equinox of date) to a CIRS R.A., measured from the CIO. | |
short | astro_planet (double jd_tt, const object *restrict ss_body, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis) |
reference_system: ICRS / GCRS observer location: geocenter position_type: geometric | |
short | astro_star (double jd_tt, const cat_entry *restrict star, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec) |
reference_system: ICRS / GCRS observer location: geocenter position_type: geometric | |
int | bary2obs (const double *pos, const double *pos_obs, double *out, double *restrict lighttime) |
Moves the origin of coordinates from the barycenter of the solar system to the observer (or the geocenter); i.e., this function accounts for parallax (annual+geocentric or just annual). | |
int | cal_date (double tjd, short *restrict year, short *restrict month, short *restrict day, double *restrict hour) |
This function will compute a broken down date on the astronomical calendar given the Julian day input. | |
short | cel2ter (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, enum novas_equatorial_class coordType, double xp, double yp, const double *in, double *out) |
short | cel_pole (double jd_tt, enum novas_pole_offset_type type, double dpole1, double dpole2) |
short | cio_array (double jd_tdb, long n_pts, ra_of_cio *restrict cio) |
short | cio_basis (double jd_tdb, double ra_cio, enum novas_cio_location_type loc_type, enum novas_accuracy accuracy, double *restrict x, double *restrict y, double *restrict z) |
short | cio_location (double jd_tdb, enum novas_accuracy accuracy, double *restrict ra_cio, short *restrict loc_type) |
short | cio_ra (double jd_tt, enum novas_accuracy accuracy, double *restrict ra_cio) |
Computes the true right ascension of the celestial intermediate origin (CIO) vs the equinox of date on the true equator of date for a given TT Julian date. | |
double | cirs_to_app_ra (double jd_tt, enum novas_accuracy accuracy, double ra) |
Converts a CIRS right ascension coordinate (measured from the CIO) to an apparent R.A. | |
int | cirs_to_gcrs (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from the Celestial Intermediate Reference System (CIRS) frame at the given epoch to the ICRS / GCRS. | |
int | cirs_to_itrs (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out) |
Rotates a position vector from the dynamical CIRS frame of date to the Earth-fixed ITRS frame (IAU 2000 standard method). | |
int | cirs_to_tod (double jd_tt, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from the Celestial Intermediate Reference System (CIRS) at the given epoch to the True of Date (TOD) reference system. | |
double | d_light (const double *pos_src, const double *pos_body) |
Returns the difference in light-time, for a star, between the barycenter of the solar system and the observer (or the geocenter) (Usage A). | |
int | e_tilt (double jd_tdb, enum novas_accuracy accuracy, double *restrict mobl, double *restrict tobl, double *restrict ee, double *restrict dpsi, double *restrict deps) |
(primarily for internal use) Computes quantities related to the orientation of the Earth's rotation axis at the specified Julian date. | |
int | ecl2equ (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, double elon, double elat, double *restrict ra, double *restrict dec) |
Convert ecliptic longitude and latitude to right ascension and declination. | |
short | ecl2equ_vec (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, const double *in, double *out) |
Converts an ecliptic position vector to an equatorial position vector. | |
short | equ2ecl (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, double ra, double dec, double *restrict elon, double *restrict elat) |
Convert right ascension and declination to ecliptic longitude and latitude. | |
short | equ2ecl_vec (double jd_tt, enum novas_equator_type coord_sys, enum novas_accuracy accuracy, const double *in, double *out) |
Converts an equatorial position vector to an ecliptic position vector. | |
int | equ2gal (double ra, double dec, double *restrict glon, double *restrict glat) |
Converts ICRS right ascension and declination to galactic longitude and latitude. | |
int | equ2hor (double jd_ut1, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const on_surface *restrict location, double ra, double dec, enum novas_refraction_model ref_option, double *restrict zd, double *restrict az, double *restrict rar, double *restrict decr) |
double | era (double jd_ut1_high, double jd_ut1_low) |
Returns the value of the Earth Rotation Angle (θ) for a given UT1 Julian date. | |
int | frame_tie (const double *in, enum novas_frametie_direction direction, double *out) |
Transforms a vector from the dynamical reference system to the International Celestial Reference System (ICRS), or vice versa. | |
int | fund_args (double t, novas_delaunay_args *restrict a) |
Compute the fundamental (a.k.a. | |
int | gal2equ (double glon, double glat, double *restrict ra, double *restrict dec) |
Converts galactic longitude and latitude to ICRS right ascension and declination. | |
short | gcrs2equ (double jd_tt, enum novas_dynamical_type sys, enum novas_accuracy accuracy, double rag, double decg, double *restrict ra, double *restrict dec) |
Converts GCRS right ascension and declination to coordinates with respect to the equator of date (mean or true). | |
int | gcrs_to_cirs (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from the ICRS / GCRS to the Celestial Intermediate Reference System (CIRS) frame at the given epoch. | |
int | gcrs_to_j2000 (const double *in, double *out) |
Changes ICRS / GCRS coordinates to J2000 coordinates. | |
int | gcrs_to_mod (double jd_tdb, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from the ICRS / GCRS to the Mean of Date (MOD) reference frame at the given epoch. | |
int | gcrs_to_tod (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from the ICRS / GCRS to the True of Date (TOD) reference frame at the given epoch. | |
short | geo_posvel (double jd_tt, double ut1_to_tt, enum novas_accuracy accuracy, const observer *restrict obs, double *restrict pos, double *restrict vel) |
Computes the geocentric GCRS position and velocity of an observer. | |
novas_nutation_provider | get_nutation_lp_provider () |
Returns the function configured for low-precision IAU 2000 nutation calculations instead of the default nu2000k(). | |
double | get_ut1_to_tt (int leap_seconds, double dut1) |
Returns the TT - UT1 time difference given the leap seconds and the actual UT1 - UTC time difference as measured and published by IERS. | |
double | get_utc_to_tt (int leap_seconds) |
Returns the difference between Terrestrial Time (TT) and Universal Coordinated Time (UTC) | |
short | grav_def (double jd_tdb, enum novas_observer_place unused, enum novas_accuracy accuracy, const double *pos_src, const double *pos_obs, double *out) |
(primarily for internal use) Computes the total gravitational deflection of light for the observed object due to the major gravitating bodies in the solar system. | |
int | grav_planets (const double *pos_src, const double *pos_obs, const novas_planet_bundle *restrict planets, double *out) |
(primarily for internal use) Computes the total gravitational deflection of light for the observed object due to the specified gravitating bodies in the solar system. | |
double | grav_redshift (double M_kg, double r_m) |
Returns the gravitational redshift (z) for light emitted near a massive spherical body at some distance from its center, and observed at some very large (infinite) distance away. | |
int | grav_undef (double jd_tdb, enum novas_accuracy accuracy, const double *pos_app, const double *pos_obs, double *out) |
(primarily for internal use) Computes the gravitationally undeflected position of an observed source position due to the major gravitating bodies in the solar system. | |
int | grav_undo_planets (const double *pos_app, const double *pos_obs, const novas_planet_bundle *restrict planets, double *out) |
(primarily for internal use) Computes the gravitationally undeflected position of an observed source position due to the specified Solar-system bodies. | |
int | grav_vec (const double *pos_src, const double *pos_obs, const double *pos_body, double rmass, double *out) |
(primarily for internal use) Corrects position vector for the deflection of light in the gravitational field of anarbitrary body. | |
int | hor_to_itrs (const on_surface *restrict location, double az, double za, double *restrict itrs) |
Converts astrometric (unrefracted) azimuth and zenith angles at the specified observer location to a unit position vector in the Earth-fixed ITRS frame. | |
int | iau2000a (double jd_tt_high, double jd_tt_low, double *restrict dpsi, double *restrict deps) |
Computes the IAU 2000A high-precision nutation series for the specified date, to 0.1 μas accuracy. | |
int | iau2000b (double jd_tt_high, double jd_tt_low, double *restrict dpsi, double *restrict deps) |
Computes the forced nutation of the non-rigid Earth based at reduced precision. | |
double | ira_equinox (double jd_tdb, enum novas_equinox_type equinox, enum novas_accuracy accuracy) |
Compute the intermediate right ascension of the equinox at the input Julian date, using an analytical expression for the accumulated precession in right ascension. | |
int | itrs_to_cirs (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out) |
Rotates a position vector from the Earth-fixed ITRS frame to the dynamical CIRS frame of date (IAU 2000 standard method). | |
int | itrs_to_hor (const on_surface *restrict location, const double *restrict itrs, double *restrict az, double *restrict za) |
Converts a position vector in the Earth-fixed ITRS frame to astrometric (unrefracted) azimuth and zenith angles at the specified observer location. | |
int | itrs_to_tod (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out) |
Rotates a position vector from the Earth-fixed ITRS frame to the dynamical True of Date (TOD) frame of date (pre IAU 2000 method). | |
int | j2000_to_gcrs (const double *in, double *out) |
Change J2000 coordinates to ICRS / GCRS coordinates. | |
int | j2000_to_tod (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from J2000 coordinates to the True of Date (TOD) reference frame at the given epoch. | |
double | julian_date (short year, short month, short day, double hour) |
Returns the Julian day for a given astronomical calendar date. | |
short | light_time (double jd_tdb, const object *restrict body, const double *pos_obs, double tlight0, enum novas_accuracy accuracy, double *pos_src_obs, double *restrict tlight) |
Computes the geocentric position of a solar system body, as antedated for light-time. | |
int | light_time2 (double jd_tdb, const object *restrict body, const double *restrict pos_obs, double tlight0, enum novas_accuracy accuracy, double *p_src_obs, double *restrict v_ssb, double *restrict tlight) |
Computes the geocentric position and velocity of a solar system body, as antedated for light-time. | |
int | limb_angle (const double *pos_src, const double *pos_obs, double *restrict limb_ang, double *restrict nadir_ang) |
Determines the angle of an object above or below the Earth's limb (horizon). | |
short | local_planet (double jd_tt, const object *restrict ss_body, double ut1_to_tt, const on_surface *restrict position, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis) |
reference_system: ICRS/GCRS observer location: topocentric (on Earth) position_type: apparent | |
short | local_star (double jd_tt, double ut1_to_tt, const cat_entry *restrict star, const on_surface *restrict position, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec) |
reference_system: ICRS / GCRS observer location: topocentric (on Earth) position_type: apparent | |
int | make_airborne_observer (const on_surface *location, const double *vel, observer *obs) |
Populates an 'observer' data structure for an observer moving relative to the surface of Earth, such as an airborne observer. | |
short | make_cat_entry (const char *restrict name, const char *restrict catalog, long cat_num, double ra, double dec, double pm_ra, double pm_dec, double parallax, double rad_vel, cat_entry *source) |
Fully populates the data structure for a catalog source, such as a star, with all parameters provided at once. | |
int | make_cat_object (const cat_entry *star, object *source) |
Populates and object data structure with the data for a catalog source. | |
int | make_cat_object_sys (const cat_entry *star, const char *restrict system, object *source) |
Populates and object data structure with the data for a catalog source for a given system of catalog coordinates. | |
int | make_ephem_object (const char *name, long num, object *body) |
Sets a celestial object to be a Solar-system ephemeris body. | |
int | make_gps_observer (double latitude, double longitude, double height, observer *obs) |
Initializes an observer data structure for a ground-based observer with the specified GPS / WGS84 location, and sets mean (annual) weather parameters based on that location. | |
int | make_gps_site (double latitude, double longitude, double height, on_surface *site) |
Initializes an observing site with the specified GPS / WGS84 location, and sets mean (annual) weather parameters based on that location. | |
int | make_in_space (const double *sc_pos, const double *sc_vel, in_space *loc) |
Populates an 'in_space' data structure, for an observer situated on a near-Earth spacecraft, with the provided position and velocity components. | |
int | make_itrf_observer (double latitude, double longitude, double height, observer *obs) |
Initializes an observer data structure for a ground-based observer with the specified International Terrestrial Reference Frame (ITRF) / GRS80 location, and sets mean (annual) weather parameters based on that location. | |
int | make_itrf_site (double latitude, double longitude, double height, on_surface *site) |
Initializes an observing site with the specified International Terrestrial Reference Frame (ITRF) / GRS80 location, and sets mean (annual) weather parameters based on that location. | |
short | make_object (enum novas_object_type, long number, const char *name, const cat_entry *star, object *source) |
short | make_observer (enum novas_observer_place, const on_surface *loc_surface, const in_space *loc_space, observer *obs) |
int | make_observer_at_geocenter (observer *restrict obs) |
Populates an 'observer' data structure for a hypothetical observer located at Earth's geocenter. | |
int | make_observer_at_site (const on_surface *restrict site, observer *restrict obs) |
Initializes an observer data structure for a ground-based observer at the specified observing site, and sets mean (annual) weather parameters based on that location. | |
int | make_observer_in_space (const double *sc_pos, const double *sc_vel, observer *obs) |
Populates an 'observer' data structure, for an observer situated on a near-Earth spacecraft, with the specified geocentric position and velocity vectors. | |
int | make_observer_on_surface (double latitude, double longitude, double height, double temperature, double pressure, observer *restrict obs) |
int | make_on_surface (double latitude, double longitude, double height, double temperature, double pressure, on_surface *restrict loc) |
int | make_planet (enum novas_planet num, object *restrict planet) |
Sets a celestial object to be a major planet, or the Sun, Moon, Solar-system Barycenter, etc. | |
int | make_redshifted_cat_entry (const char *name, double ra, double dec, double z, cat_entry *source) |
Populates a celestial object data structure with the parameters for a redhifted catalog source, such as a distant quasar or galaxy. | |
int | make_redshifted_object (const char *name, double ra, double dec, double z, object *source) |
Populates a celestial object data structure with the ICRS astrometric parameters for a redhifted catalog source, such as a distant quasar or galaxy. | |
int | make_redshifted_object_sys (const char *name, double ra, double dec, const char *restrict system, double z, object *source) |
Populates a celestial object data structure with the parameters for a redhifted catalog source, such as a distant quasar or galaxy, for a given system of catalog coordinates. | |
int | make_solar_system_observer (const double *sc_pos, const double *sc_vel, observer *obs) |
Populates an 'observer' data structure, for an observer situated on a near-Earth spacecraft, with the specified geocentric position and velocity vectors. | |
int | make_xyz_site (const double *restrict xyz, on_surface *restrict site) |
Initializes an observing site with the specified Cartesian geocentric xyz location, and sets mean (annual) weather parameters based on that location. | |
double | mean_obliq (double jd_tdb) |
Computes the mean obliquity of the ecliptic. | |
short | mean_star (double jd_tt, double tra, double tdec, enum novas_accuracy accuracy, double *restrict ira, double *restrict idec) |
reference system: TOD → ICRS observer location: geocenter → SSB position_type: apparent → geometric | |
int | mod_to_gcrs (double jd_tdb, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from Mean of Date (MOD) reference frame at the given epoch to the ICRS / GCRS. | |
int | novas_app_to_geom (const novas_frame *restrict frame, enum novas_reference_system sys, double ra, double dec, double dist, double *restrict geom_icrs) |
Converts an observed apparent sky position of a source to an ICRS geometric position, by undoing the gravitational deflection and aberration corrections. | |
int | novas_app_to_hor (const novas_frame *restrict frame, enum novas_reference_system sys, double ra, double dec, RefractionModel ref_model, double *restrict az, double *restrict el) |
Converts an observed apparent position vector in the specified coordinate system to local horizontal coordinates in the specified observer frame. | |
int | novas_cartesian_to_geodetic (const double *restrict xyz, enum novas_reference_ellipsoid ellipsoid, double *restrict lon, double *restrict lat, double *restrict alt) |
Converts geocentric Cartesian site coordinates to geodetic coordinates on the given reference ellipsoid. | |
void | novas_case_sensitive (int value) |
Enables or disables case-sensitive processing of the object name. | |
int | novas_change_observer (const novas_frame *orig, const observer *obs, novas_frame *out) |
Change the observer location for an observing frame. | |
double | novas_clock_skew (const novas_frame *frame, enum novas_timescale timescale) |
Returns the instantaneous incremental rate at which the observer's clock (i.e. | |
double | novas_date (const char *restrict date) |
Returns a Julian date (in non-specific timescale) corresponding the specified input string date / time. | |
double | novas_date_scale (const char *restrict date, enum novas_timescale *restrict scale) |
Returns a Julian date and the timescale corresponding the specified input string date/time and timescale marker. | |
int | novas_day_of_week (double tjd) |
Returns the one-based ISO 8601 day-of-week index of a given Julian Date. | |
int | novas_day_of_year (double tjd, enum novas_calendar_type calendar, int *restrict year) |
Returns the one-based day index in the calendar year for a given Julian Date. | |
void | novas_debug (enum novas_debug_mode mode) |
Enables or disables reporting errors and traces to the standard error stream. | |
double | novas_diff_tcb (const novas_timespec *t1, const novas_timespec *t2) |
Returns the Barycentric Coordinate Time (TCB) based time difference (t1 - t2) in days between two astronomical time specifications. | |
double | novas_diff_tcg (const novas_timespec *t1, const novas_timespec *t2) |
Returns the Geocentric Coordinate Time (TCG) based time difference (t1 - t2) in days between two astronomical time specifications. | |
double | novas_diff_time (const novas_timespec *t1, const novas_timespec *t2) |
Returns the Terrestrial Time (TT) based time difference (t1 - t2) in days between two astronomical time specifications. | |
int | novas_diurnal_eop (double gmst, const novas_delaunay_args *restrict delaunay, double *restrict xp, double *restrict yp, double *restrict dut1) |
Calculate corrections to the Earth orientation parameters (EOP) due to short term (diurnal and semidiurnal) libration and the ocean tides. | |
int | novas_diurnal_eop_at_time (const novas_timespec *restrict time, double *restrict dxp, double *restrict dyp, double *restrict dut1) |
Calculate corrections to the Earth orientation parameters (EOP) due to short term (diurnal and semidiurnal) libration and the ocean tides at a given astromtric time. | |
int | novas_diurnal_libration (double gmst, const novas_delaunay_args *restrict delaunay, double *restrict xp, double *restrict yp, double *restrict dut1) |
Calculate diurnal and semi-diurnal libration corrections to the Earth orientation parameters (EOP) for the non-rigid Earth. | |
int | novas_diurnal_ocean_tides (double gmst, const novas_delaunay_args *restrict delaunay, double *restrict xp, double *restrict yp, double *restrict dut1) |
Calculate corrections to the Earth orientation parameters (EOP) due to the ocean tides. | |
double | novas_dms_degrees (const char *restrict dms) |
Returns the decimal degrees for a DMS string specification. | |
int | novas_e2h_offset (double dra, double ddec, double pa, double *restrict daz, double *restrict del) |
Converts coordinate offsets, from the local equatorial system to local horizontal offsets. | |
double | novas_epa (double ha, double dec, double lat) |
Returns the Parallactic Angle (PA) calculated for an RA/Dec location of the sky at a given sidereal time. | |
double | novas_epoch (const char *restrict system) |
Returns the Julian day corresponding to an astronomical coordinate epoch. | |
double | novas_equ_sep (double ra1, double dec1, double ra2, double dec2) |
Returns the angular separation of two equatorial locations on a sphere. | |
int | novas_equ_track (const object *restrict source, const novas_frame *restrict frame, double dt, novas_track *restrict track) |
Calculates equatorial tracking position and motion (first and second time derivatives) for the specified source in the given observing frame. | |
double | novas_frame_lst (const novas_frame *restrict frame) |
Returns the Local (apparent) Sidereal Time for an observing frame of an Earth-bound observer. | |
double | novas_gast (double jd_ut1, double ut1_to_tt, enum novas_accuracy accuracy) |
Returns the Greenwich Apparent Sidereal Time (GAST) for a given UT1 date. | |
int | novas_geodetic_to_cartesian (double lon, double lat, double alt, enum novas_reference_ellipsoid ellipsoid, double *xyz) |
Converts geodetic site coordinates to geocentric Cartesian coordinates, using the specified reference ellipsoid. | |
int | novas_geodetic_transform_site (enum novas_reference_ellipsoid from_ellipsoid, const on_surface *in, enum novas_reference_ellipsoid to_ellipsoid, on_surface *out) |
Transforms a geodetic location from one reference ellipsoid to another. | |
int | novas_geom_posvel (const object *restrict source, const novas_frame *restrict frame, enum novas_reference_system sys, double *restrict pos, double *restrict vel) |
Calculates the geometric position and velocity vectors, relative to the observer, for a source in the given observing frame, in the specified coordinate system of choice. | |
int | novas_geom_to_app (const novas_frame *restrict frame, const double *restrict pos, enum novas_reference_system sys, sky_pos *restrict out) |
Converts an geometric position in ICRS to an apparent position on sky, by applying appropriate corrections for aberration and gravitational deflection for the observer's frame. | |
enum novas_debug_mode | novas_get_debug_mode () |
Returns the current, thread-local, mode for reporting errors encountered (and traces). | |
double | novas_get_split_time (const novas_timespec *restrict time, enum novas_timescale timescale, long *restrict ijd) |
Returns the fractional Julian date of an astronomical time in the specified timescale, as an integer and fractional part. | |
double | novas_get_time (const novas_timespec *restrict time, enum novas_timescale timescale) |
Returns the fractional Julian date of an astronomical time in the specified timescale. | |
time_t | novas_get_unix_time (const novas_timespec *restrict time, long *restrict nanos) |
Returns the UNIX time for an astronomical time instant. | |
double | novas_gmst (double jd_ut1, double ut1_to_tt) |
Returns the Greenwich Mean Sidereal Time (GMST) for a given UT1 date, using eq. | |
int | novas_h2e_offset (double daz, double del, double pa, double *restrict dra, double *restrict ddec) |
Converts coordinate offsets, from the local horizontal system to local equatorial offsets. | |
double | novas_hms_hours (const char *restrict hms) |
Returns the decimal hours for a HMS string specification. | |
int | novas_hor_to_app (const novas_frame *restrict frame, double az, double el, RefractionModel ref_model, enum novas_reference_system sys, double *restrict ra, double *restrict dec) |
Converts an observed azimuth and elevation coordinate to right ascension (R.A.) and declination coordinates expressed in the coordinate system of choice. | |
int | novas_hor_track (const object *restrict source, const novas_frame *restrict frame, RefractionModel ref_model, novas_track *restrict track) |
Calculates horizontal tracking position and motion (first and second time derivatives) for the specified source in the given observing frame. | |
double | novas_hpa (double az, double el, double lat) |
Returns the Parallactic Angle (PA) calculated for a horizontal Az/El location of the sky. | |
int | novas_init_cat_entry (cat_entry *restrict source, const char *restrict name, double ra, double dec) |
Initializes a catalog entry, such as a star or a distant galaxy, with a name and coordinates. | |
double | novas_inv_refract (RefractionModel model, double jd_tt, const on_surface *restrict loc, enum novas_refraction_type type, double el0) |
Computes the reverse atmospheric refraction for a given refraction model. | |
int | novas_invert_transform (const novas_transform *transform, novas_transform *inverse) |
Inverts a novas coordinate transformation matrix. | |
int | novas_iso_timestamp (const novas_timespec *restrict time, char *restrict dst, int maxlen) |
Prints a UTC-based ISO timestamp to millisecond precision to the specified string buffer. | |
int | novas_itrf_transform (int from_year, const double *restrict from_coords, const double *restrict from_rates, int to_year, double *to_coords, double *to_rates) |
Converts ITRF coordinates between different realizations of the ITRF coordinate system. | |
int | novas_itrf_transform_eop (int from_year, double from_xp, double from_yp, double from_dut1, int to_year, double *restrict to_xp, double *restrict to_yp, double *restrict to_dut1) |
Transforms Earth orientation parameters (xp, yp, dUT1) from one ITRF realization to another. | |
int | novas_itrf_transform_site (int from_year, const on_surface *in, int to_year, on_surface *out) |
Transforms a geodetic location between two International Terrestrial Reference Frame (ITRF) realizations. | |
double | novas_jd_from_date (enum novas_calendar_type calendar, int year, int month, int day, double hour) |
Returns the Julian day for a given calendar date. | |
int | novas_jd_to_date (double tjd, enum novas_calendar_type calendar, int *restrict year, int *restrict month, int *restrict day, double *restrict hour) |
This function will compute a broken down date on the specified calendar for given the Julian day input. | |
int | novas_los_to_xyz (const double *los, double lon, double lat, double *xyz) |
Converts a 3D line-of-sight vector (δφ, δθ δr) to a rectangular equatorial (δx, δy, δz) vector. | |
double | novas_lsr_to_ssb_vel (double epoch, double ra, double dec, double vLSR) |
Returns a Solar System Baricentric (SSB) radial velocity for a radial velocity that is referenced to the Local Standard of Rest (LSR). | |
int | novas_make_frame (enum novas_accuracy accuracy, const observer *obs, const novas_timespec *time, double xp, double yp, novas_frame *frame) |
Sets up a observing frame for a specific observer location, time of observation, and accuracy requirement. | |
int | novas_make_transform (const novas_frame *frame, enum novas_reference_system from_system, enum novas_reference_system to_system, novas_transform *transform) |
Calculates a transformation matrix that can be used to convert positions and velocities from one coordinate reference system to another. | |
double | novas_mean_clock_skew (const novas_frame *frame, enum novas_timescale timescale) |
Returns the averaged incremental rate at which the observer's clock (i.e. | |
double | novas_norm_ang (double angle) |
Returns the normalized angle in the [0:2π) range. | |
double | novas_object_sep (const object *source1, const object *source2, const novas_frame *restrict frame) |
Returns the angular separation of two objects from the observer's point of view. | |
int | novas_offset_time (const novas_timespec *time, double seconds, novas_timespec *out) |
Increments the astrometric time by a given amount. | |
double | novas_optical_refraction (double jd_tt, const on_surface *loc, enum novas_refraction_type type, double el) |
Returns an optical refraction correction using the weather parameters defined for the observer location. | |
double | novas_parse_date (const char *restrict date, char **restrict tail) |
Parses an astronomical date/time string into a Julian date specification. | |
double | novas_parse_date_format (enum novas_calendar_type calendar, enum novas_date_format format, const char *restrict date, char **restrict tail) |
Parses a calendar date/time string, expressed in the specified type of calendar, into a Julian day (JD). | |
double | novas_parse_degrees (const char *restrict str, char **restrict tail) |
Parses an angle in degrees from a string that contains either a decimal degrees or else a broken-down DMS representation. | |
double | novas_parse_dms (const char *restrict str, char **restrict tail) |
Parses the decimal degrees for a DMS string specification. | |
double | novas_parse_hms (const char *restrict str, char **restrict tail) |
Parses the decimal hours for a HMS string specification. | |
double | novas_parse_hours (const char *restrict str, char **restrict tail) |
Parses a time or time-like angle from a string that contains either a decimal hours or else a broken-down HMS representation. | |
double | novas_parse_iso_date (const char *restrict date, char **restrict tail) |
Parses an ISO 8601 timestamp, converting it to a Julian day. | |
enum novas_timescale | novas_parse_timescale (const char *restrict str, char **restrict tail) |
Parses the timescale from a string containing a standard abbreviation (case insensitive), and returns the updated parse position after the timescale specification (if any). | |
int | novas_print_dms (double degrees, enum novas_separator_type sep, int decimals, char *restrict buf, int len) |
Prints an angle in degrees as [-]ddd:mm:ss[.S...] into the specified buffer, with up to nanosecond precision. | |
int | novas_print_hms (double hours, enum novas_separator_type sep, int decimals, char *restrict buf, int len) |
Prints a time in hours as hh:mm:ss[.S...] into the specified buffer, with up to nanosecond precision. | |
int | novas_print_timescale (enum novas_timescale scale, char *restrict buf) |
Prints the standard string representation of the timescale to the specified buffer. | |
double | novas_radio_refraction (double jd_tt, const on_surface *loc, enum novas_refraction_type type, double el) |
Atmospheric refraction model for radio wavelengths (Berman & Rockwell 1976). | |
int | novas_refract_wavelength (double microns) |
Sets the observing wavelength for which refraction is to be calculated when using a wavelength-depenendent model, such as novas_wave_refraction(). | |
double | novas_rises_above (double el, const object *restrict source, const novas_frame *restrict frame, RefractionModel ref_model) |
Returns the UTC date at which a distant source appears to rise above the specified elevation angle. | |
double | novas_sep (double lon1, double lat1, double lon2, double lat2) |
Returns the angular separation of two locations on a sphere. | |
int | novas_set_catalog (cat_entry *restrict source, const char *restrict catalog, long num) |
Sets the optional catalog information for a sidereal source. | |
int | novas_set_current_time (int leap, double dut1, novas_timespec *restrict time) |
Sets the time eith the UNIX time obtained from the system clock. | |
int | novas_set_default_weather (on_surface *site) |
Sets default weather parameters based on an approximate global model to the mean annualized temperatures, based on Feulner et al. | |
int | novas_set_distance (cat_entry *source, double parsecs) |
Sets the distance for a catalog source. | |
int | novas_set_lsr_vel (cat_entry *source, double epoch, double v_kms) |
Sets a radial velocity for a catalog source, defined w.r.t. | |
int | novas_set_parallax (cat_entry *source, double mas) |
Sets the parallax (a measure of distance) for a catalog source. | |
int | novas_set_proper_motion (cat_entry *source, double pm_ra, double pm_dec) |
Sets the proper-motion for a (Galactic) catalog source. | |
int | novas_set_redshift (cat_entry *source, double z) |
Sets a redhift for a catalog source, as a relativistic measure of velocity. | |
int | novas_set_split_time (enum novas_timescale timescale, long ijd, double fjd, int leap, double dut1, novas_timespec *restrict time) |
Sets an astronomical time to the split Julian Date value, defined in the specified timescale. | |
int | novas_set_ssb_vel (cat_entry *source, double v_kms) |
Sets a radial velocity for a catalog source, defined w.r.t. | |
int | novas_set_str_time (enum novas_timescale timescale, const char *restrict str, int leap, double dut1, novas_timespec *restrict time) |
Sets astronomical time in a specific timescale using a string specification. | |
int | novas_set_time (enum novas_timescale timescale, double jd, int leap, double dut1, novas_timespec *restrict time) |
Sets an astronomical time to the fractional Julian Date value, defined in the specified timescale. | |
int | novas_set_unix_time (time_t unix_time, long nanos, int leap, double dut1, novas_timespec *restrict time) |
Sets an astronomical time to a UNIX time value. | |
double | novas_sets_below (double el, const object *restrict source, const novas_frame *restrict frame, RefractionModel ref_model) |
Returns the UTC date at which a distant source appears to set below the specified elevation angle. | |
int | novas_sky_pos (const object *restrict object, const novas_frame *restrict frame, enum novas_reference_system sys, sky_pos *restrict out) |
Calculates an apparent location on sky for the source. | |
double | novas_ssb_to_lsr_vel (double epoch, double ra, double dec, double vLSR) |
Returns a radial-velocity referenced to the Local Standard of Rest (LSR) for a given Solar-System Barycentric (SSB) radial velocity. | |
double | novas_standard_refraction (double jd_tt, const on_surface *loc, enum novas_refraction_type type, double el) |
Returns an optical refraction correction for a standard atmosphere. | |
double | novas_str_degrees (const char *restrict dms) |
Returns an angle parsed from a string that contains either a decimal degrees or else a broken-down DMS representation. | |
double | novas_str_hours (const char *restrict hms) |
Returns a time or time-like angleparsed from a string that contains either a decimal hours or else a broken-down HMS representation. | |
double | novas_time_gst (const novas_timespec *restrict time, enum novas_accuracy accuracy) |
Returns the Greenwich (apparent) Sidereal Time (GST/GaST) for a given astronomical time specification. | |
double | novas_time_lst (const novas_timespec *restrict time, double lon, enum novas_accuracy accuracy) |
Returns the Local (apparent) Sidereal Time (LST/LaST) for a given astronomical time specification and observer location. | |
enum novas_timescale | novas_timescale_for_string (const char *restrict str) |
Returns the timescale constant for a string that denotes the timescale in with a standard abbreviation (case insensitive). | |
int | novas_timestamp (const novas_timespec *restrict time, enum novas_timescale scale, char *restrict dst, int maxlen) |
Prints an astronomical timestamp to millisecond precision in the specified timescale to the specified string buffer. | |
int | novas_track_pos (const novas_track *track, const novas_timespec *time, double *restrict lon, double *restrict lat, double *restrict dist, double *restrict z) |
Calculates a projected position, distance, and redshift for a source, given its near-term trajectory on sky, in the system for which the track was calculated. | |
int | novas_transform_sky_pos (const sky_pos *in, const novas_transform *restrict transform, sky_pos *out) |
Transforms a position or velocity 3-vector from one coordinate reference system to another. | |
int | novas_transform_vector (const double *in, const novas_transform *restrict transform, double *out) |
Transforms a position or velocity 3-vector from one coordinate reference system to another. | |
double | novas_transit_time (const object *restrict source, const novas_frame *restrict frame) |
Returns the UTC date at which a source transits the local meridian. | |
int | novas_uvw_to_xyz (const double *uvw, double ha, double dec, double *xyz) |
Converts equatorial u,v,w projected (absolute or relative) coordinates to rectangular telescope x,y,z coordinates (in ITRS) to for a specified line of sight. | |
double | novas_v2z (double vel) |
Converts a radial recession velocity to a redshift value (z = Δλ / λrest). | |
double | novas_wave_refraction (double jd_tt, const on_surface *loc, enum novas_refraction_type type, double el) |
The wavelength-dependent IAU atmospheric refraction model, based on the SOFA iauRefco() function, in compliance to the 'SOFA Software License' terms of the original source. | |
int | novas_xyz_to_los (const double *xyz, double lon, double lat, double *los) |
Converts a 3D rectangular equatorial (δx, δy, δz) vector to a polar (δφ, δθ δr) vector along a line-of-sight. | |
int | novas_xyz_to_uvw (const double *xyz, double ha, double dec, double *uvw) |
Converts rectangular telescope x,y,z (absolute or relative) coordinates (in ITRS) to equatorial u,v,w projected coordinates for a specified line of sight. | |
double | novas_z2v (double z) |
Converts a redshift value (z = Δλ / λrest) to a radial velocity (i.e. | |
double | novas_z_add (double z1, double z2) |
Compounds two redshift corrections, e.g. | |
double | novas_z_inv (double z) |
Returns the inverse of a redshift value, that is the redshift for a body moving with the same velocity as the original but in the opposite direction. | |
int | nu2000k (double jd_tt_high, double jd_tt_low, double *restrict dpsi, double *restrict deps) |
Computes the forced nutation of the non-rigid Earth: Model NU2000K. | |
int | nutation (double jd_tdb, enum novas_nutation_direction direction, enum novas_accuracy accuracy, const double *in, double *out) |
Nutates equatorial rectangular coordinates from mean equator and equinox of epoch to true equator and equinox of epoch. | |
int | nutation_angles (double t, enum novas_accuracy accuracy, double *restrict dpsi, double *restrict deps) |
Returns the IAU2000 / 2006 values for nutation in longitude and nutation in obliquity for a given TDB Julian date and the desired level of accuracy. | |
int | obs_planets (double jd_tdb, enum novas_accuracy accuracy, const double *restrict pos_obs, int pl_mask, novas_planet_bundle *restrict planets) |
Calculates the positions and velocities for the Solar-system bodies, e.g. | |
int | obs_posvel (double jd_tdb, double ut1_to_tt, enum novas_accuracy accuracy, const observer *restrict obs, const double *restrict geo_pos, const double *restrict geo_vel, double *restrict pos, double *restrict vel) |
Calculates the ICRS position and velocity of the observer relative to the Solar System Barycenter (SSB). | |
short | place (double jd_tt, const object *restrict source, const observer *restrict location, double ut1_to_tt, enum novas_reference_system coord_sys, enum novas_accuracy accuracy, sky_pos *restrict output) |
Computes the apparent direction of a celestial object at a specified time and in a specified coordinate system and for a given observer location. | |
int | place_cirs (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos) |
reference_system: CIRS observer location: geocenter position_type: apparent | |
int | place_gcrs (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos) |
reference_system: ICRS/GCRS observer location: geocenter position_type: apparent | |
int | place_icrs (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos) |
reference_system: ICRS / GCRS observer location: geocenter position_type: geometric | |
int | place_j2000 (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos) |
reference_system: J2000 observer location: geocenter position_type: apparent | |
int | place_mod (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos) |
reference_system: MOD observer location: geocenter position_type: apparent | |
int | place_star (double jd_tt, const cat_entry *restrict star, const observer *restrict obs, double ut1_to_tt, enum novas_reference_system system, enum novas_accuracy accuracy, sky_pos *restrict pos) |
Computes the apparent place of a star at the specified date, given its catalog mean place, proper motion, parallax, and radial velocity. | |
int | place_tod (double jd_tt, const object *restrict source, enum novas_accuracy accuracy, sky_pos *restrict pos) |
reference_system: TOD observer location: geocenter position_type: apparent | |
double | planet_lon (double t, enum novas_planet planet) |
Returns the planetary longitude, for Mercury through Neptune, w.r.t. | |
short | precession (double jd_tdb_in, const double *in, double jd_tdb_out, double *out) |
Precesses equatorial rectangular coordinates from one epoch to another using the IAU2006 (P03) precession model of Capitaine et al. | |
int | proper_motion (double jd_tdb_in, const double *pos, const double *restrict vel, double jd_tdb_out, double *out) |
Applies proper motion, including foreshortening effects, to a star's position. | |
int | rad_vel (const object *restrict source, const double *restrict pos_src, const double *vel_src, const double *vel_obs, double d_obs_geo, double d_obs_sun, double d_src_sun, double *restrict rv) |
Predicts the radial velocity of the observed object as it would be measured by spectroscopic means. | |
double | rad_vel2 (const object *restrict source, const double *pos_emit, const double *vel_src, const double *pos_det, const double *vel_obs, double d_obs_geo, double d_obs_sun, double d_src_sun) |
Predicts the radial velocity of the observed object as it would be measured by spectroscopic means. | |
int | radec2vector (double ra, double dec, double dist, double *restrict pos) |
Converts equatorial spherical coordinates to a vector (equatorial rectangular coordinates). | |
int | radec_planet (double jd_tt, const object *restrict ss_body, const observer *restrict obs, double ut1_to_tt, enum novas_reference_system sys, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis, double *restrict rv) |
Computes the place of a solar system body at the specified time for an observer in the specified coordinate system. | |
int | radec_star (double jd_tt, const cat_entry *restrict star, const observer *restrict obs, double ut1_to_tt, enum novas_reference_system sys, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict rv) |
Computes the place of a star at date 'jd_tt', for an observer in the specified coordinate system, given the star's ICRS catalog place, proper motion, parallax, and radial velocity. | |
double | redshift_vrad (double vrad, double z) |
Applies an incremental redshift correction to a radial velocity. | |
double | refract (const on_surface *restrict location, enum novas_refraction_model model, double zd_obs) |
Computes atmospheric optical refraction for an observed (already refracted!) zenith distance through the atmosphere. | |
double | refract_astro (const on_surface *restrict location, enum novas_refraction_model model, double zd_astro) |
Computes atmospheric optical refraction for a source at an astrometric zenith distance (e.g. | |
int | set_cio_locator_file (const char *restrict filename) |
int | set_nutation_lp_provider (novas_nutation_provider func) |
Set the function to use for low-precision IAU 2000 nutation calculations instead of the default nu2000k(). | |
short | sidereal_time (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_equinox_type gst_type, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, double *restrict gst) |
int | spin (double angle, const double *in, double *out) |
Transforms a vector from one coordinate system to another with same origin and axes rotated about the z-axis. | |
int | starvectors (const cat_entry *restrict star, double *restrict pos, double *restrict motion) |
Converts angular quantities for stars to vectors. | |
int | tdb2tt (double jd_tdb, double *restrict jd_tt, double *restrict secdiff) |
Computes the Terrestrial Time (TT) based Julian date corresponding to a Barycentric Dynamical Time (TDB) Julian date, and retuns th TDB-TT time difference also. | |
short | ter2cel (double jd_ut1_high, double jd_ut1_low, double ut1_to_tt, enum novas_earth_rotation_measure erot, enum novas_accuracy accuracy, enum novas_equatorial_class coordType, double xp, double yp, const double *in, double *out) |
int | terra (const on_surface *restrict location, double gast, double *restrict pos, double *restrict vel) |
Computes the position and velocity vectors of a terrestrial observer with respect to the center of the Earth, based on the GRS80 reference ellipsoid, used for the International Terrestrial Reference Frame (ITRF) and its realizations. | |
int | tod_to_cirs (double jd_tt, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from the True of Date (TOD) reference system to the Celestial Intermediate Reference System (CIRS) at the given epoch to the . | |
int | tod_to_gcrs (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from True of Date (TOD) reference frame at the given epoch to the ICRS / GCRS. | |
int | tod_to_itrs (double jd_tt_high, double jd_tt_low, double ut1_to_tt, enum novas_accuracy accuracy, double xp, double yp, const double *in, double *out) |
Rotates a position vector from the dynamical True of Date (TOD) frame of date the Earth-fixed ITRS frame (pre IAU 2000 method). | |
int | tod_to_j2000 (double jd_tdb, enum novas_accuracy accuracy, const double *in, double *out) |
Transforms a rectangular equatorial (x, y, z) vector from True of Date (TOD) reference frame at the given epoch to the J2000 coordinates. | |
short | topo_planet (double jd_tt, const object *restrict ss_body, double ut1_to_tt, const on_surface *restrict position, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis) |
reference_system: TOD observer location: topocentric (on Earth) position_type: apparent | |
short | topo_star (double jd_tt, double ut1_to_tt, const cat_entry *restrict star, const on_surface *restrict position, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec) |
reference_system: TOD observer location: topocentric (on Earth) position_type: apparent | |
short | transform_cat (enum novas_transform_type, double jd_tt_in, const cat_entry *in, double jd_tt_out, const char *out_id, cat_entry *out) |
Transform a star's catalog quantities for a change the coordinate system and/or the date for which the positions are calculated. | |
int | transform_hip (const cat_entry *hipparcos, cat_entry *hip_2000) |
Convert Hipparcos catalog data at epoch J1991.25 to epoch J2000.0, for use within NOVAS. | |
double | tt2tdb (double jd_tt) |
Returns the TDB - TT time difference in seconds for a given TT or TDB date, with a maximum error of 10 μs for dates between 1600 and 2200. | |
double | tt2tdb_fp (double jd_tt, double limit) |
Returns the TDB-TT time difference with flexible precision. | |
double | tt2tdb_hp (double jd_tt) |
Returns the TDB-TT time difference with high precision. | |
double | unredshift_vrad (double vrad, double z) |
Undoes an incremental redshift correction that was applied to radial velocity. | |
short | vector2radec (const double *restrict pos, double *restrict ra, double *restrict dec) |
Converts an vector in equatorial rectangular coordinates to equatorial spherical coordinates. | |
short | virtual_planet (double jd_tt, const object *restrict ss_body, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec, double *restrict dis) |
reference_system: ICRS / GCRS observer location: geocenter position_type: apparent | |
short | virtual_star (double jd_tt, const cat_entry *restrict star, enum novas_accuracy accuracy, double *restrict ra, double *restrict dec) |
reference_system: ICRS / GCRS observer location: geocenter position_type: apparent | |
int | wobble (double jd_tt, enum novas_wobble_direction direction, double xp, double yp, const double *in, double *out) |
Corrects a vector in the ITRS (rotating Earth-fixed system) for polar motion, and also corrects the longitude origin (by a tiny amount) to the Terrestrial Intermediate Origin (TIO). | |
SuperNOVAS types, definitions, and function prototypes.
SuperNOVAS astrometry software based on the Naval Observatory Vector Astrometry Software (NOVAS). It has been modified to fix outstanding issues, make it easier to use, and provide a ton of new features.
Based on the NOVAS C Edition, Version 3.1:
U. S. Naval Observatory
Astronomical Applications Dept.
Washington, DC
http://www.usno.navy.mil/USNO/astronomical-applications
#define BARYC NOVAS_BARYCENTER |
Deprecated.
#define HELIOC NOVAS_HELIOCENTER |
Deprecated.
#define NOVAS_DELAUNAY_ARGS_INIT |
#define NOVAS_EQUATOR_TYPES (NOVAS_GCRS_EQUATOR + 1) |
The number of equator types defined in enum novas_equator_type
.
#define NOVAS_MATRIX_IDENTITY |
#define NOVAS_MATRIX_INIT |
#define NOVAS_OBJECT_TYPES (NOVAS_ORBITAL_OBJECT + 1) |
The number of object types distinguished by NOVAS.
#define NOVAS_OBSERVER_PLACES (NOVAS_SOLAR_SYSTEM_OBSERVER + 1) |
The number of observer place types supported.
#define NOVAS_ORIGIN_TYPES (NOVAS_HELIOCENTER + 1) |
The number of different ICSR origins available in NOVAS.
#define NOVAS_PLANET_BUNDLE_INIT |
#define NOVAS_PLANET_INIT | ( | num, | |
name ) |
object
initializer macro for major planets, the Sun, Moon, and barycenters.
num | An enum novas_planet number |
name | The designated planet name |
#define NOVAS_PLANETS (NOVAS_PLUTO_BARYCENTER + 1) |
The number of major planets defined in NOVAS.
#define NOVAS_REFERENCE_SYSTEMS (NOVAS_ITRS + 1) |
The number of basic coordinate reference systems in NOVAS.
#define NOVAS_REFRACTION_MODELS (NOVAS_WAVE_REFRACTION + 1) |
The number of built-in refraction models available in SuperNOVAS.
#define NOVAS_SOLAR_RADIUS 696340000.0 |
[m] Solar radius (photosphere)
#define NOVAS_TIMESCALES (NOVAS_UT1 + 1) |
#define NOVAS_TRANSFORM_TYPES (ICRS_TO_J2000 + 1) |
The number of coordinate transfor types in NOVAS.
#define NOVAS_WOBBLE_DIRECTIONS (WOBBLE_PEF_TO_ITRS + 1) |
System in which CIO location is calculated.
Enumerator | |
---|---|
CIO_VS_GCRS | The location of the CIO relative to the GCRS frame. |
CIO_VS_EQUINOX | The location of the CIO relative to the true equinox in the dynamical frame. |
Constants that determine the type of rotation measure to use.
The class of celestial coordinates used as parameters for ter2cel() and cel2ter().
Enumerator | |
---|---|
NOVAS_REFERENCE_CLASS | Celestial coordinates are in GCRS. |
NOVAS_DYNAMICAL_CLASS | Celestial coordinates are apparent values (CIRS or TOD) |
The convention in which the celestial pole offsets are defined for polar wobble.
Enumerator | |
---|---|
POLE_OFFSETS_DPSI_DEPS | Offsets are Δdψ, Δdε pairs (pre IAU 2006 precession-nutation model). |
POLE_OFFSETS_X_Y | Offsets are dx, dy pairs (IAU 2006 precession-nutation model) |
int aberration | ( | const double * | pos, |
const double * | vobs, | ||
double | lighttime, | ||
double * | out ) |
Corrects position vector for aberration of light.
Algorithm includes relativistic terms.
NOTES:
place()
to account for aberration when calculating the position of the source. REFERENCES:
pos | [AU] Position vector of source relative to observer | |
vobs | [AU/day] Velocity vector of observer, relative to the solar system barycenter. | |
lighttime | [day] Light time from object to Earth (if known). Or set to 0, and this function will compute it as needed. | |
[out] | out | [AU] Position vector, referred to origin at center of mass of the Earth, corrected for aberration. It can be the same vector as one of the inputs. |
References novas_vlen().
int bary2obs | ( | const double * | pos, |
const double * | pos_obs, | ||
double * | out, | ||
double *restrict | lighttime ) |
Moves the origin of coordinates from the barycenter of the solar system to the observer (or the geocenter); i.e., this function accounts for parallax (annual+geocentric or just annual).
REFERENCES:
pos | [AU] Position vector, referred to origin at solar system barycenter. | |
pos_obs | [AU] Position vector of observer (or the geocenter), with respect to origin at solar system barycenter. | |
[out] | out | [AU] Position vector, referred to origin at center of mass of the Earth. It may be NULL if not required, or be the same vector as either of the inputs. |
[out] | lighttime | [day] Light time from object to Earth. It may be NULL if not required. |
References novas_vlen().
short cel2ter | ( | double | jd_ut1_high, |
double | jd_ut1_low, | ||
double | ut1_to_tt, | ||
enum novas_earth_rotation_measure | erot, | ||
enum novas_accuracy | accuracy, | ||
enum novas_equatorial_class | coordType, | ||
double | xp, | ||
double | yp, | ||
const double * | in, | ||
double * | out ) |
Rotates a vector from the celestial to the terrestrial system. Specifically, it transforms a vector in the GCRS, or the dynamical (CIRS or TOD) frames to the ITRS (a rotating Earth-fixed system) by applying rotations for the GCRS-to-dynamical frame tie, precession, nutation, Earth rotation, and polar motion.
If 'system' is NOVAS_CIRS then method EROT_ERA must be used. Similarly, if 'system' is NOVAS_TOD then method must be EROT_ERA. Otherwise an error 3 is returned.
If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. | |
jd_ut1_low | [day] Low-order part of UT1 Julian date. | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
erot | Unused. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
coordType | Input coordinate class, NOVAS_REFERENCE_CLASS (0) or NOVAS_DYNAMICAL_CLASS (1). Use the former if the input coordinates are in the GCRS, and the latter if they are CIRS or TOD (the 'erot' parameter selects which dynamical system the input is specified in.) | |
xp | [arcsec] Conventionally-defined X coordinate of celestial intermediate pole with respect to ITRS pole. If you have defined pole offsets the old (pre IAU2000) way, via cel_pole() , then use 0 here. | |
yp | [arcsec] Conventionally-defined Y coordinate of celestial intermediate pole with respect to ITRS pole. If you have defined pole offsets the old (pre IAU2000) way, via cel_pole() , then use 0 here. | |
in | Input position vector, geocentric equatorial rectangular coordinates in the specified input coordinate system (GCRS if 'class' is NOVAS_REFERENCE_CLASS; or else either CIRS if 'erot' is EROT_ERA, or TOD if 'erot' is EROT_GST). | |
[out] | out | ITRS position vector, geocentric equatorial rectangular coordinates (terrestrial system). It can be the same vector as the input. |
References era(), EROT_ERA, EROT_GST, gcrs_to_cirs(), gcrs_to_tod(), NOVAS_DYNAMICAL_CLASS, NOVAS_FULL_ACCURACY, novas_gast(), NOVAS_REDUCED_ACCURACY, spin(), wobble(), WOBBLE_PEF_TO_ITRS, and WOBBLE_TIRS_TO_ITRS.
short cel_pole | ( | double | jd_tt, |
enum novas_pole_offset_type | type, | ||
double | dpole1, | ||
double | dpole2 ) |
novas_app_to_hor()
/ novas_hor_to_app()
or else wobble()
.Specifies the unmodeled celestial pole offsets for high-precision applications to be applied to the True of Date (TOD) equator, in the old, pre IAU 2006 methodology. Nonetheless, these offsets should be specified relative to the IAU2006 precession / nutation model to provide a correction to the modeled (precessed and nutated) position of Earth's pole, such those derived from observations and published by IERS.
The call sets the global variables PSI_COR
and EPS_COR
, for subsequent calls to e_tilt()
. As such, it should be called to specify pole offsets prior to legacy NOVAS C equinox-specific calls. The global values of PSI_COR
and EPS_COR
specified via this function will be effective until explicitly changed again.
NOTES:
novas_app_to_hor(),
novas_hor_to_app(), or wobble()
). novas_make_frame()
, the offsets will be applied for the TIRS / ITRS conversion only, and not to the TOD equator per se. REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. Used only if 'type' is POLE_OFFSETS_X_Y (2), to transform dx and dy to the equivalent Δδψ and Δδε values. |
type | POLE_OFFSETS_DPSI_DEPS (1) if the offsets are Δδψ, Δδε relative to the IAU 20006 precession/nutation model; or POLE_OFFSETS_X_Y (2) if they are dx, dy offsets relative to the IAU 2000 / 2006 precession-nutation model. |
dpole1 | [mas] Value of celestial pole offset in first coordinate, (Δδψ for or dx) in milliarcseconds, relative to the IAU2006 precession/nutation model. |
dpole2 | [mas] Value of celestial pole offset in second coordinate, (Δδε or dy) in milliarcseconds, relative to the IAU2006 precession/nutation model. |
References POLE_OFFSETS_DPSI_DEPS, and POLE_OFFSETS_X_Y.
short cio_array | ( | double | jd_tdb, |
long | n_pts, | ||
ra_of_cio *restrict | cio ) |
Given an input TDB Julian date and the number of data points desired, this function returns a set of Julian dates and corresponding values of the GCRS right ascension of the celestial intermediate origin (CIO). The range of dates is centered on the requested date.
NOTES:
The CIO locator file that is bundled was prepared with the original IAU2000A nutation model, not with the newer R06 (a.k.a. IAU2006) nutation model, resulting in an error up to the few tens of micro-arcseconds level for dates between 1900 and 2100, and larger errors further away from the current epoch.
CIO_RA.TXT
or cio_ra.bin
). The current version no longer does, and instead generates the requested data on the fly. jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
n_pts | Number of Julian dates and right ascension values requested (not less than 2 or more than NOVAS_CIO_CACHE_SIZE). | |
[out] | cio | A time series (array) of the right ascension of the Celestial Intermediate Origin (CIO) with respect to the GCRS. |
References ra_of_cio::jd_tdb, NOVAS_ARCSEC, NOVAS_HOURANGLE, and ra_of_cio::ra_cio.
short cio_basis | ( | double | jd_tdb, |
double | ra_cio, | ||
enum novas_cio_location_type | loc_type, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | x, | ||
double *restrict | y, | ||
double *restrict | z ) |
Computes the CIRS basis vectors, with respect to the GCRS (geocentric ICRS), of the celestial intermediate system defined by the celestial intermediate pole (CIP) (in the z direction) and the celestial intermediate origin (CIO) (in the x direction).
This function effectively constructs the CIRS to GCRS transformation matrix C in eq. (3) of the reference.
NOTES:
REFERENCES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
ra_cio | [h] Right ascension of the CIO at epoch (hours). | |
loc_type | CIO_VS_GCRS (1) if the cio location is relative to the GCRS or else CIO_VS_EQUINOX (2) if relative to the true equinox of date. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | x | Unit 3-vector toward the CIO, equatorial rectangular coordinates, referred to the GCRS. |
[out] | y | Unit 3-vector toward the y-direction, equatorial rectangular coordinates, referred to the GCRS. |
[out] | z | Unit 3-vector toward north celestial pole (CIP), equatorial rectangular coordinates, referred to the GCRS. |
References CIO_VS_EQUINOX, CIO_VS_GCRS, NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, novas_vlen(), and tod_to_gcrs().
short cio_location | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
double *restrict | ra_cio, | ||
short *restrict | loc_type ) |
cio_ra()
to get the CIO location w.r.t. the equinox of date (on the same dynamical equator), or equivalently ira_equinox()
to return the negated value of the same.Returns the location of the celestial intermediate origin (CIO) for a given Julian date, as a right ascension with respect to the true equinox of date. The CIO is always located on the true equator (= intermediate equator) of date.
NOTES:
jd_tdb | [day] Barycentric Dynamic Time (TDB) based Julian date | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra_cio | [h] Right ascension of the CIO, in hours, or NAN if returning with an error. |
[out] | loc_type | Pointer in which to return the reference system in which right ascension is given, which is either CIO_VS_GCRS (1) if the location was obtained via interpolation of the available data file, or else CIO_VS_EQUINOX (2) if it was calculated locally. It is set to -1 if returning with an error. |
References CIO_VS_EQUINOX, ira_equinox(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, and NOVAS_TRUE_EQUINOX.
double d_light | ( | const double * | pos_src, |
const double * | pos_body ) |
Returns the difference in light-time, for a star, between the barycenter of the solar system and the observer (or the geocenter) (Usage A).
Alternatively (Usage B), this function returns the light-time from the observer (or the geocenter) to a point on a light ray that is closest to a specific solar system body. For this purpose, 'pos_src' is the position vector toward observed object, with respect to origin at observer (or the geocenter); 'pos_body' is the position vector of solar system body, with respect to origin at observer (or the geocenter), components in AU; and the returned value is the light time to point on line defined by 'pos' that is closest to solar system body (positive if light passes body before hitting observer, i.e., if 'pos_body' is within 90 degrees of 'pos_src').
NOTES:
novas_geom_posvel()
, novas_sky_pos(), or place()
pos_src | Position vector towards observed object, with respect to the SSB (Usage A), or relative to the observer / geocenter (Usage B). |
pos_body | [AU] Position of observer relative to SSB (Usage A), or position of intermediate solar-system body with respect to the observer / geocenter (Usage B). |
References novas_vlen().
int equ2hor | ( | double | jd_ut1, |
double | ut1_to_tt, | ||
enum novas_accuracy | accuracy, | ||
double | xp, | ||
double | yp, | ||
const on_surface *restrict | location, | ||
double | ra, | ||
double | dec, | ||
enum novas_refraction_model | ref_option, | ||
double *restrict | zd, | ||
double *restrict | az, | ||
double *restrict | rar, | ||
double *restrict | decr ) |
novas_app_to_hor()
instead, or else the more explicit (less ambiguous) tod_to_itrs()
followed by itrs_to_hor()
, and possibly following it with an atmospheric refraction correction if appropriate.Transforms topocentric (TOD) apparent right ascension and declination to zenith distance and azimuth. This method should not be used to convert CIRS apparent coordinates (IAU 2000 standard) – for those you should use cirs_to_itrs() followed by itrs_to_hor() instead.
It uses a method that properly accounts for polar motion, which is significant at the sub-arcsecond level. This function can also adjust coordinates for atmospheric refraction.
NOTES:
REFERENCES:
jd_ut1 | [day] UT1 based Julian date | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
xp | [arcsec] Conventionally-defined x coordinate of celestial intermediate pole with respect to ITRS reference pole, e.g. from IERS Bulletin A. If you have defined pole offsets to be incorporated into the TOD input coordinates (pre-IAU2000 method) via cel_pole() , then you should set this to 0. | |
yp | [arcsec] Conventionally-defined y coordinate of celestial intermediate pole with respect to ITRS reference pole, e.g. from IERS Bulletin A. If you have defined pole offsets to be incorporated into the TOD input coordinates (pre-IAU2000 method) via cel_pole() , then you should set this to 0. | |
location | Geodetic (ITRF / GRS80) observer location on Earth | |
ra | [h] Topocentric apparent (TOD) right ascension of object of interest, referred to true equator and equinox of date. | |
dec | [deg] Topocentric apparent (TOD) declination of object of interest, referred to true equator and equinox of date. | |
ref_option | Refraction model to use. E.g., NOVAS_STANDARD_ATMOSPHERE (1), or NOVAS_WEATHER_AT_LOCATION (2) if to use the weather. | |
[out] | zd | [deg] Topocentric zenith distance in degrees (unrefracted). |
[out] | az | [deg] Topocentric azimuth (measured east from north) in degrees. |
[out] | rar | [h] Topocentric right ascension of object of interest, in hours, referred to true equator and equinox of date, affected by refraction if 'ref_option' is non-zero. (It may be NULL if not required) |
[out] | decr | [deg] Topocentric declination of object of interest, in degrees, referred to true equator and equinox of date. (It may be NULL if not required) |
References EROT_GST, NOVAS_DYNAMICAL_CLASS, refract_astro(), and ter2cel().
short grav_def | ( | double | jd_tdb, |
enum novas_observer_place | unused, | ||
enum novas_accuracy | accuracy, | ||
const double * | pos_src, | ||
const double * | pos_obs, | ||
double * | out ) |
(primarily for internal use) Computes the total gravitational deflection of light for the observed object due to the major gravitating bodies in the solar system.
This function valid for an observed body within the solar system as well as for a star.
If 'accuracy' is NOVAS_FULL_ACCURACY (0), the deflections due to the Sun, Jupiter, Saturn, and Earth are calculated. Otherwise, only the deflection due to the Sun is calculated. In either case, deflection for a given body is ignored if the observer is within ~1500 km of the center of the gravitating body.
The number of bodies used at full and reduced accuracy can be set by changing the grav_bodies_reduced_accuracy
, and grav_bodies_full_accuracy
lobal variables.
NOTES:
novas_sky_pos()
novas_app_to_geom(),
novas_geom_to_app() and
place()` to calculate gravitational deflections as appropriate for positioning sources precisely. The gravitational deflection due to planets requires a planet calculator function to be configured, e.g. via set_planet_provider_hp(). REFERENCES:
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
unused | The type of observer frame (no longer used) | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1). In full accuracy mode, it will calculate the deflection due to the Sun, Jupiter, Saturn and Earth. In reduced accuracy mode, only the deflection due to the Sun is calculated. | |
pos_src | [AU] Position 3-vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU. | |
pos_obs | [AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU. | |
[out] | out | [AU] Position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, corrected for gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs. |
References grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, grav_planets(), NOVAS_FULL_ACCURACY, and obs_planets().
int grav_planets | ( | const double * | pos_src, |
const double * | pos_obs, | ||
const novas_planet_bundle *restrict | planets, | ||
double * | out ) |
(primarily for internal use) Computes the total gravitational deflection of light for the observed object due to the specified gravitating bodies in the solar system.
This function is valid for an observed body within the solar system as well as for a star.
NOTES:
REFERENCES:
pos_src | [AU] Position 3-vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU. | |
pos_obs | [AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU. | |
planets | Apparent planet data containing positions and velocities for the major gravitating bodies in the solar-system. | |
[out] | out | [AU] Position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, corrected for gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs. |
References d_light(), grav_vec(), NOVAS_PLANETS, NOVAS_RMASS_INIT, and novas_vlen().
int grav_undef | ( | double | jd_tdb, |
enum novas_accuracy | accuracy, | ||
const double * | pos_app, | ||
const double * | pos_obs, | ||
double * | out ) |
(primarily for internal use) Computes the gravitationally undeflected position of an observed source position due to the major gravitating bodies in the solar system.
This function is valid for an observed body within the solar system as well as for a star.
If 'accuracy' is set to zero (full accuracy), three bodies (Sun, Jupiter, and Saturn) are used in the calculation. If the reduced-accuracy option is set, only the Sun is used in the calculation. In both cases, if the observer is not at the geocenter, the deflection due to the Earth is included.
he number of bodies used at full and reduced accuracy can be set by changing the grav_bodies_reduced_accuracy
, and grav_bodies_full_accuracy
lobal variables.
REFERENCES:
jd_tdb | [day] Barycentric Dynamical Time (TDB) based Julian date | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
pos_app | [AU] Apparent position 3-vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU. | |
pos_obs | [AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU. | |
[out] | out | [AU] Nominal position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, without gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs. |
References grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, grav_undo_planets(), NOVAS_FULL_ACCURACY, and obs_planets().
int grav_undo_planets | ( | const double * | pos_app, |
const double * | pos_obs, | ||
const novas_planet_bundle *restrict | planets, | ||
double * | out ) |
(primarily for internal use) Computes the gravitationally undeflected position of an observed source position due to the specified Solar-system bodies.
REFERENCES:
pos_app | [AU] Apparent position 3-vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, components in AU. | |
pos_obs | [AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes, components in AU. | |
planets | Apparent planet data containing positions and velocities for the major gravitating bodies in the solar-system. | |
[out] | out | [AU] Nominal position vector of observed object, with respect to origin at observer (or the geocenter), referred to ICRS axes, without gravitational deflection, components in AU. It can be the same vector as the input, but not the same as pos_obs. |
References grav_planets(), novas_inv_max_iter, and novas_vlen().
int grav_vec | ( | const double * | pos_src, |
const double * | pos_obs, | ||
const double * | pos_body, | ||
double | rmass, | ||
double * | out ) |
(primarily for internal use) Corrects position vector for the deflection of light in the gravitational field of anarbitrary body.
This function valid for an observed body within the solar system as well as for a star.
NOTES:
REFERENCES:
pos_src | [AU] Position 3-vector of observed object, with respect to origin at observer (or the geocenter), components in AU. | |
pos_obs | [AU] Position vector of gravitating body, with respect to origin at solar system barycenter, components in AU. | |
pos_body | [AU] Position 3-vector of gravitating body, with respect to origin at solar system barycenter, components in AU. | |
rmass | [1/Msun] Reciprocal mass of gravitating body in solar mass units, that is, Sun mass / body mass. | |
[out] | out | [AU] Position 3-vector of observed object, with respect to origin at observer (or the geocenter), corrected for gravitational deflection, components in AU. It can the same vector as the input. |
References novas_vlen().
short local_planet | ( | double | jd_tt, |
const object *restrict | ss_body, | ||
double | ut1_to_tt, | ||
const on_surface *restrict | position, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ra, | ||
double *restrict | dec, | ||
double *restrict | dis ) |
reference_system: ICRS/GCRS
observer location: topocentric (on Earth)
position_type: apparent
Computes the local apparent place of a solar system body, in the GCRS. This is the same as calling place() for the body for the same observer location and NOVAS_GCRS as the reference system, except the different set of return values used.
The calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Local apparent right ascension, referred to the GCRS (it may be NULL if not required). |
[out] | dec | [deg] Local apparent right ascension, referred to the GCRS (it may be NULL if not required). |
[out] | dis | [AU] Apparent distance from Earth to the body at 'jd_tt' (it may be NULL if not required). |
References make_observer_at_site(), NOVAS_GCRS, and radec_planet().
short local_star | ( | double | jd_tt, |
double | ut1_to_tt, | ||
const cat_entry *restrict | star, | ||
const on_surface *restrict | position, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ra, | ||
double *restrict | dec ) |
reference_system: ICRS / GCRS
observer location: topocentric (on Earth)
position_type: apparent
Computes the local apparent place of a star at date 'jd_tt', in the GCRS, given its catalog mean place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with the same observer location NOVAS_GCRS for an object that specifies the star.
The calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Local apparent right ascension, referred to the GCRS (it may be NULL if not required). |
[out] | dec | [deg] Local apparent right ascension, referred to the GCRS (it may be NULL if not required). |
References make_observer_at_site(), NOVAS_GCRS, and radec_star().
short make_object | ( | enum novas_object_type | type, |
long | number, | ||
const char * | name, | ||
const cat_entry * | star, | ||
object * | source ) |
make_cat_object()
, make_redshifted_object()
, make_planet()
, make_ephem_object()
, or make_orbital_object()
.Populates an object data structure using the parameters provided. By default (for compatibility with NOVAS C) source names are converted to upper-case internally. You can however enable case-sensitive processing by calling novas_case_sensitive() before.
NOTES:
orbit
field (added in v1.2) with zeroes to remain ABI compatible with versions <1.2. type | The type of object. NOVAS_PLANET (0), NOVAS_EPHEM_OBJECT (1) or NOVAS_CATALOG_OBJECT (2), or NOVAS_ORBITAL_OBJECT (3). | |
number | The novas ID number (for solar-system bodies only, otherwise ignored) | |
name | The name of the object (case insensitive). It should be shorter than SIZE_OF_OBJ_NAME or else an error will be returned. The name is converted to upper internally unless novas_case_sensitive() was called before to change that. | |
star | Pointer to structure to populate with the catalog data for a celestial object located outside the solar system. Used only if type is NOVAS_CATALOG_OBJECT, otherwise ignored and can be NULL. | |
[out] | source | Pointer to the celestial object data structure to be populated. |
References object::name, NOVAS_CATALOG_OBJECT, NOVAS_OBJECT_TYPES, NOVAS_PLANET, NOVAS_PLANETS, object::number, object::star, and object::type.
short make_observer | ( | enum novas_observer_place | where, |
const on_surface * | loc_surface, | ||
const in_space * | loc_space, | ||
observer * | obs ) |
make_itrf_observer()
, make_gps_observer()
, make_observer_at_site()
, make_airborne_observer()
make_solar_system_observer(), or
make_observer_at_geocenter()`. This function will be available for the foreseeable future also.Populates an 'observer' data structure given the parameters. The output data structure may be used an the the inputs to NOVAS-C functions, such as make_frame()
or place()
.
where | The location type of the observer | |
loc_surface | Pointer to data structure that defines a location on Earth's surface. Used only if 'where' is NOVAS_OBSERVER_ON_EARTH, otherwise can be NULL. | |
loc_space | Pointer to data structure that defines a near-Earth location in space. Used only if 'where' is NOVAS_OBSERVER_IN_EARTH_ORBIT, otherwise can be NULL. | |
[out] | obs | Pointer to observer data structure to populate. |
References observer::near_earth, NOVAS_AIRBORNE_OBSERVER, NOVAS_OBSERVER_AT_GEOCENTER, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_SOLAR_SYSTEM_OBSERVER, observer::on_surf, in_space::sc_vel, and observer::where.
int make_on_surface | ( | double | latitude, |
double | longitude, | ||
double | height, | ||
double | temperature, | ||
double | pressure, | ||
on_surface *restrict | loc ) |
make_itrf_site()
or make_gps_site()
instead (both of which set default mean annual weather parameters for approximate refraction correction), and optionally set actual weather data afterwards, based on the measurements available. This function will be available for the foreseeable future also.Populates an 'on_surface' data structure, for an observer on the surface of the Earth, with the given parameters.
Note, that because this is an original NOVAS C routine, it does not have an argument to set a humidity value (e.g. for radio refraction). As such, the humidity is set to a a default mean annual value for the location. To set an actual humidity, set the output structure's field after calling this funcion.
NOTES
humidity
) that was not yet part of the on_surface
structure in v1.0. As such, linking SuperNOVAS v1.1 or later with application code compiled for SuperNOVAS v1.0 can result in memory corruption or segmentation fault when this function is called. To be safe, make sure your application has been (re)compiled against SuperNOVAS v1.1 or later. novas_itrf_transform()
, possibly after novas_geodetic_to_cartesian()
with NOVAS_GRS80_ELLIPSOID
as necessary for polar ITRF coordinates. novas_cartesian_to_geodetic()
with NOVAS_GRS80_ELLIPSOID
as the reference ellipsoid parameter. make_gps_site()
instead, and then set weather parameters afterwards as necessary. latitude | [deg] Geodetic (ITRF / GRS80) latitude; north positive. | |
longitude | [deg] Geodetic (ITRF / GRS80) longitude; east positive. | |
height | [m] Geodetic (ITRF / GSR80) altitude above sea level of the observer. | |
temperature | [C] Temperature (degrees Celsius) [-120:70]. | |
pressure | [mbar] Atmospheric pressure [0:1200]. | |
[out] | loc | Pointer to Earth location data structure to populate. |
References make_itrf_site().
short mean_star | ( | double | jd_tt, |
double | tra, | ||
double | tdec, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ira, | ||
double *restrict | idec ) |
reference system: TOD → ICRS
observer location: geocenter → SSB
position_type: apparent → geometric
Computes the barycentric ICRS position of a star, given its True of Date (TOD) apparent place, as seen by an observer at the geocenter, at date 'jd_tt'. Proper motion, parallax and radial velocity are assumed to be zero.
Equivalently, in the new frame-based approach, you might use novas_app_to_geom()
with NOVAS_TOD as the reference system, and distance set to 0, and then convert the geometric position to the output ira
, idec
coordinates using vector2radec()
. The advantage of using novas_app_to_geom()
is that you are not restricted to using TOD input coordinates, but rather the apparent RA/Dec may be specified in any reference system.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
[in] | tra | [h] Geocentric apparent (TOD) right ascension, referred to true equator and equinox of date. |
[in] | tdec | [deg] Geocentric apparent (TOD) declination, referred to true equator and equinox of date. |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ira | [h] Barycentric geometric ICRS right ascension, or NAN if returning with an error code. It may be a pointer to the input value. |
[out] | idec | [deg] Barycentric geometric ICRS declination, or NAN if returning with an error code. It may be a pointer to the input value. |
References app_star(), CAT_ENTRY_INIT, cat_entry::dec, novas_inv_max_iter, precession(), cat_entry::ra, starvectors(), and vector2radec().
short place | ( | double | jd_tt, |
const object *restrict | source, | ||
const observer *restrict | location, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | coord_sys, | ||
enum novas_accuracy | accuracy, | ||
sky_pos *restrict | output ) |
Computes the apparent direction of a celestial object at a specified time and in a specified coordinate system and for a given observer location.
While coord_sys
defines the coordinate referfence system, it is location->where
sets the origin of the reference place relative to which positions and velocities are reported.
For all but ICRS coordinate outputs, the calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
In case of a dynamical equatorial system (such as CIRS or TOD) and an Earth-based observer, the polar wobble parameters set via a prior call to cel_pole() together with he ut1_to_tt argument decide whether the resulting 'topocentric' output frame is Pseudo Earth Fixed (PEF; if cel_pole() was not set and DUT1 is 0) or ITRS (actual rotating Earth; if cel_pole() was set and ut1_to_tt includes the DUT1 component).
NOTES:
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
source | Pointer to a celestrial object data structure. Catalog objects musy have ICRS coordinates. You can use transform_cat() to convert other catalog systems to ICRS as necessary. | |
location | The observer location, relative to which the output positions and velocities are to be calculated | |
ut1_to_tt | [s] TT - UT1 time difference. Used only when 'location->where' is NOVAS_OBSERVER_ON_EARTH (1) or NOVAS_OBSERVER_IN_EARTH_ORBIT (2). | |
coord_sys | The coordinate system that defines the orientation of the celestial pole. If it is NOVAS_ICRS (3), a geometric position and radial velocity is returned. For all other systems, the returned position is the apparent position including aberration and gravitational deflection corrections, and the radial velocity is in the direction the eflected light was emitted from the source. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | output | Data structure to populate with the result. |
References aberration(), bary2obs(), d_light(), ephemeris(), era(), gcrs_to_cirs(), gcrs_to_j2000(), gcrs_to_mod(), gcrs_to_tod(), grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, grav_planets(), light_time2(), make_observer_at_geocenter(), NOVAS_BARYCENTER, NOVAS_CATALOG_OBJECT, NOVAS_CIRS, NOVAS_EARTH_INIT, NOVAS_FULL_ACCURACY, NOVAS_ICRS, NOVAS_ITRS, NOVAS_J2000, NOVAS_MOD, NOVAS_REDUCED_ACCURACY, NOVAS_SUN_INIT, NOVAS_TIRS, NOVAS_TOD, novas_vlen(), obs_planets(), obs_posvel(), proper_motion(), rad_vel2(), spin(), starvectors(), tt2tdb(), and vector2radec().
int place_icrs | ( | double | jd_tt, |
const object *restrict | source, | ||
enum novas_accuracy | accuracy, | ||
sky_pos *restrict | pos ) |
reference_system: ICRS / GCRS
observer location: geocenter
position_type: geometric
Computes the International Celestial Reference System (ICRS) position of a source. (from the geocenter). Unlike place_gcrs()
, this version does not include aberration or gravitational deflection corrections.
jd_tt | [day] Terrestrial Time (TT) based Julian date of observation. | |
source | Catalog source or solar_system body. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | Structure to populate with the calculated geometric ICRS position data for a fictitous geocentric observer. (Unlike place_gcrs(), the calculated coordinates do not account for aberration or gravitational deflection). |
References NOVAS_ICRS, and place().
int place_star | ( | double | jd_tt, |
const cat_entry *restrict | star, | ||
const observer *restrict | obs, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | system, | ||
enum novas_accuracy | accuracy, | ||
sky_pos *restrict | pos ) |
Computes the apparent place of a star at the specified date, given its catalog mean place, proper motion, parallax, and radial velocity.
See place()
for more information.
It is effectively the same as calling place(), except for the cat_entry
type target argument.
For all but ICRS coordinate outputs, the calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
obs | Observer location (NULL defaults to geocentric) | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
system | The type of coordinate reference system in which coordinates are to be returned. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | pos | The position and radial velocity of of the catalog source in the specified coordinate system and relative to the specified observer location (if applicable) |
References NOVAS_CATALOG_OBJECT, place(), object::star, and object::type.
double planet_lon | ( | double | t, |
enum novas_planet | planet ) |
Returns the planetary longitude, for Mercury through Neptune, w.r.t.
mean dynamical ecliptic and equinox of J2000, with high order terms omitted (Simon et al. 1994, 5.8.1-5.8.8).
REFERENCES:
t | [cy] Julian centuries since J2000 |
planet | Novas planet id, e.g. NOVAS_MARS. |
planet
id is out of range.References NOVAS_NEPTUNE, and TWOPI.
int radec_planet | ( | double | jd_tt, |
const object *restrict | ss_body, | ||
const observer *restrict | obs, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | sys, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ra, | ||
double *restrict | dec, | ||
double *restrict | dis, | ||
double *restrict | rv ) |
Computes the place of a solar system body at the specified time for an observer in the specified coordinate system.
This is the same as calling place() with the same arguments, except the different set of return values used.
Notwithstanding the different set of return values, this is the same as calling place_planet() with the same arguments.
For all but ICRS coordinate outputs, the calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
obs | Observer location. It may be NULL if not relevant. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
sys | Coordinate reference system in which to produce output values | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric apparent right ascension in hours, referred to the true equator and equinox of date, or NAN when returning with an error code. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric apparent declination in degrees referred to the true equator and equinox of date, or NAN when returning with an error code. (It may be NULL if not required) |
[out] | dis | [AU] True distance from Earth to the body at 'jd_tt' in AU, or NAN when returning with an error code. (It may be NULL if not needed). |
[out] | rv | [AU/day] radial velocity relative ot observer, or NAN when returning with an error code. (It may be NULL if not required) |
References sky_pos::dec, sky_pos::dis, NOVAS_EPHEM_OBJECT, NOVAS_ORBITAL_OBJECT, NOVAS_PLANET, place(), sky_pos::ra, sky_pos::rv, and SKY_POS_INIT.
int radec_star | ( | double | jd_tt, |
const cat_entry *restrict | star, | ||
const observer *restrict | obs, | ||
double | ut1_to_tt, | ||
enum novas_reference_system | sys, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ra, | ||
double *restrict | dec, | ||
double *restrict | rv ) |
Computes the place of a star at date 'jd_tt', for an observer in the specified coordinate system, given the star's ICRS catalog place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with the same arguments.
For all but ICRS coordinate outputs, the calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
obs | Observer location. It may be NULL if not relevant. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
sys | Coordinate reference system in which to produce output values | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric right ascension in hours, referred to true equator and equinox of date 'jd_tt' or NAN when returning with an error code. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric declination in degrees, referred to true equator and equinox of date 'jd_tt' or NAN when returning with an error code. (It may be NULL if not required) |
[out] | rv | [AU/day] radial velocity relative ot observer, or NAN when returning with an error code. (It may be NULL if not required) |
References sky_pos::dec, place_star(), sky_pos::ra, sky_pos::rv, and SKY_POS_INIT.
int set_cio_locator_file | ( | const char *restrict | filename | ) |
It used to set the CIO interpolaton data file to use to interpolate CIO locations vs the GCRS. As of version 1.5, this call does exactly nothing. It simply returns 0.
filename | (unused) Used to be the path (preferably absolute path) CIO_RA.TXT or else to the binary cio_ra.bin data, or NULL to disable using a CIO locator file altogether. |
short sidereal_time | ( | double | jd_ut1_high, |
double | jd_ut1_low, | ||
double | ut1_to_tt, | ||
enum novas_equinox_type | gst_type, | ||
enum novas_earth_rotation_measure | erot, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | gst ) |
Computes the Greenwich sidereal time, either mean or apparent, at the specified Julian date. The Julian date can be broken into two parts if convenient, but for the highest precision, set 'jd_high' to be the integral part of the Julian date, and set 'jd_low' to be the fractional part.
NOTES:
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. | |
jd_ut1_low | [day] Low-order part of UT1 Julian date. (You can leave it at zero if 'jd_high' specified the date with sufficient precision) | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
gst_type | NOVAS_MEAN_EQUINOX (0) or NOVAS_TRUE_EQUINOX (1), depending on whether wanting mean or apparent GST, respectively. | |
erot | Unused as of 1.5.0. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | gst | [h] Greenwich (mean or apparent) sidereal time, in hours [0:24]. (In case the returned error code is >1 the gst value will be set to NAN.) |
References NOVAS_FULL_ACCURACY, novas_gast(), novas_gmst(), NOVAS_REDUCED_ACCURACY, and NOVAS_TRUE_EQUINOX.
short ter2cel | ( | double | jd_ut1_high, |
double | jd_ut1_low, | ||
double | ut1_to_tt, | ||
enum novas_earth_rotation_measure | erot, | ||
enum novas_accuracy | accuracy, | ||
enum novas_equatorial_class | coordType, | ||
double | xp, | ||
double | yp, | ||
const double * | in, | ||
double * | out ) |
Rotates a vector from the terrestrial to the celestial system. Specifically, it transforms a vector in the ITRS (rotating earth-fixed system) to the True of Date (TOD), CIRS, or GCRS (a local space-fixed system) by applying rotations for polar motion, Earth rotation (for TOD); and nutation, precession, and the dynamical-to-GCRS frame tie (for GCRS).
If 'system' is NOVAS_CIRS then method EROT_ERA must be used. Similarly, if 'system' is NOVAS_TOD then method must be EROT_ERA. Otherwise an error 3 is returned.
If both 'xp' and 'yp' are set to 0 no polar motion is included in the transformation.
REFERENCES:
jd_ut1_high | [day] High-order part of UT1 Julian date. | |
jd_ut1_low | [day] Low-order part of UT1 Julian date. | |
ut1_to_tt | [s] TT - UT1 Time difference in seconds | |
erot | Unused. | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
coordType | Output coordinate class NOVAS_REFERENCE_CLASS (0, or any value other than 1) or NOVAS_DYNAMICAL_CLASS (1). Use the former if the output coordinates are to be in the GCRS, and the latter if they are to be in CIRS or TOD (the 'erot' parameter selects which dynamical system to use for the output.) | |
xp | [arcsec] Conventionally-defined X coordinate of celestial intermediate pole with respect to ITRS pole. If you have defined pole offsets the old (pre IAU2000) way, via cel_pole() , then use 0 here. | |
yp | [arcsec] Conventionally-defined Y coordinate of celestial intermediate pole with respect to ITRS pole. If you have defined pole offsets the old (pre IAU2000) way, via cel_pole() , then use 0 here. | |
in | Position vector, geocentric equatorial rectangular coordinates, referred to ITRS axes (terrestrial system) in the normal case where 'option' is NOVAS_GCRS (0). | |
[out] | out | Position vector, equatorial rectangular coordinates in the specified output system (GCRS if 'class' is NOVAS_REFERENCE_CLASS; or else either CIRS if 'erot' is EROT_ERA, or TOD if 'erot' is EROT_GST). It may be the same vector as the input. |
References cirs_to_gcrs(), era(), EROT_ERA, EROT_GST, NOVAS_DYNAMICAL_CLASS, NOVAS_FULL_ACCURACY, novas_gast(), NOVAS_REDUCED_ACCURACY, spin(), tod_to_gcrs(), wobble(), WOBBLE_ITRS_TO_PEF, and WOBBLE_ITRS_TO_TIRS.
short topo_planet | ( | double | jd_tt, |
const object *restrict | ss_body, | ||
double | ut1_to_tt, | ||
const on_surface *restrict | position, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ra, | ||
double *restrict | dec, | ||
double *restrict | dis ) |
reference_system: TOD
observer location: topocentric (on Earth)
position_type: apparent
Computes the topocentric apparent place of a solar system body at the specified time. This is the same as calling place() for the body for the same observer location and NOVAS_TOD as the reference system, except the different set of return values used.
The calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terretrial Time (TT) based Julian date. | |
ss_body | Pointer to structure containing the body designation for the solar system body. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric apparent right ascension, referred to the true equator and equinox of date. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric apparent declination, referred to the true equator and equinox of date. (It may be NULL if not required) |
[out] | dis | [AU] Apparent distance from Earth to the body at 'jd_tt' (may be NULL if not needed). |
References make_observer_at_site(), NOVAS_TOD, and radec_planet().
short topo_star | ( | double | jd_tt, |
double | ut1_to_tt, | ||
const cat_entry *restrict | star, | ||
const on_surface *restrict | position, | ||
enum novas_accuracy | accuracy, | ||
double *restrict | ra, | ||
double *restrict | dec ) |
reference_system: TOD
observer location: topocentric (on Earth)
position_type: apparent
Computes the topocentric (True of Date; TOD) apparent place of a star at date 'jd_tt', given its ICRS catalog place, proper motion, parallax, and radial velocity.
Notwithstanding the different set of return values, this is the same as calling place_star() with the same observer location and NOVAS_TOD for an object that specifies the star.
The calculated positions and velocities include aberration corrections for the moving frame of the observer as well as gravitational deflection due to the Sun and Earth and other major gravitating bodies in the Solar system, provided planet positions are available via a novas_planet_provider function.
REFERENCES:
jd_tt | [day] Terrestrial Time (TT) based Julian date. | |
ut1_to_tt | [s] Difference TT-UT1 at 'jd_tt', in seconds of time. | |
star | Pointer to catalog entry structure containing catalog data for the object in the ICRS. | |
position | Position of the observer | |
accuracy | NOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1) | |
[out] | ra | [h] Topocentric appatent right ascension, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required) |
[out] | dec | [deg] Topocentric apparent declination, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required) |
References make_observer_at_site(), NOVAS_TOD, and radec_star().