SuperNOVAS v1.5
The NOVAS C library, made better
Loading...
Searching...
No Matches
Apparent equatorial positions on sky

Data Structures

struct  novas_observable
 Spherical and spectral coordinate set. More...
 
struct  novas_track
 The spherical and spectral tracking position of a source, and its first and second time derivatives. More...
 
struct  sky_pos
 Celestial object's place on the sky; contains the output from place() More...
 

Macros

#define DEFAULT_GRAV_BODIES_FULL_ACCURACY   ( DEFAULT_GRAV_BODIES_REDUCED_ACCURACY | (1 << NOVAS_JUPITER) | (1 << NOVAS_SATURN) )
 Default set of gravitating bodies to use for deflection calculations in full accuracy mode.
 
#define DEFAULT_GRAV_BODIES_REDUCED_ACCURACY   ( (1 << NOVAS_SUN) )
 Default set of gravitating bodies to use for deflection calculations in reduced accuracy mode.
 
#define NOVAS_OBSERVABLE_INIT
 Empty initializer for novas_observable.
 
#define NOVAS_TRACK_INIT
 Empty initializer for novas_track.
 
#define SKY_POS_INIT
 Initializer for a NOVAS sky_pos structure.
 

Functions

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

 
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).
 
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_approx_sky_pos (enum novas_planet id, const novas_frame *restrict frame, enum novas_reference_system sys, sky_pos *restrict out)
 Calculates an approximate apparent location on sky for a major planet, Sun, Moon, Earth-Moon Barycenter (EMB) – typically to arcmin level accuracy – using Keplerian orbital elements.
 
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.
 
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_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_moon_phase (double jd_tdb)
 Calculates the Moon's phase at a given time.
 
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_sep (double lon1, double lat1, double lon2, double lat2)
 Returns the angular separation of two locations on a sphere.
 
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 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_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_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

 
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

 

Variables

int grav_bodies_full_accuracy
 Current set of gravitating bodies to use for deflection calculations in full accuracy mode.
 
int grav_bodies_reduced_accuracy
 Current set of gravitating bodies to use for deflection calculations in reduced accuracy mode.
 

Detailed Description

Macro Definition Documentation

◆ DEFAULT_GRAV_BODIES_FULL_ACCURACY

#define DEFAULT_GRAV_BODIES_FULL_ACCURACY   ( DEFAULT_GRAV_BODIES_REDUCED_ACCURACY | (1 << NOVAS_JUPITER) | (1 << NOVAS_SATURN) )

Default set of gravitating bodies to use for deflection calculations in full accuracy mode.

Since
1.1
Author
Attila Kovacs
See also
grav_bodies_full_accuracy, grav_planets(), grav_undo_planets()

◆ DEFAULT_GRAV_BODIES_REDUCED_ACCURACY

#define DEFAULT_GRAV_BODIES_REDUCED_ACCURACY   ( (1 << NOVAS_SUN) )

Default set of gravitating bodies to use for deflection calculations in reduced accuracy mode.

(only apply gravitational deflection for the Sun.)

Since
1.1
Author
Attila Kovacs
See also
grav_bodies_reduced_accuracy, grav_planets(), grav_undo_planets()

◆ NOVAS_OBSERVABLE_INIT

#define NOVAS_OBSERVABLE_INIT

Empty initializer for novas_observable.

Since
1.3
See also
novas_observable

◆ NOVAS_TRACK_INIT

#define NOVAS_TRACK_INIT

Empty initializer for novas_track.

Since
1.3
See also
novas_track

◆ SKY_POS_INIT

#define SKY_POS_INIT

Initializer for a NOVAS sky_pos structure.

Since
1.1.1
Author
Attila Kovacs
See also
sky_pos

Function Documentation

◆ app_planet()

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

Computes the True-of-Date (TOD) apparent place of a solar system body as would be seen by an observer at the geocenter. This is the same as calling place() for the body with NOVAS_TOD as the system and a geocentric observer, 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:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Explanatory Supplement to the Astronomical Almanac (1992),Chapter 3.
Parameters
jd_tt[day] Terretrial Time (TT) based Julian date.
ss_bodyPointer to structure containing the body designation for the solar system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]ra[h] Apparent (TOD) right ascension for a fictitous geocentric observer, referred to true equator and equinox of date 'jd_tt'. (It may be NULL if not required)
[out]dec[deg] Apparent (TOD) declination for a fictitous geocentric observer, referred to true equator and equinox of date 'jd_tt'. (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 needed).
Returns
0 if successful, or -1 if the object argument is NULL, or else 1 if the value of 'type' in structure 'ss_body' is invalid, or 10 + the error from place().
See also
place_tod(), astro_planet(), local_planet(), topo_planet(), virtual_planet(), app_star()
novas_sky_pos(), make_observer_at_geocenter(), NOVAS_TOD

