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

Functions

int aberration (const double *pos, const double *vobs, double lighttime, double *out)
 
int bary2obs (const double *pos, const double *pos_obs, double *out, double *restrict lighttime)
 
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)
 
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)
 
int make_airborne_observer (const on_surface *location, const double *vel, observer *obs)
 
int make_gps_observer (double latitude, double longitude, double height, observer *obs)
 
int make_gps_site (double latitude, double longitude, double height, on_surface *site)
 
int make_in_space (const double *sc_pos, const double *sc_vel, in_space *loc)
 
int make_itrf_observer (double latitude, double longitude, double height, observer *obs)
 
int make_itrf_site (double latitude, double longitude, double height, on_surface *site)
 
short make_observer (enum novas_observer_place where, const on_surface *loc_surface, const in_space *loc_space, observer *obs)
 
int make_observer_at_geocenter (observer *restrict obs)
 
int make_observer_at_site (const on_surface *restrict site, observer *restrict obs)
 
int make_observer_in_space (const double *sc_pos, const double *sc_vel, observer *obs)
 
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_solar_system_observer (const double *sc_pos, const double *sc_vel, observer *obs)
 
int make_xyz_site (const double *restrict xyz, on_surface *restrict site)
 
int novas_e2h_offset (double dra, double ddec, double pa, double *restrict daz, double *restrict del)
 
double novas_epa (double ha, double dec, double lat)
 
int novas_h2e_offset (double daz, double del, double pa, double *restrict dra, double *restrict ddec)
 
double novas_hpa (double az, double el, double lat)
 
int novas_los_to_xyz (const double *los, double lon, double lat, double *xyz)
 
int novas_set_default_weather (on_surface *site)
 
int novas_uvw_to_xyz (const double *uvw, double ha, double dec, double *xyz)
 
int novas_xyz_to_los (const double *xyz, double lon, double lat, double *los)
 
int novas_xyz_to_uvw (const double *xyz, double ha, double dec, double *uvw)
 
int obs_planets (double jd_tdb, enum novas_accuracy accuracy, const double *restrict pos_obs, int pl_mask, novas_planet_bundle *restrict planets)
 
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)
 

Detailed Description

Date
Created on Mar 6, 2025
Author
Attila Kovacs and G. Kaplan

Function that define an astronomical observer location or are related to observer location.

Function Documentation

◆ aberration()

int aberration ( const double * pos,
const double * vobs,
double lighttime,
double * out )

Corrects position vector for aberration of light. Algorithm includes relativistic terms.

NOTES:

  1. This function is called by place() to account for aberration when calculating the position of the source.

REFERENCES:

  1. Murray, C. A. (1981) Mon. Notices Royal Ast. Society 195, 639-648.
  2. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
Parameters
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.
Returns
0 if successful, or -1 if any of the vector arguments are NULL.
See also
frame_aberration()

References novas_vlen().

◆ bary2obs()

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:

  1. Kaplan, G. H. et. al. (1989). Astron. Journ. 97, 1197-1210.
Parameters
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.
Returns
0 if successful, or -1 if any of the essential pointer arguments is NULL.
See also
novas_make_frame(), light_time2()

References novas_vlen().

◆ light_time()

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.

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date
bodyPointer to structure containing the designation for the solar system body
pos_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes.
tlight0[day] First approximation to light-time (can be set to 0.0 if not readily available – it will then be computed as needed).
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]pos_src_obs[AU] Position 3-vector of body, with respect to origin at observer (or the geocenter), referred to ICRS axes. It can be the same vector as either of the inputs.
[out]tlight[day] Calculated light time
Returns
0 if successful, -1 if any of the poiinter arguments is NULL, 1 if the algorithm failed to converge after 10 iterations, or 10 + the error from ephemeris().
See also
light_time2(), novas_make_frame()

References light_time2().

◆ light_time2()

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. It is effectively the same as the original NOVAS C light_time(), except that this returns the antedated source velocity vector also.

