SuperNOVAS v1.5
The NOVAS C library, made better
Loading...
Searching...
No Matches
frames.c File Reference

High-level and efficient astrometric calculations using observing frames. More...

Functions

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_change_observer (const novas_frame *orig, const observer *obs, novas_frame *out)
 Change the observer location for an observing frame.
 
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.
 
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.
 
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.
 
int novas_invert_transform (const novas_transform *transform, novas_transform *inverse)
 Inverts a novas coordinate transformation matrix.
 
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_moon_angle (const object *restrict source, const novas_frame *restrict frame)
 Returns the apparent angular distance of a source from the Moon from the observer's point of view.
 
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.
 
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_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_solar_illum (const object *restrict source, const novas_frame *restrict frame)
 Returns the Solar illumination fraction of a source, assuming a spherical geometry for the observed body.
 
double novas_sun_angle (const object *restrict source, const novas_frame *restrict frame)
 Returns the apparent angular distance of a source from the Sun from the observer's point of view.
 
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.
 

Detailed Description

High-level and efficient astrometric calculations using observing frames.

An observing frame represents an observer location at a specific astronomical time (instant), which can be re-used again and again to calculate or transform positions of celestial sources in a a range of astronomical coordinate systems.

To use frames, you start with novas_make_frame() with an astrometric time and an observer location, and optionally Earth-orientation parameters if precision below the arcsecond level is needed. E.g.:

observer obs; // observer location
novas_timespec time; // astrometric time
novas_frame frame; // observing frame
// Specify an observer, with GPS coordinates
make_gps_observer(gps_lat, gps_lon, gps_alt, &obs);
// Specify time, e.g. current time with leap seconds and UT1-UTC time difference:
novas_set_current_time(leap_seconds, dut1, &time);
// Set up a frame with the desired accuracy and Earth orientation parameters xp, yp.
novas_make_frame(NOVAS_REDUCED_ACCURACY, &obs, &time, xp, yp, &frame);
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 require...
Definition frames.c:515
@ NOVAS_REDUCED_ACCURACY
Calculate with truncated terms.
Definition novas.h:543
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 loc...
Definition observer.c:277
A set of parameters that uniquely define the place and time of observation.
Definition novas.h:1375
A structure, which defines a precise instant of time that can be extpressed in any of the astronomica...
Definition novas.h:1285
Observer location.
Definition novas.h:1182
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.
Definition timescale.c:889

Once a frame is defined, you can perform various astrometric calculation within it efficiently. For example, calculate apparent positions, in the coordinate system of choice, to predict where objects would be seen by the observer:

object source = ...; // observed source
sky_pos app; // apparent position to calculate.
// Calculate True-of-Date apparent coordinates for a source.
novas_sky_pos(&source, &frame, NOVAS_TOD, &app);
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.
Definition frames.c:821
@ NOVAS_TOD
True equinox Of Date: dynamical system of the 'true' equator, with its origin at the 'true' equinox (...
Definition novas.h:449
Celestial object's place on the sky; contains the output from place()
Definition novas.h:1211

If your observer was defined on Earth (ground-based or airborne), you might continue to calculate Az/El positions using novas_app_to_hor(), or you might just want to know the same position in another reference system also:

novas_transform T; // Transformation between reference systems
sky_pos tirs; // calculated apparent position in TIRS
// Transform from True-of-Date to TIRS coordinates
// Transform the above position to TIRS...
novas_transform_sky_pos(&app, T, &tirs);
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 coor...
Definition frames.c:1234
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.
Definition frames.c:1437
@ NOVAS_TIRS
Terrestrial Intermediate Reference System.
Definition novas.h:473
A transformation between two astronomical coordinate systems for the same observer location and time.
Definition novas.h:1432

Or, you can calculate geometric positions, which are not corrected for aberration or gravitational deflection (but still corrected for light travel time in case of Solar-system sources):

object source = ...; // observed source
double pos[3], vel[3]; // geometric position and velocity vectors to calculate.
// Calculate geometric positions / velocities, say in ICRS
novas_geom_posvel(&source, &frame, NOVAS_ICRS, pos, vel);
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...
Definition frames.c:707
@ NOVAS_ICRS
International Celestial Reference system.
Definition novas.h:457

You can also transform geometric positions / velocities to other reference frames using the same novas_transform as above.

Or, perhaps you are interested when the source will rise (or transit, or set) next, for a ground-based or airborne observer:

object source = ...; // observed source
// UTC-based Julian date when source rises above 20 degrees given a refraction model.
double jd_rise = novas_rises_above(20.0, &source, &frame, novas_standard_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.
Definition frames.c:1654
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.
Definition refract.c:258

Or, you might simply want to know how far your source is from the Sun or Moon (or another source) on the sky:

object source = ...; // observed source
// Calculate the angular distance, in degrees, of the source from the Sun.
double angle_deg = novas_sun_angle(&source, &frame);
double novas_sun_angle(const object *restrict source, const novas_frame *restrict frame)
Returns the apparent angular distance of a source from the Sun from the observer's point of view.
Definition frames.c:1823

These are just some of the common use-case scenarios. There is even more possibilities with frames...

Date
Created on Jun 23, 2024
Author
Attila Kovacs
Since
1.1
See also
timescale.c, observer.c, target.c, ephemeris.c

Function Documentation

◆ novas_app_to_geom()

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.

Parameters
frameThe observer frame, defining the location and time of observation
sysThe reference system in which the observed position is specified.
ra[h] Observed ICRS right-ascension of the source
dec[deg] Observed ICRS declination of the source
dist[AU] Observed distance from observer. A value of <=0 will translate to 1015 AU (around 5 Gpc).
[out]geom_icrs[AU] The corresponding geometric position for the source, in ICRS.
Returns
0 if successful, or else an error from grav_undef2(), or -1 (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_geom_to_app()
novas_hor_to_app(), novas_geom_to_hor(), novas_transform_vector()

References grav_undo_planets(), NOVAS_CIRS, NOVAS_ITRS, NOVAS_J2000, NOVAS_MOD, NOVAS_REFERENCE_SYSTEMS, NOVAS_TIRS, NOVAS_TOD, radec2vector(), spin(), wobble(), and WOBBLE_ITRS_TO_TIRS.

◆ novas_app_to_hor()

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.

The observer must be located on the surface of Earth, or else the call will return with an error. The caller may optionally supply a refraction model of choice to calculate an appropriate elevation angle that includes a refraction correction for Earth's atmosphere. If no such model is provided the calculated elevation will be the astrometric elevation without a refraction correction.

Parameters
frameObserver frame, defining the time and place of observation (on Earth).
sysAstronomical coordinate system in which the observed position is given.
ra[h] Observed apparent right ascension (R.A.) coordinate
dec[deg] Observed apparent declination coordinate
ref_modelAn appropriate refraction model, or NULL to calculate unrefracted elevation. Depending on the refraction model, you might want to make sure that the weather parameters were set when the observing frame was defined.
[out]az[deg] Calculated azimuth angle. It may be NULL if not required.
[out]el[deg] Calculated elevation angle. It may be NULL if not required.
Returns
0 if successful, or else an error from tod_to_itrs() or cirs_to_itrs(), or -1 (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_hor_to_app()
novas_app_to_geom(), novas_standard_refraction(), novas_optical_refraction(), novas_radio_refraction(), novas_wave_refraction()

References novas_timespec::fjd_tt, novas_timespec::ijd_tt, itrs_to_hor(), NOVAS_AIRBORNE_OBSERVER, NOVAS_CIRS, NOVAS_GCRS, NOVAS_ICRS, NOVAS_ITRS, NOVAS_J2000, NOVAS_MOD, NOVAS_OBSERVER_ON_EARTH, NOVAS_REFRACT_ASTROMETRIC, NOVAS_TIRS, NOVAS_TOD, radec2vector(), spin(), wobble(), WOBBLE_PEF_TO_ITRS, and WOBBLE_TIRS_TO_ITRS.

◆ novas_change_observer()

int novas_change_observer ( const novas_frame * orig,
const observer * obs,
novas_frame * out )

Change the observer location for an observing frame.

Parameters
origPointer to original observing frame
obsNew observer location
[out]outObserving frame to populate with a original frame data and new observer location. It can be the same as the input.
Returns
0 if successfule or else an an error code from geo_posvel() (errno will also indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_make_frame()

References novas_frame::accuracy, grav_bodies_full_accuracy, grav_bodies_reduced_accuracy, NOVAS_FULL_ACCURACY, novas_get_time(), NOVAS_TDB, obs_planets(), novas_frame::obs_pos, novas_frame::observer, novas_frame::planets, novas_frame::state, and novas_frame::time.

◆ novas_equ_track()

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.

The position and its derivatives are calculated via the more precise IAU2006 method, and CIRS.

Parameters
sourceObserved source
frameObserving frame, defining the observer location and astronomical time of observation.
dt[s] Time step used for calculating derivatives.
[out]trackOutput tracking parameters to populate
Returns
0 if successful, or else -1 if any of the pointer arguments are NULL, or else an error code from novas_sky_pos().
Since
1.3
Author
Attila Kovacs
See also
novas_hor_track(), novas_track_pos()

References sky_pos::dec, sky_pos::dis, novas_timespec::fjd_tt, ira_equinox(), NOVAS_CIRS, novas_make_frame(), novas_sky_pos(), NOVAS_TRUE_EQUINOX, novas_v2z(), sky_pos::ra, sky_pos::rv, and SKY_POS_INIT.

◆ novas_frame_lst()

double novas_frame_lst ( const novas_frame *restrict frame)

Returns the Local (apparent) Sidereal Time for an observing frame of an Earth-bound observer.

Parameters
frameObserver frame, defining the location and time of observation
Returns
[h] The LST for an Earth-bound observer [0.0–24.0), or NAN otherwise. If NAN is returned errno will indicate the type of error.
Since
1.3
Author
Attila Kovacs
See also
novas_time_lst()

References NOVAS_AIRBORNE_OBSERVER, and NOVAS_OBSERVER_ON_EARTH.

◆ novas_geom_posvel()

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.

The geometric position includes proper motion, and for solar-system bodies it is antedated for light travel time, so it effectively represents the geometric position as seen by the observer. However, the geometric does not include aberration correction, nor gravitational deflection.

If you want apparent positions, which account for aberration and gravitational deflection, use novas_skypos() instead.

You can also use novas_transform_vector() to convert the output position and velocity vectors to a dfferent coordinate system of choice afterwards if you want the results expressed in more than one coordinate system.

It implements the same geometric transformations as place() but at a reduced computational cost. See place() for references.

NOTES:

  1. As of SuperNOVAS v1.3, the returned velocity vector is a proper observer-based velocity measure. In prior releases, and in NOVAS C 3.1, this was inconsistent, with pseudo LSR-based measures being returned for catalog sources.
Parameters
sourcePointer to a celestial source data structure that is observed. Catalog sources should have coordinates and properties in ICRS. You can use transform_cat() to convert catalog entries to ICRS as necessary.
frameObserver frame, defining the location and time of observation
sysThe coordinate system in which to return positions and velocities.
[out]pos[AU] Calculated geometric position vector of the source relative to the observer location, in the designated coordinate system. It may be NULL if not required.
[out]vel[AU/day] The calculated velocity vector of the source relative to the observer in the designated coordinate system. It must be distinct from the pos output vector, and may be NULL if not required.
Returns
0 if successful, or else -1 if any of the arguments is invalid, 50–70 error is 50 + error from light_time2().
Since
1.1
Author
Attila Kovacs
See also
novas_geom_to_app(), novas_sky_pos(),
novas_transform_vector(), novas_make_frame()

References bary2obs(), d_light(), light_time2(), NOVAS_CATALOG_OBJECT, NOVAS_FULL_ACCURACY, novas_get_time(), NOVAS_JD_J2000, NOVAS_PLANET, NOVAS_REDUCED_ACCURACY, NOVAS_TDB, proper_motion(), and starvectors().

◆ novas_geom_to_app()

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.

Unlike place() the output reports the distance calculated from the parallax for sidereal sources. The radial velocity of the output is set to NAN (if needed use novas_sky_pos() instead).

Parameters
frameThe observer frame, defining the location and time of observation
pos[AU] Geometric position of source in ICRS coordinates
sysThe coordinate system in which to return the apparent sky location
[out]outPointer to the data structure which is populated with the calculated apparent location in the designated coordinate system. It may be the same pounter as the input position.
Returns
0 if successful, or an error from grav_def2(), or else -1 (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_app_to_geom()
novas_app_to_hor(), novas_sky_pos(), novas_geom_posvel()

References grav_planets(), NOVAS_FULL_ACCURACY, NOVAS_REDUCED_ACCURACY, novas_vlen(), and vector2radec().

◆ novas_hor_to_app()

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.

The observer must be located on the surface of Earth, or else the call will return with an error. The caller may optionally supply a refraction model of choice to calculate an appropriate elevation angle that includes a refraction correction for Earth's atmosphere. If no such model is provided, the provided elevation value will be assumed to be an astrometric elevation without a refraction correction.

Parameters
frameObserver frame, defining the time and place of observation (on Earth).
az[deg] Observed azimuth angle. It may be NULL if not required.
el[deg] Observed elevation angle. It may be NULL if not required.
ref_modelAn appropriate refraction model, or NULL to assume unrefracted elevation. Depending on the refraction model, you might want to make sure that the weather parameters were set when the observing frame was defined.
sysAstronomical coordinate system in which the output is R.A. and declination values are to be calculated.
[out]ra[h] Calculated apparent right ascension (R.A.) coordinate
[out]dec[deg] Calculated apparent declination coordinate
Returns
0 if successful, or else an error from itrs_to_tod() or itrs_to_cirs(), or -1 (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_app_to_hor()
novas_app_to_geom(), novas_standard_refraction(), novas_optical_refraction(), novas_radio_refraction(), novas_wave_refraction()

References novas_timespec::fjd_tt, hor_to_itrs(), novas_timespec::ijd_tt, NOVAS_AIRBORNE_OBSERVER, NOVAS_CIRS, NOVAS_GCRS, NOVAS_ICRS, NOVAS_ITRS, NOVAS_J2000, NOVAS_MOD, NOVAS_OBSERVER_ON_EARTH, NOVAS_REFRACT_OBSERVED, NOVAS_TIRS, NOVAS_TOD, spin(), vector2radec(), wobble(), WOBBLE_ITRS_TO_PEF, and WOBBLE_ITRS_TO_TIRS.

◆ novas_hor_track()

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.

The position and its derivatives are calculated via the more precise IAU2006 method, and CIRS, and then converted to local horizontal coordinates using the specified refraction model (if any).

Parameters
sourceObserved source
frameObserving frame, defining the observer location and astronomical time of observation.
ref_modelRefraction model to use, or NULL for an unrefracted track.
[out]trackOutput tracking parameters to populate
Returns
0 if successful, or else -1 if any of the pointer arguments are NULL, or else an error code from novas_sky_pos(), or from novas_app_hor().
Since
1.3
Author
Attila Kovacs
See also
novas_equ_track(), novas_track_pos()

References sky_pos::dec, sky_pos::dis, novas_timespec::fjd_tt, ira_equinox(), NOVAS_AIRBORNE_OBSERVER, novas_app_to_hor(), NOVAS_CIRS, novas_make_frame(), NOVAS_OBSERVER_ON_EARTH, novas_sky_pos(), NOVAS_TOD, NOVAS_TRUE_EQUINOX, novas_v2z(), sky_pos::ra, sky_pos::rv, and SKY_POS_INIT.

◆ novas_invert_transform()

int novas_invert_transform ( const novas_transform * transform,
novas_transform * inverse )

Inverts a novas coordinate transformation matrix.

Parameters
transformPointer to a coordinate transformation matrix.
[out]inversePointer to a coordinate transformation matrix to populate with the inverse transform. It may be the same as the input.
Returns
0 if successful, or else -1 if the was an error (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_make_transform()

References novas_transform::matrix.

◆ novas_make_frame()

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.

The frame is initialized using the currently configured planet ephemeris provider function (see set_planet_provider() and set_planet_provider_hp()), and in case of reduced accuracy mode, the currently configured IAU nutation model provider (see set_nutation_lp_provider()).

Note, that to construct full accuracy frames, you will need a high-precision ephemeris provider for the major planets (not just the default Earth/Sun), as without it, gravitational bending around massive plannets cannot be accounted for, and therefore μas accuracy cannot be ensured, in general. Attempting to construct a high-accuracy frame without a high-precision ephemeris provider for the major planets will result in an error in the 10–40 range from the required ephemeris() call.

NOTES:

  1. This function expects the Earth polar wobble parameters to be defined on a per-frame basis and will not use the legacy global (undated) orientation parameters set via cel_pole().
  2. The Earth orientation parameters xp, yp should be provided in the same ITRF realization as the observer location for an Earth-based observer. You can use novas_itrf_transform_eop() to convert the EOP values as necessary.
Parameters
accuracyAccuracy requirement, NOVAS_FULL_ACCURACY (0) for the utmost precision or NOVAS_REDUCED_ACCURACY (1) if ~1 mas accuracy is sufficient.
obsObserver location
timeTime of observation
xp[mas] Earth orientation parameter, polar offset in x, e.g. from the IERS Bulletins, and possibly corrected for diurnal and semi-diurnal variations, e.g. via novas_diurnal_eop(). (The global, undated value set by cel_pole() is not not used here.) You can use 0.0 if sub-arcsecond accuracy is not required.
yp[mas] Earth orientation parameter, polar offset in y, e.g. from the IERS Bulletins, and possibly corrected for diurnal and semi-diurnal variations, e.g. via novas_diurnal_eop(). (The global, undated value set by cel_pole() is not not used here.) You can use 0.0 if sub-arcsecond accuracy is not required.
[out]framePointer to the observing frame to configure.
Returns
0 if successful, 10–40: error is 10 + the error from ephemeris(), 40–50: error is 40 + the error from geo_posvel(), or else -1 if there was an error (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_change_observer(), novas_sky_pos(), novas_geom_posvel(), novas_make_transform()
set_planet_provider(), set_planet_provider_hp(), set_nutation_lp_provider(), novas_diurnal_eop(), novas_itrf_transform_eop()

References novas_frame::accuracy, novas_frame::deps0, novas_frame::dpsi0, novas_frame::dx, novas_frame::dy, novas_frame::earth_pos, novas_frame::earth_vel, novas_frame::ee, ephemeris(), era(), novas_frame::era, novas_timespec::fjd_tt, novas_frame::gst, novas_timespec::ijd_tt, mean_obliq(), novas_frame::mobl, NOVAS_BARYCENTER, novas_change_observer(), NOVAS_EARTH_INIT, novas_get_split_time(), NOVAS_JD_J2000, NOVAS_OBSERVER_PLACES, NOVAS_REDUCED_ACCURACY, NOVAS_SUN_INIT, NOVAS_UT1, nutation_angles(), novas_frame::state, novas_frame::sun_pos, novas_frame::sun_vel, novas_frame::time, novas_frame::tobl, tt2tdb(), and observer::where.

◆ novas_make_transform()

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.

Parameters
frameObserver frame, defining the location and time of observation
from_systemOriginal coordinate reference system
to_systemNew coordinate reference system
[out]transformPointer to the transform data structure to populate.
Returns
0 if successful, or else -1 if there was an error (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_transform_vector(), novas_transform_sky_pos(), novas_invert_transform()
novas_geom_posvel(), novas_app_to_geom()

References novas_frame::era, novas_transform::frame, novas_transform::from_system, novas_frame::gcrs_to_cirs, novas_frame::icrs_to_j2000, novas_matrix::M, novas_transform::matrix, NOVAS_CIRS, NOVAS_GCRS, NOVAS_ICRS, NOVAS_ITRS, NOVAS_J2000, NOVAS_MOD, NOVAS_REFERENCE_SYSTEMS, NOVAS_TIRS, NOVAS_TOD, novas_frame::nutation, novas_frame::precession, and novas_transform::to_system.

◆ novas_moon_angle()

double novas_moon_angle ( const object *restrict source,
const novas_frame *restrict frame )

Returns the apparent angular distance of a source from the Moon from the observer's point of view.

Parameters
sourceAn observed source
frameObserving frame, defining the observer location and astronomical time of observation.
Returns
[deg] Apparent angular distance between the source an the Moon, from the observer's point of view.
Since
1.3
Author
Attila Kovacs
See also
novas_sun_angle()

References NOVAS_MOON_INIT, and novas_object_sep().

◆ novas_object_sep()

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.

The calculated separation includes light-time corrections, aberration and gravitational deflection for both sources, and thus represents a precise observed separation between the two sources.

Parameters
source1An observed source
source2Another observed source
frameObserving frame, defining the observer location and astronomical time of observation.
Returns
[deg] Apparent angular separation between the two observed source from the observer's point-of-view.
Since
1.3
Author
Attila Kovacs
See also
novas_sep(), novas_sun_angle(), novas_moon_angle()

References sky_pos::dec, sky_pos::dis, novas_equ_sep(), NOVAS_GCRS, novas_sky_pos(), sky_pos::ra, and SKY_POS_INIT.

◆ novas_rises_above()

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.

The calculated time will account for the (slow) motion for Solar-system bodies, and optionally for atmospheric refraction also.

NOTES:

  1. The current implementation is not suitable for calculating the nearest successive rise times for near-Earth objects, at or within the geostationary orbit.
  2. This function calculates the time when the center (not the limb!) of the source rises above the specified elevation threshold. Something to keep in mind for calculating Sun/Moon rise times.
Parameters
el[deg] Elevation angle.
sourceObserved source
frameObserving frame, defining the observer location and astronomical time of observation.
ref_modelRefraction model, or NULL to calculate unrefracted rise time.
Returns
[day] UTC-based Julian date at which the object rises above the specified elevation next after the specified date, or else NAN if the source stays above or below the given elevation for the entire 24-hour period.
Since
1.3
Author
Attila Kovacs
See also
novas_sets_below(), novas_transit_time()

◆ novas_sets_below()

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.

The calculated time will account for the (slow) motion of Solar-system bodies, and optionally for atmopsheric refraction also.

NOTES:

  1. The current implementation is not suitable for calculating the nearest successive set times for near-Earth objects, at or within the geostationary orbit.
  2. This function calculates the time when the center (not the limb!) of the source sets below the specified elevation threshold. Something to keep in mind for calculating Sun/Moon rise times.
Parameters
el[deg] Elevation angle.
sourceObserved source
frameObserving frame, defining the observer location and astronomical time of observation.
ref_modelRefraction model, or NULL to calculate unrefracted setting time.
Returns
[day] UTC-based Julian date at which the object sets below the specified elevation next after the specified date, or else NAN if the source stays above or below the given elevation for the entire 24-hour day..
Since
1.3
Author
Attila Kovacs
See also
novas_rises_above(), novas_transit_time()

◆ novas_sky_pos()

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.

The position takes into account the proper motion (for sidereal source), or is antedated for light-travel time (for Solar-System bodies). It also applies an appropriate aberration correction and gravitational deflection of the light.

To calculate corresponding local horizontal coordinates, you can pass the output RA/Dec coordinates to novas_app_to_hor(). Or to calculate apparent coordinates in other systems, you may pass the result to novas_transform_sy_pos() after.

And if you want geometric positions instead (not corrected for aberration or gravitational deflection), you may want to use novas_geom_posvel() instead.

The approximate 'inverse' of this function is novas_app_to_geom().

This function implements the same aberration and gravitational deflection corrections as place(), but at reduced computational cost. See place() for references. Unlike place(), however, the output always reports the distance calculated from the parallax for sidereal sources. Note also, that while place() does not apply aberration and gravitational deflection corrections when sys is NOVAS_ICRS (3), this routine will apply those corrections consistently for all coordinate systems (and you can use novas_geom_posvel() instead to get positions without aberration or deflection in any system).

NOTES:

  1. As of SuperNOVAS v1.3, the returned radial velocity component is a proper observer-based spectroscopic measure. In prior releases, and in NOVAS C 3.1, this was inconsistent, with LSR-based measures being returned for catalog sources.
Parameters
objectPointer to a celestial object data structure that is observed. Catalog sources should have coordinates and properties in ICRS. You can use transform_cat() to convert catalog entries to ICRS as necessary.
frameThe observer frame, defining the location and time of observation.
sysThe coordinate system in which to return the apparent sky location.
[out]outPointer to the data structure which is populated with the calculated apparent location in the designated coordinate system.
Returns
0 if successful, 50–70 error is 50 + error from light_time2(), 70–80 error is 70 + error from grav_def(), or else -1 (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_geom_to_app(), novas_app_to_hor(), novas_make_frame()

References grav_planets(), NOVAS_CATALOG_OBJECT, NOVAS_FULL_ACCURACY, novas_geom_posvel(), novas_geom_to_app(), NOVAS_ICRS, NOVAS_REDUCED_ACCURACY, novas_vlen(), rad_vel2(), and object::type.

◆ novas_solar_illum()

double novas_solar_illum ( const object *restrict source,
const novas_frame *restrict frame )

Returns the Solar illumination fraction of a source, assuming a spherical geometry for the observed body.

Parameters
sourceObserved source. Usually a Solar-system source. (For other source types, 1.0 is returned by default.)
frameObserving frame, defining the observer location and astronomical time of observation.
Returns
Solar illumination fraction [0.0:1.0] of a spherical body observed at the source location from the given observer location, or NAN if there was an error (errno will indicate the type of error).
Since
1.3
Author
Attila Kovacs
See also
novas_solar_power()

References NOVAS_CATALOG_OBJECT, novas_geom_posvel(), NOVAS_ICRS, and novas_vlen().

◆ novas_sun_angle()

double novas_sun_angle ( const object *restrict source,
const novas_frame *restrict frame )

Returns the apparent angular distance of a source from the Sun from the observer's point of view.

Parameters
sourceAn observed source
frameObserving frame, defining the observer location and astronomical time of observation.
Returns
[deg] the apparent angular distance between the source an the Sun, from the observer's point of view.
Since
1.3
Author
Attila Kovacs
See also
novas_moon_angle()

References novas_object_sep(), and NOVAS_SUN_INIT.

◆ novas_track_pos()

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.

Thus if the input track was obtained with novas_hor_track() it will calculate momentary azimuth and elevation (Az/El) for the specified proximate time. And, if you used novas_equ_track() then it will give you theinstantaneous R.A. and declination for that time.

Using the quadratic trajectories via a novas_track to project positions can be much faster than the repeated full recalculation of the source position, but still quite accurate over a suffciently brief period.

In SuperNOVAS terminology a 'track' is a 2nd order Taylor series expansion of the observed position and redshift in time. For most but the fastest moving sources, horizontal (Az/El) tracks are sufficiently precise for telescope positioning on minute timescales, whereas depending on the type of source, equatorial tracks can be precise for up to days, especially for sidereal (non-Solar-system) sources.

Parameters
trackTracking position and motion (first and second derivatives), obtained e.g. in horizontal (Az/El) coordinates with novas_hor_track() or in equatorial coordinates with novas_equ_track().
timeAstrometric time of observation for which to calculate projected positions, around the time for which the track was obtained.
[out]lon[deg] projected observed Eastward longitude in the tracking coordinate system. Thus, if the track was calculated with novas_hor_track() it will return the projected azimuth coordinate for the given time. And, if you used novas_equ_track(), then this will be R.A. (in degrees, nevertheless).
[out]lat[deg] projected observed latitude in the tracking coordinate system. Thus, if the track was calculated with novas_hor_track() it will return the projected elevation angle for the given time. And, if you used novas_equ_track(), then this will be declination.
[out]dist[AU] projected apparent distance to source from observer.
[out]zprojected observed redshift (z = Δλ / λrest, e.g. for spectroscopy).
Returns
0 if successful, or else -1 if either input pointer is NULL (errno is set to EINVAL).
Since
1.3
Author
Attila Kovacs
See also
novas_equ_track(), novas_hor_track()
novas_z2v()

References novas_track::accel, novas_observable::dist, novas_observable::lat, novas_observable::lon, novas_diff_time(), novas_track::pos, novas_track::rate, novas_track::time, and novas_observable::z.

◆ novas_transform_sky_pos()

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.

Parameters
inInput apparent position on sky in the original coordinate reference system
transformPointer to a coordinate transformation matrix
[out]outOutput apparent position on sky in the new coordinate reference system. It may be the same as the input.
Returns
0 if successful, or else -1 if there was an error (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_sky_pos(), novas_make_transform(), novas_transform_vector()

References sky_pos::dec, sky_pos::r_hat, sky_pos::ra, and vector2radec().

◆ novas_transform_vector()

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.

Parameters
inInput 3-vector in the original coordinate reference system
transformPointer to a coordinate transformation matrix
[out]outOutput 3-vector in the new coordinate reference system. It may be the same as the input.
Returns
0 if successful, or else -1 if there was an error (errno will indicate the type of error).
Since
1.1
Author
Attila Kovacs
See also
novas_make_transform(), novas_transform_skypos()

◆ novas_transit_time()

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.

The calculated time will account for the (slow) motion of Solar-system bodies.

NOTES:

  1. The current implementation is not suitable for calculating the nearest successive transit times for near-Earth objects, at or within the geostationary orbit.
Parameters
sourceObserved source
frameObserving frame, defining the observer location and astronomical time of observation.
Returns
[day] UTC-based Julian date at which the object transits the local meridian next after the specified date, or NAN if either input pointer is NULL.
Since
1.3
Author
Attila Kovacs
See also
novas_rises_above(), novas_sets_below()