References NOVAS_TOD, and radec_planet().

◆ app_star()

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

Computes the True-of-Date (TOD) apparent place of a star, as would be seen by an observer at the geocenter, referenced to the dynamical equator of date, 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 a NULL observer location and NOVAS_TOD as the system 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:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Explanatory Supplement to the Astronomical Almanac (1992), Chapter 3.
Parameters
jd_tt[day] Terretrial Time (TT) based Julian date.
starPointer to catalog entry structure containing catalog data for the object in the ICRS.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]ra[h] Apparent (TOD) right ascension for a fictitous geocentric observer, referred to true equator and equinox of date 'jd_tt' (it may be NULL if not required).
[out]dec[deg] Apparent (TOD) declination in degrees for a fictitous geocentric observer, referred to true equator and equinox of date 'jd_tt' (it may be NULL if not required).
Returns
0 if successful, -1 if a required pointer argument is NULL, or else an the error from make_object(), or 20 + the error from place().
See also
place_tod(), novas_sky_pos(), place_star(), astro_star(), local_star(), topo_star(), virtual_star(), app_planet()
novas_sky_pos(), make_observer_at_geocenter(), NOVAS_TOD

References NOVAS_TOD, and radec_star().

◆ limb_angle()

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).

The geometric limb is computed, assuming the Earth to be an airless sphere (no refraction or oblateness is included). The observer can be on or above the Earth. For an observer on the surface of the Earth, this function returns the approximate unrefracted elevation.

Parameters
pos_src[AU] Position 3-vector of observed object, with respect to origin at geocenter, components in AU.
pos_obs[AU] Position 3-vector of observer, with respect to origin at geocenter, components in AU.
[out]limb_ang[deg] Angle of observed object above (+) or below (-) limb in degrees, or NAN if reurning with an error. It may be NULL if not required.
[out]nadir_angNadir angle of observed object as a fraction of apparent radius of limb: lt;1.0 if below the limb; 1.0 on the limb; or >1.0 if above the limb. Returns NAN in case of an error return. It may be NULL if not required.
Returns
0 if successful, or -1 if either of the input vectors is NULL or if either input position is a null vector (at the geocenter).

References M_PI, and novas_vlen().

◆ 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_approx_sky_pos()

int novas_approx_sky_pos ( enum novas_planet id,
const novas_frame *restrict frame,
enum novas_reference_system sys,
sky_pos *restrict out )

Calculates an approximate apparent location on sky for a major planet, Sun, Moon, Earth-Moon Barycenter (EMB) – typically to arcmin level accuracy – using Keplerian orbital elements.

The returned position is antedated for light-travel time (for Solar-System bodies). It also applies an appropriate aberration correction (but not gravitational deflection).

The orbitals can provide planet positions to arcmin-level precision for the rocky inner planets, and to a fraction of a degree precision for the gas and ice giants and Pluto. The accuracies for Uranus, Neptune, and Pluto are significantly improved (to the arcmin level) if used in the time range of 1800 AD to 2050 AD. For a more detailed summary of the typical accuracies, see either of the top two references below.

For accurate positions, you should use planetary ephemerides (such as the JPL ephemerides via the CALCEPH or CSPICE plugins) and novas_sky_pos() instead.

While this function is generally similar to creating an orbital object with an orbit initialized with novas_make_planet_orbit() or novas_make_moon_orbit(), and then calling novas_sky_pos(), there are a few important differences to note:

  1. This function calculates Earth and Moon positions about the Keplerian orbital position of the Earth-Moon Barycenter (EMB). In constrast, novas_make_planet_orbit() does not provide orbitals for the Earth directly, and make_moot_orbit() references the Moon's orbital to the Earth position returned by the currently configured planet calculator function (see set_planet_provider()).
  2. This function ignores gravitational deflection. It makes little sense to bother about corrections that are orders of magnitude below the accuracy of the orbital positions obtained.

REFERENCES:

  1. E.M. Standish and J.G. Williams 1992.
  2. https://ssd.jpl.nasa.gov/planets/approx_pos.html
  3. Chapront, J. et al., 2002, A&A 387, 700–709
  4. Chapront-Touze, M, and Chapront, J. 1983, Astronomy and Astrophysics (ISSN 0004-6361), vol. 124, no. 1, July 1983, p. 50-62.