NOTES:

  1. This function is called by novas_geom_posvel() or place() to calculate observed positions, radial velocity, and distance for the time when the observed light originated from the source.
Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date
bodyPointer to structure containing the designation for the solar system body
pos_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes.
tlight0[day] First approximation to light-time (can be set to 0.0 if not readily avasilable – so it will be calculated as needed).
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
[out]p_src_obs[AU] Position 3-vector of body, relative to observer, referred to ICRS axes, components in AU.
[out]v_ssb[AU/day] Velocity 3-vector of body, with respect to the Solar-system barycenter, referred to ICRS axes.
[out]tlight[day] Calculated light time, or NAN when returning with an error code.
Returns
0 if successful, -1 if any of the pointer arguments is NULL or if the output vectors are the same or if they are the same as pos_obs, 1 if the algorithm failed to converge after 10 iterations, or 10 + the error from ephemeris().
See also
light_time()
novas_sky_pos()
Since
1.0
Author
Attila Kovacs

References bary2obs(), ephemeris(), NOVAS_BARYCENTER, NOVAS_FULL_ACCURACY, and novas_inv_max_iter.

◆ make_airborne_observer()

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. Airborne observers have an earth fixed momentary location, defined by longitude, latitude, and altitude, the same was as for a stationary observer on Earth, but are moving relative to the surface, such as in an aircraft or balloon observatory.

Parameters
locationCurrent geodetic location, e.g. as populated with make_gps_site() or similar.
vel[km/s] Surface velocity.
[out]obsPointer to data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_itrf_site(), make_gps_site(), make_xyz_site()
make_gps_observer(), make_itrf_observer(), make_observer_at_site(), make_observer_in_space(), make_solar_system_observer(), make_observer_at geocenter(), novas_make_frame()
Since
1.1
Author
Attila Kovacs

References IN_SPACE_INIT, make_observer(), NOVAS_AIRBORNE_OBSERVER, and in_space::sc_vel.

◆ make_gps_observer()

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.

For the highest (μas / mm level) precision, you probably should use an ITRF location instead of a GPS based location.

Parameters
latitude[deg] Geodetic (GPS / WGS84) latitude north positive.
longitude[deg] Geodetic (GPS / WGS84) longitude east positive.
height[m] Geodetic (GPS / WGS84) altitude above sea level of the observer.
[out]obsPointer to the data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL (errno set to EINVAL), or if the latitude is outside of the [-90:90] range (errno set to ERANGE).
Since
1.5
Author
Attila Kovacs
See also
make_itrf_observer(), make_observer_at_site(), make_airborne_observer(), make_observer_in_space(), make_observer_at_geocenter(), make_solar_system_observer()
make_gps_site(), novas_set_default_weather(), novas_make_frame(), novas_geodetic_transform_site()

References make_gps_site(), and make_observer_at_site().

◆ make_gps_site()

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.

For the highest (μas / mm level) precision, you probably should use an ITRF location instead of a GPS based location.

Parameters
latitude[deg] Geodetic (GPS / WGS84) latitude; north positive.
longitude[deg] Geodetic (GPS / WGS84) longitude; east positive.
height[m] Geodetic (GPS / WGS84) altitude above sea level of the observer.
[out]sitePointer to the data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL (errno set to EINVAL), or if the latitude is outside of the [-90:90] range (errno set to ERANGE).
Since
1.5
Author
Attila Kovacs
See also
make_gps_site(), make_xyz_site()
make_itrf_observer(), make_observer_at_site(), novas_set_default_weather()

References make_itrf_site(), novas_cartesian_to_geodetic(), novas_geodetic_to_cartesian(), NOVAS_GRS80_ELLIPSOID, and NOVAS_WGS84_ELLIPSOID.

◆ make_in_space()

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. Both input vectors are assumed with respect to true equator and equinox of date.

Parameters
sc_pos[km] Geocentric (x, y, z) position vector. NULL defaults to the origin
sc_vel[km/s] Geocentric (x, y, z) velocity vector. NULL defaults to zero speed.
[out]locPointer to earth-orbit location data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_observer_in_space(), IN_SPACE_INIT

References in_space::sc_pos, and in_space::sc_vel.

◆ make_itrf_observer()

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.

For the highest precision (μas level) applications you should make sure that the location provided here and the Earth-orientation parameters (EOP) used (in setting novas_timescale and in novas_make_frame()) are provided in the same ITRF realization. You can use novas_itrf_transform_eop() to change the ITRF realization for the EOP values, if necessary.

Parameters
latitude[deg] Geodetic (ITRF / GRS80) latitude; north positive.
longitude[deg] Geodetic (ITRF / GRS80) longitude; east positive.
height[m] Geodetic (ITRF / GRS80) altitude above sea level of the observer.
[out]obsPointer to the data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL (errno set to EINVAL), or if the latitude is outside of the [-90:90] range (errno set to ERANGE).
Since
1.5
Author
Attila Kovacs
See also
make_gps_observer(), make_observer_at_site(), make_airborne_observer(), make_observer_in_space(), make_observer_at_geocenter(), make_solar_system_observer()
make_itrf_site(), novas_itrf_transform_site(), novas_set_default_weather(), novas_make_frame()

References make_itrf_site(), and make_observer_at_site().

◆ make_itrf_site()

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.

For the highest precision (μas level) applications you should make sure that the location provided here and the Earth-orientation parameters (EOP) used (in setting novas_timescale and in novas_make_frame()) are provided in the same ITRF realization. You can use novas_itrf_transform_eop() to change the ITRF realization for the EOP values, if necessary.

Parameters
latitude[deg] Geodetic (ITRF / GRS80) latitude; north positive.
longitude[deg] Geodetic (ITRF / GRS80) longitude; east positive.
height[m] Geodetic (ITRF / GRS80) altitude above sea level of the observer.
[out]sitePointer to the data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL (errno set to EINVAL), or if the latitude is outside of the [-90:90] range (errno set to ERANGE).
Since
1.5
Author
Attila Kovacs
See also
make_gps_site(), make_xyz_site()
make_itrf_observer(), make_observer_at_site(), novas_set_default_weather(), novas_itrf_transform_site(), novas_itrf_transform_eop()

References on_surface::height, on_surface::latitude, on_surface::longitude, and novas_set_default_weather().

◆ make_observer()

short make_observer ( enum novas_observer_place where,
const on_surface * loc_surface,
const in_space * loc_space,
observer * obs )
Deprecated
It is recommended that you use one of the more specific ways of initializing the observer data structure, e.g. with make_itrf_observer(), make_gps_observer(), make_observer_at_site(), make_airborne_observer() make_solar_system_observer(), ormake_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().

Parameters
whereThe location type of the observer
loc_surfacePointer 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_spacePointer 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]obsPointer to observer data structure to populate.
Returns
0 if successful, -1 if a required argument is NULL, or 1 if the 'where' argument is invalid.
See also
make_observer_at_geocenter(), make_itrf_observer(), make_gps_observer(), make_airborne_observer(), make_observer_in_space(), make_solar_system_observer()

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.

◆ make_observer_at_geocenter()

int make_observer_at_geocenter ( observer *restrict obs)

Populates an 'observer' data structure for a hypothetical observer located at Earth's geocenter. The output data structure may be used an the the inputs to NOVAS-C functions, such as make_frame() or place().

Parameters
[out]obsPointer to data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_gps_observer(), make_itrf_observer(), make_observer_at_site(), make_airborne_observer(), make_observer_in_space(), make_solar_system_observer()
novas_make_frame()

References make_observer(), and NOVAS_OBSERVER_AT_GEOCENTER.

◆ make_observer_at_site()

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.

Parameters
sitePointer to observing site data structure.
[out]obsPointer to the data structure to populate.
Returns
0 if successful, or -1 if the either argument is NULL (errno set to EINVAL).
Since
1.5
Author
Attila Kovacs
See also
make_itrf_observer(), make_gps_observer(), make_airborne_observer(), make_observer_in_space(), make_observer_at_geocenter(), make_solar_system_observer()
make_itrf_site(), make_gps_site(), make_xyz_site(), novas_make_frame()

References NOVAS_OBSERVER_ON_EARTH.

◆ make_observer_in_space()

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. Both input vectors are with respect to true equator and equinox of date. The output data structure may be used an the the inputs to NOVAS-C functions, such as make_frame() or place().

Parameters
sc_pos[km] Geocentric (x, y, z) position vector.
sc_vel[km/s] Geocentric (x, y, z) velocity vector.
[out]obsPointer to the data structure to populate
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_gps_observer(), make_itrf_observer(), make_observer_at_site(), make_airborne_observer(), make_observer_at_geocenter(), make_solar_system_observer()
novas_make_frame()