Parameters
idNOVAS major planet ID. All major planets, plus the Sun, Moon, Earth-Moon Barycenter (EMB), and Pluto system Barycenter are supported. (For Pluto, the Pluto System Barycenter values are returned.)
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 approximate apparent location in the designated coordinate system.
Returns
0 if successful, or else -1 in case of an error (errno will indicate the type of error).
Since
1.4
Author
Attila Kovacs
See also
novas_sky_pos(), novas_app_to_hor()
novas_make_frame()

References make_planet(), novas_approx_heliocentric(), novas_geom_to_app(), novas_get_time(), NOVAS_TDB, novas_vlen(), and rad_vel2().

◆ novas_equ_sep()

double novas_equ_sep ( double ra1,
double dec1,
double ra2,
double dec2 )

Returns the angular separation of two equatorial locations on a sphere.

Parameters
ra1[h] right ascension of first location
dec1[deg] declination of first location
ra2[h] right ascension of second location
dec2[deg] declination of second location
Returns
[deg] the angular separation of the two locations.
Since
1.3
Author
Attila Kovacs
See also
novas_sep(), novas_sun_angle(), novas_moon_angle()

References novas_sep().

◆ 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_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_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_moon_phase()

double novas_moon_phase ( double jd_tdb)

Calculates the Moon's phase at a given time.

It uses orbital models for Earth (E.M. Standish and J.G. Williams 1992), and for the Moon (Chapront, J. et al., 2002), and takes into account the slightly eccentric nature of both orbits.

NOTES:

  1. The Moon's phase here follows the definition by the Astronomical Almanac, as the excess ecliptic longitude of the Moon over that of the Sun seen from the geocenter.
  2. There are other definitions of the phase too, depending on which you might find slightly different answers, but regardless of the details most phase calculations should match to within a few degrees.

REFERENCES:

  1. The Explanatory Supplement to the Astronomical Almanac, University Science Books, 3rd ed., p. 507
  2. E.M. Standish and J.G. Williams 1992.
  3. https://ssd.jpl.nasa.gov/planets/approx_pos.html
  4. Chapront, J. et al., 2002, A&A 387, 700–709
  5. Chapront-Touze, M, and Chapront, J. 1983, Astronomy and Astrophysics (ISSN 0004-6361), vol. 124, no. 1, July 1983, p. 50-62.
Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian Date.
Returns
[deg] The Moon's phase, or more precisely the ecliptic longitude difference between the Sun and the Moon, as seen from the geocenter. 0: New Moon, 90: 1st quarter, +/- 180 Full Moon, -90: 3rd quarter or NAN if the solution failed to converge (errno will be set to ECANCELED), or if the JD date is outside the range of the orbital model (errno set to EINVAL).
Since
1.4
Author
Attila Kovacs
See also
novas_next_moon_phase(), novas_make_moon_orbit(), novas_solar_illum()

References NOVAS_EMB, novas_make_moon_orbit(), novas_make_planet_orbit(), NOVAS_ORBIT_INIT, novas_orbit_native_posvel(), and vector2radec().

◆ 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_sep()

double novas_sep ( double lon1,
double lat1,
double lon2,
double lat2 )

Returns the angular separation of two locations on a sphere.

Parameters
lon1[deg] longitude of first location
lat1[deg] latitude of first location
lon2[deg] longitude of second location
lat2[deg] latitude of second location
Returns
[deg] the angular separation of the two locations.
Since
1.3
Author
Attila Kovacs
See also
novas_equ_sep(), novas_sun_angle(), novas_moon_angle()

◆ 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().

◆ place_cirs()

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

Computes the Celestial Intermediate Reference System (CIRS) dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

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.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated apparent CIRS position data for a fictitous geocentric observer.
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
Since
1.0
Author
Attila Kovacs
See also
place_tod(), place_gcrs()
novas_sky_pos(), NOVAS_CIRS, make_observer_at_geocenter()

References NOVAS_CIRS, and place().

◆ place_gcrs()

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

Computes the Geocentric Celestial Reference System (GCRS) position of a source (as 'seen' from the geocenter) at the given time of observation. Unlike place_icrs(), this includes aberration for the moving frame of the geocenter as well as gravitational deflections calculated for a virtual observer located at the geocenter. See place() for more information.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated apparent GCRS position data for a fictitous geocentric observer.
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
Since
1.0
Author
Attila Kovacs
See also
place_icrs(), place_cirs(), place_tod(), virtual_star(), virtual_planet()
novas_sky_pos(), NOVAS_GCRS, make_observer_at_geocenter(),

References NOVAS_GCRS, and place().

◆ place_j2000()

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