References make_in_space(), make_observer(), and NOVAS_OBSERVER_IN_EARTH_ORBIT.

◆ make_observer_on_surface()

int make_observer_on_surface ( double latitude,
double longitude,
double height,
double temperature,
double pressure,
observer *restrict obs )
Deprecated
This old NOVAS C function has a few too many caveats. It is recommended that you use make_itrf_observer(), make_gps_observer(), or make_observer_at_site() instead (all 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.

Initializes an observer data structure with the specified ITRF / GRS80 location, and he specified local pressure and temperature. This old NOVAS C function does not set humidity. Thus, if humidity is needed for refraction correction, then you should set it explicitly after this call.

NOTES:

  1. You can convert coordinates among ITRF realization using novas_itrf_transform(), possibly after novas_geodetic_to_cartesian() with NOVAS_GRS80_ELLIPSOID as necessary for polar ITRF coordinates.
  2. You can convert ITRF Cartesian xyz locations to geodetic locations by using novas_cartesian_to_geodetic() with NOVAS_GRS80_ELLIPSOID as the reference ellipsoid parameter.
  3. If you have longitude, latitude, and height defined as GPS (WGS84) values, you might want to use make_gps_observer() intead
Parameters
latitude[deg] Geodetic (ITRF / GRS80) latitude; north positive.
longitude[deg] Geodetic (ITRF / GRS80) longitude; east positive.
height[m] Geodetic (ITRF / GRS80) altitude above sea level of the observer.
temperature[C] Temperature (degrees Celsius).
pressure[mbar] Atmospheric pressure.
[out]obsPointer to the data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL, or if the latitude is outside of the [-90:90] range, or if the temperature or pressure values are impossible for an Earth based observer (errno set to ERANGE).
See also
make_itrf_observer(), make_gps_observer(), make_observer_at_site()
novas_set_weather(), novas_cartesian_to_geodetic(), novas_geodetic_to_cartesian(), NOVAS_GRS80_ELLIPSOID, novas_make_frame()

References make_observer(), make_on_surface(), and NOVAS_OBSERVER_ON_EARTH.

◆ make_on_surface()

int make_on_surface ( double latitude,
double longitude,
double height,
double temperature,
double pressure,
on_surface *restrict loc )
Deprecated
This old NOVAS C function has a few too many caveats. It is recommended that you use 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

  1. This implementation breaks strict v1.0 ABI compatibility since it writes to (initializes) a field (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.
  2. You can convert coordinates among ITRF realization using novas_itrf_transform(), possibly after novas_geodetic_to_cartesian() with NOVAS_GRS80_ELLIPSOID as necessary for polar ITRF coordinates.
  3. You can convert Cartesian xyz locations to geodetic locations by using novas_cartesian_to_geodetic() with NOVAS_GRS80_ELLIPSOID as the reference ellipsoid parameter.
  4. If you have longitude, latitude, and height defined as GPS (WGS84) values, you might want to use make_gps_site() instead, and then set weather parameters afterwards as necessary.
Parameters
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]locPointer to Earth location data structure to populate.
Returns
0 if successful, or -1 if the output argument is NULL (errno set to EINVAL), or if the latitude is outside of the [-90:90] range, or if the temperature or pressure values are impossible for an Earth based observer (errno set to ERANGE).
See also
make_itrf_site(), make_gps_site(), make_xyz_site()
novas_set_default_weather(), novas_cartesian_to_geodetic(), ON_SURFACE_INIT, ON_SURFACE_LOC, NOVAS_GRS80_ELLIPSOID

References make_itrf_site().

◆ make_solar_system_observer()

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. Solar-system observers are similar to observers in Earth-orbit but their momentary position and velocity is defined relative to the Solar System Barycenter, instead of the geocenter.

Parameters
sc_pos[AU] Solar-system barycentric (x, y, z) position vector in ICRS.
sc_vel[AU/day] Solar-system barycentric (x, y, z) velocity vector in ICRS.
[out]obsPointer to the data structure to populate
Returns
0 if successful, or -1 if the output argument is NULL.
See also
make_gps_observer(), make_itrf_observer(), make_observer_at_site(), make_airborne_observer(), make_observer_in_space(), make_solar_system_observer()
novas_make_frame()
Since
1.1
Author
Attila Kovacs

References make_in_space(), make_observer(), and NOVAS_SOLAR_SYSTEM_OBSERVER.

◆ make_xyz_site()

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.

For the highest precision (μas level) applications you should make sure that the site coordinates and the Earth-orientation parameters (EOP) used (in setting novas_timescale and in novas_make_frame()) are provided in the same ITRF realization. You can use novas_itrf_transform() or novas_itrf_transform_eop() to change the ITRF realization for the site coordinates and/or the EOP values, if necessary.

Parameters
xyz[m] Cartesian geocentric position.
[out]sitePointer to the data structure to populate.
Returns
0 if successful, or -1 if either of the arguments is NULL (errno set to EINVAL), or if the latitude is outside of the [-90:90] range (errno set to ERANGE).
Since
1.5
Author
Attila Kovacs
See also
make_gps_site(), make_itrf_site(), make_xyz_site()
make_observer_at_site(), novas_itrf_transform(), novas_itrf_transform_eop(), novas_set_default_weather()

References make_itrf_site(), novas_cartesian_to_geodetic(), and NOVAS_GRS80_ELLIPSOID.

◆ novas_e2h_offset()

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. Converting between local flat projections and spherical coordinates usually requires a WCS projection.

REFERENCES:

  1. Calabretta, M.R., & Greisen, E.W., (2002), Astronomy & Astrophysics, 395, 1077-1122.
Parameters
dra[arcsec] Projected ffset position in the apparent true-of-date R.A. direction. E.g. The projected offset between two RA coordinates at a same reference declination, is δRA = (RA2 - RA1) * cos(Dec0).
ddec[arcsec] Projected offset position in the apparent true-of-date declination direction.
pa[deg] Parallactic Angle
[out]daz[arcsec] Output offset position in the local azimuth direction. It can be a pointer to one of the input coordinates, or NULL if not required.
[out]del[arcsec] Output offset position in the local elevation direction. It can be a pointer to one of the input coordinates, or NULL if not required.
Returns
0
Since
1.3
Author
Attila Kovacs
See also
novas_h2e_offset(), novas_epa()

References novas_h2e_offset().

◆ novas_epa()

double novas_epa ( double ha,
double dec,
double lat )

Returns the equatorial Parallactic Angle (PA) calculated for an R.A./Dec location of the sky at a given sidereal time. The PA is the angle between the local horizontal coordinate directions and the local true-of-date equatorial coordinate directions, at the given location and time. The polar wobble is not included in the calculation.

The Parallactic Angle is sometimes referrred to as the Vertical Position Angle (VPA). Both define the same quantity.

Parameters
ha[h] Hour angle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source.
dec[deg] Apparent (true-of-date) declination of observed source
lat[deg] Geodetic latitude of observer
Returns
[deg] Parallactic Angle (PA). I.e., the clockwise position angle of the elevation direction w.r.t. the declination axis in the equatorial system. Same as the clockwise position angle of the declination direction w.r.t. the elevation axis, in the horizontal system.
Since
1.3
Author
Attila Kovacs
See also
novas_hpa(), novas_e2h_offset()

◆ novas_h2e_offset()

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. Converting between local flat projections and spherical coordinates usually requires a WCS projection.

REFERENCES:

  1. Calabretta, M.R., & Greisen, E.W., (2002), Astronomy & Astrophysics, 395, 1077-1122.
Parameters
daz[arcsec] Projected offset position in the azimuth direction. The projected offset between two azimuth positions at the same reference elevation is δAz = (Az2 - Az1) * cos(El0).
del[arcsec] projected offset position in the elevation direction
pa[deg] Parallactic Angle
[out]dra[arcsec] Output offset position in the local true-of-date R.A. direction. It can be a pointer to one of the input coordinates, or NULL if not required.
[out]ddec[arcsec] Output offset position in the local true-of-date declination direction. It can be a pointer to one of the input coordinates, or NULL if not required.
Returns
0
Since
1.3
Author
Attila Kovacs
See also
novas_e2h_offset(), novas_hpa()

◆ novas_hpa()

double novas_hpa ( double az,
double el,
double lat )

Returns the horizontal Parallactic Angle (PA) calculated for a gorizontal Az/El location of the sky. The PA is the angle between the local horizontal coordinate directions and the local true-of-date equatorial coordinate directions at the given location. The polar wobble is not included in the calculation.

The Parallactic Angle is sometimes referrred to as the Vertical Position Angle (VPA). Both define the same quantity.

Parameters
az[deg] Azimuth angle
el[deg] Elevation angle
lat[deg] Geodetic latitude of observer
Returns
[deg] Parallactic Angle (PA). I.e., the clockwise position angle of the declination direction w.r.t. the elevation axis in the horizontal system. Same as the the clockwise position angle of the elevation direction w.r.t. the declination axis in the equatorial system.
Since
1.3
Author
Attila Kovacs
See also
novas_epa(), novas_h2e_offset()

◆ novas_los_to_xyz()

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.

Parameters
los[arb.u.] Line-of-sight 3-vector (δφ, δθ δr).
lon[deg] Line-of-sight longitude.
lat[deg] Line-of-sight latitude.
[out]xyz[arb.u.] Output rectangular equatorial 3-vector (δx, δy, δz), in the same units as the input. It may be the same vector as the input.
Returns
0 if successful, or else -1 if either vector argument is NULL (errno will be set to EINVAL).
Since
1.3
Author
Attila Kovacs
See also
novas_xyz_to_los(), novas_uvw_to_xyz()

◆ novas_set_default_weather()

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. (2013), and scaling relations with altitude (up to 12 km).

Humidity is set to 70% at sea level (which is a typical value globally), and adjusted to decrease with altitude linearly at a rate of 7.5% per km up to 8000 meters. Above that a quadratic model is assumed, peaking at 45% at 14 km – based on the measured distribution by Mendez-Astudillo et al. (2021). Finally above 20.8 km, zero humidity is assumed.

These parameters are all very approximate, but in the absence of measured data, they represent a best guess default model of sorts.

REFERENCES:

  1. Feulner, G., Rahmstorf, S., Levermann, A., and Volkwardt, S. (2013), Journal of Climate 26, 7136
  2. Mendez-Astudillo, J., et al. (2021), Urban Heat Island (UHI) Mitigation (pp.43-59), DOI:10.1007/978-981-33-4050-3_3
Parameters
[in,out]siteSite containing geodetic loation as input, and poulated with typical mean weather parameters for the output.
Returns
0 if successful, or else -1 if the site is NULL (errno is set to EINVAL).
Since
1.5
Author
Attila Kovacs
See also
make_itrf_site(), make_gps_site(), make_xyz_site(), make_itrf_observer(), make_gps_observer(), NOVAS_STANDARD_ATMOSPHERE

References on_surface::height, on_surface::humidity, on_surface::latitude, on_surface::pressure, and on_surface::temperature.

◆ novas_uvw_to_xyz()

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.

u,v,w are Cartesian coordinates (u,v) along the local equatorial R.A. and declination directions as seen from a direction on the sky (w). As such, they are effectively ITRS-based line-of-sight (LOS) coordinates.

x,y,z are Cartesian coordinates w.r.t the Greenwich meridian in the ITRS frame. The directions are x: long=0, lat=0; y: long=90, lat=0; z: lat=90.

Parameters
xyz[arb.u.] Absolute or relative u,v,w coordinates (double[3]).
ha[h] Hourangle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source.
dec[deg] Apparent (true-of-date) declination of source
[out]uvw[arb.u.] Converted x,y,z coordinates (double[3]) in the same unit as uvw. It may be the same vector as the input.
Returns
0 if successful, or else -1 if either vector argument is NULL (errno will be set to EINVAL)
Since
1.3
Author
Attila Kovacs
See also
novas_xyz_to_uvw()

References novas_los_to_xyz().

◆ novas_xyz_to_los()

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.

Parameters
xyz[arb.u.] Rectangular equatorial 3-vector (δx, δy, δz).
lon[deg] Line-of-sight longitude.
lat[deg] Line-of-sight latitude.
[out]los[arb.u.] Output line-of-sight 3-vector (δφ, δθ δr), in the same units as the input. It may be the same vector as the input.
Returns
0 if successful, or else -1 if either vector argument is NULL (errno will be set to EINVAL).
Since
1.3
Author
Attila Kovacs
See also
novas_los_to_xyz(), novas_xyz_to_uvw()

◆ novas_xyz_to_uvw()

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.

x,y,z are Cartesian coordinates w.r.t the Greenwich meridian, in the ITRS frame. The directions are x: long=0, lat=0; y: long=90, lat=0; z: lat=90.

u,v,w are Cartesian coordinates (u,v) along the local equatorial R.A. and declination directions as seen from a direction on the sky (w). As such, they are effectively ITRS-based line-of-sight (LOS) coordinates.

Parameters
xyz[arb.u.] Absolute or relative x,y,z coordinates (double[3]).
ha[h] Hourangle (LST - RA) i.e., the difference between the Local (apparent) Sidereal Time and the apparent (true-of-date) Right Ascension of observed source.
dec[deg] Apparent (true-of-date) declination of source
[out]uvw[arb.u.] Converted u,v,w coordinates (double[3]) in same units as xyz. It may be the same vector as the input.
Returns
0 if successful, or else -1 if either vector argument is NULL (errno will be set to EINVAL)
Since
1.3
Author
Attila Kovacs
See also
novas_uvw_to_xyz()

References novas_xyz_to_los().

◆ obs_planets()

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. for use for gravitational deflection calculations. The planet positions are calculated relative to the observer location, while velocities are w.r.t. the SSB. Both positions and velocities are antedated for light travel time, so they accurately reflect the apparent position (and barycentric motion) of the bodies from the observer's perspective.

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date
accuracyNOVAS_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_obs[AU] Position 3-vector of observer (or the geocenter), with respect to origin at solar system barycenter, referred to ICRS axes.
pl_maskBitwise (1 << planet-number) mask indicating which planets to request data for. See enum novas_planet for the enumeration of planet numbers.
[out]planetsPointer to apparent planet data to populate. The planets with non-zero mask bits will have have positions and velocities calculated. See enum novas_planet for the enumeration of planet numbers.
Returns
0 if successful, -1 if any of the pointer arguments is NULL or if the output vector is the same as pos_obs, or the error from ephemeris().
See also
enum novas_planet, grav_planets(), grav_undo_planets(), set_planet_provider(), set_planet_provider_hp()
Since
1.1
Author
Attila Kovacs

References light_time2(), make_planet(), novas_debug(), NOVAS_DEBUG_EXTRA, NOVAS_DEBUG_OFF, novas_get_debug_mode(), NOVAS_PLANETS, and NOVAS_SUN.

◆ obs_posvel()

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

Parameters
jd_tdb[day] Barycentric Dynamical Time (TDB) based Julian date.
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).
accuracyNOVAS_FULL_ACCURACY (0) or NOVAS_REDUCED_ACCURACY (1)
obsThe observer location, relative to which the output positions and velocities are to be calculated
geo_pos[AU] ICRS position vector of the geocenter w.r.t. the Solar System Barycenter (SSB). If either geo_pos or geo_vel is NULL, it will be calculated when needed.
geo_vel[AU/day] ICRS velocity vector of the geocenter w.r.t. the Solar System Barycenter (SSB). If either geo_pos or geo_vel is NULL, it will be calculated when needed.
[out]pos[AU] Position 3-vector of the observer w.r.t. the Solar System Barycenter (SSB). It may be NULL if not required.
[out]vel[AU/day] Velocity 3-vector of the observer w.r.t. the Solar System Barycenter (SSB). It must be distinct from the pos output vector, and may be NULL if not required.
Returns
0 if successful, or the error from geo_posvel(), or else -1 (with errno indicating the type of error).
Author
Attila Kovacs
Since
1.3
See also
novas_make_frame()

References ephemeris(), geo_posvel(), NOVAS_AIRBORNE_OBSERVER, NOVAS_BARYCENTER, NOVAS_EARTH_INIT, NOVAS_OBSERVER_IN_EARTH_ORBIT, NOVAS_OBSERVER_ON_EARTH, NOVAS_OBSERVER_PLACES, NOVAS_SOLAR_SYSTEM_OBSERVER, and tt2tdb().