Computes the J2000 dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

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.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated apparent J2000 position data for a fictitous geocentric observer.
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
Since
1.1
Author
Attila Kovacs
See also
place_cirs(), place_gcrs(), app_star(), app_planet()
novas_sky_pos(), NOVAS_TOD, make_observer_at_geocenter()

References NOVAS_J2000, and place().

◆ place_mod()

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

Computes the Mean of Date (MOD) dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

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.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated apparent Mean-of-Date (MOD) position data for a fictitous geocentric observer.
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
Since
1.1
Author
Attila Kovacs
See also
place_cirs(), place_gcrs(), app_star(), app_planet()
novas_sky_pos(), NOVAS_TOD, make_observer_at_geocenter()

References NOVAS_MOD, and place().

◆ place_tod()

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

Computes the True of Date (TOD) dynamical position position of a source as 'seen' from the geocenter at the given time of observation. See place() for more information.

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.

Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date of observation.
sourceCatalog source or solar_system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]posStructure to populate with the calculated apparent True-of-Date (TOD) position data for a fictitous geocentric observer.
Returns
0 if successful, or -1 if any of the input pointer arguments is NULL, or else an error from place().
Since
1.0
Author
Attila Kovacs
See also
place_cirs(), place_gcrs(), app_star(), app_planet()
novas_sky_pos(), NOVAS_TOD, make_observer_at_geocenter()

References NOVAS_TOD, and place().

◆ virtual_planet()

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

Computes the virtual place of a solar system body, as would be seen by an observer at the geocenter, referenced to the GCRS. This is the same as calling place_gcrs() for the body, 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:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Explanatory Supplement to the Astronomical Almanac (1992),Chapter 3.
Parameters
jd_tt[day] Terretrial Time (TT) based Julian date.
ss_bodyPointer to structure containing the body designation for the solar system body.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]ra[h] Virtual apparent right ascension for a fictitous geocentric observer, referred to the GCRS (it may be NULL if not required).
[out]dec[deg] Virtual apparent declination, for a fictitous geocentric observer, 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 needed).
Returns
0 if successful, or -1 if the object argument is NULL, or else 1 if the value of 'type' in structure 'ss_body' is invalid, or 10 + the error from place().
See also
place_gcrs(), app_planet(), astro_planet(), local_planet(), topo_planet(), app_star()
novas_sky_pos(), make_observer_at_geocenter(), NOVAS_GCRS

References NOVAS_GCRS, and radec_planet().

◆ virtual_star()

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

Computes the virtual place of a star, as would be seen by an observer at the geocenter, referenced to GCRS, at date 'jd_tt', 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 a NULL observer location and NOVAS_GCRS as the system, or place_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:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
  2. Explanatory Supplement to the Astronomical Almanac (1992), Chapter 3.
Parameters
jd_tt[day] Terrestrial Time (TT) based Julian date.
starPointer to catalog entry structure containing catalog data for the object in the ICRS.
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]ra[h] Virtual apparent right ascension for a fictitous geocentric observer, referred to the GCRS (it may be NULL if not required).
[out]dec[deg] Virtual apparent declination for a fictitous geocentric observer, referred to the GCRS (it may be NULL if not required).
Returns
0 if successful, or -1 if a required pointer argument is NULL, or 20 + the error from place().
See also
place_star(), place_gcrs(), app_star(), astro_star(), local_star(), topo_star(), virtual_planet()
novas_sky_pos(), make_observer_at_geocenter(), NOVAS_GCRS

References NOVAS_GCRS, and radec_star().

Variable Documentation

◆ grav_bodies_full_accuracy

int grav_bodies_full_accuracy
extern

Current set of gravitating bodies to use for deflection calculations in full accuracy mode.

Each bit signifies whether a given body is to be accounted for as a gravitating body that bends light, such as the bit (1 << NOVAS_JUPITER) indicates whether or not Jupiter is considered as a deflecting body. You should also be sure that you provide ephemeris data for bodies that are designated for the deflection calculation.

Since
1.1
See also
grav_def(), grav_planets(), DEFAULT_GRAV_BODIES_FULL_ACCURACY
set_ephem_provider_hp()

◆ grav_bodies_reduced_accuracy

int grav_bodies_reduced_accuracy
extern

Current set of gravitating bodies to use for deflection calculations in reduced accuracy mode.

Each bit signifies whether a given body is to be accounted for as a gravitating body that bends light, such as the bit (1 << NOVAS_JUPITER) indicates whether or not Jupiter is considered as a deflecting body. You should also be sure that you provide ephemeris data for bodies that are designated for the deflection calculation.

Since
1.1
See also
grav_def(), grav_planets(), DEFAULT_GRAV_BODIES_REDUCED_ACCURACY
set_ephem_provider()